v0.13.0
Public Member Functions | Public Attributes | List of all members
PlaneRigidBody Struct Reference

#include <users_modules/multifield_plasticity/src/RigidBodies.hpp>

Inheritance diagram for PlaneRigidBody:
[legend]
Collaboration diagram for PlaneRigidBody:
[legend]

Public Member Functions

 PlaneRigidBody (VectorDouble c_coords, VectorDouble roller_disp, int id)
 
MoFEMErrorCode getRollerDataForTag ()
 
Tensor1< double, 3 > getNormal (Tensor1< TPack3, 3 > &t_coords, Tensor1< TPack1, 3 > &t_disp)
 
Tensor1< double, 3 > getNormal (Tensor1< TPack3, 3 > &t_coords, Tensor1< double, 3 > &t_disp)
 
Tensor1< double, 3 > getNormal (Tensor1< double, 3 > &t_coords, Tensor1< double, 3 > &t_disp)
 
Tensor2< double, 3, 3 > getDiffNormal (Tensor1< TPack3, 3 > &t_coords, Tensor1< TPack1, 3 > &t_disp, Tensor1< TPack1, 3 > &t_normal)
 
double getGap (Tensor1< TPack3, 3 > &t_coords)
 
Tensor1< double, 3 > getdGap (Tensor1< TPack3, 3 > &t_coords, Tensor1< TPack1, 3 > &t_normal)
 
template<typename T1 , typename T2 >
Tensor1< double, 3 > & getNormalImpl (Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_disp)
 
template<typename T1 >
double getGapImpl (Tensor1< T1, 3 > &t_coords)
 
template<typename T1 , typename T2 , typename T3 >
Tensor2< double, 3, 3 > getDiffNormalImpl (Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_disp, Tensor1< T3, 3 > &t_normal)
 
template<typename T1 , typename T2 >
Tensor1< double, 3 > getdGapImpl (Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_normal)
 
MoFEMErrorCode getBodyOptions ()
 
- Public Member Functions inherited from RigidBodyData
 RigidBodyData (VectorDouble c_coords, VectorDouble roller_disp, int id)
 
 RigidBodyData ()=delete
 
virtual ~RigidBodyData ()
 
MoFEMErrorCode computeRotationMatrix ()
 
Tensor1< double, 3 > getBodyOffset ()
 
template<typename T1 , typename T2 >
Tensor1< double, 3 > & getPointCoords (Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_disp)
 
MoFEMErrorCode saveBasicDataOnTag (moab::Interface &moab_debug, EntityHandle &vertex)
 

Public Attributes

PetscBool fLip
 
- Public Attributes inherited from RigidBodyData
const int iD
 
Index< 'i', 3 > i
 
Index< 'j', 3 > j
 
Index< 'k', 3 > k
 
VectorDouble originCoords
 
VectorDouble rollerDisp
 
VectorDouble BodyDispScaled
 
double gAp
 
Tensor1< double, 3 > tNormal
 
Tensor1< double, 3 > dGap
 
Tensor1< double, 3 > pointCoords
 
Tensor1< double, 3 > closestPoint
 
Tensor1< double, 3 > defaultOrientation
 
Tensor1< double, 3 > oRientation
 
Tensor2< double, 3, 3 > rotationMat
 
Tensor2< double, 3, 3 > diffNormal
 
string positionDataParamName
 
int bodyType
 
array< double, 9 > dataForTags
 

Additional Inherited Members

- Public Types inherited from RigidBodyData
using TPack3 = PackPtr< double *, 3 >
 
using TPack1 = PackPtr< double *, 1 >
 

Detailed Description

Definition at line 164 of file RigidBodies.hpp.

Constructor & Destructor Documentation

◆ PlaneRigidBody()

PlaneRigidBody::PlaneRigidBody ( VectorDouble  c_coords,
VectorDouble  roller_disp,
int  id 
)

Definition at line 166 of file RigidBodies.hpp.

167  : RigidBodyData(c_coords, roller_disp, id), fLip(PETSC_FALSE) {}
PetscBool fLip
RigidBodyData()=delete

Member Function Documentation

◆ getBodyOptions()

MoFEMErrorCode PlaneRigidBody::getBodyOptions ( )
virtual

Implements RigidBodyData.

Definition at line 238 of file RigidBodies.hpp.

238  {
240  PetscBool rflg;
241  PetscInt nb_dirs = 3;
242  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "", "");
243  string param = "-direction" + to_string(iD);
244  CHKERR PetscOptionsGetRealArray(NULL, NULL, param.c_str(), &oRientation(0),
245  &nb_dirs, &rflg);
246  param = "-flip" + to_string(iD);
247  CHKERR PetscOptionsBool(param.c_str(), "flip the plane normal", "", fLip,
248  &fLip, PETSC_NULL);
249  if (rflg)
251  ierr = PetscOptionsEnd();
252  CHKERRG(ierr);
254  }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:87
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
Tensor1< double, 3 > oRientation
Definition: RigidBodies.hpp:43
MoFEMErrorCode computeRotationMatrix()
Definition: RigidBodies.hpp:85
const int iD
Definition: RigidBodies.hpp:30

◆ getdGap()

Tensor1<double, 3> PlaneRigidBody::getdGap ( Tensor1< TPack3, 3 > &  t_coords,
Tensor1< TPack1, 3 > &  t_normal 
)
virtual

Implements RigidBodyData.

Definition at line 198 of file RigidBodies.hpp.

199  {
200  return getdGapImpl(t_coords, t_normal);
201  }
Tensor1< double, 3 > getdGapImpl(Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_normal)

◆ getdGapImpl()

template<typename T1 , typename T2 >
Tensor1<double, 3> PlaneRigidBody::getdGapImpl ( Tensor1< T1, 3 > &  t_coords,
Tensor1< T2, 3 > &  t_normal 
)

Definition at line 232 of file RigidBodies.hpp.

233  {
234  dGap(i) = tNormal(i);
235  return dGap;
236  };
Tensor1< double, 3 > tNormal
Definition: RigidBodies.hpp:38
Tensor1< double, 3 > dGap
Definition: RigidBodies.hpp:39
Index< 'i', 3 > i
Definition: RigidBodies.hpp:31

◆ getDiffNormal()

Tensor2<double, 3, 3> PlaneRigidBody::getDiffNormal ( Tensor1< TPack3, 3 > &  t_coords,
Tensor1< TPack1, 3 > &  t_disp,
Tensor1< TPack1, 3 > &  t_normal 
)
virtual

Implements RigidBodyData.

Definition at line 188 of file RigidBodies.hpp.

190  {
191  return getDiffNormalImpl(t_coords, t_disp, t_normal);
192  }
Tensor2< double, 3, 3 > getDiffNormalImpl(Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_disp, Tensor1< T3, 3 > &t_normal)

◆ getDiffNormalImpl()

template<typename T1 , typename T2 , typename T3 >
Tensor2<double, 3, 3> PlaneRigidBody::getDiffNormalImpl ( Tensor1< T1, 3 > &  t_coords,
Tensor1< T2, 3 > &  t_disp,
Tensor1< T3, 3 > &  t_normal 
)

Definition at line 224 of file RigidBodies.hpp.

226  {
227  const Tensor2<double, 3, 3> t_diff_norm(0., 0., 0., 0., 0., 0., 0., 0., 0.);
228  return t_diff_norm;
229  };

◆ getGap()

double PlaneRigidBody::getGap ( Tensor1< TPack3, 3 > &  t_coords)
virtual

Implements RigidBodyData.

Definition at line 194 of file RigidBodies.hpp.

194  {
195  return getGapImpl(t_coords);
196  }
double getGapImpl(Tensor1< T1, 3 > &t_coords)

◆ getGapImpl()

template<typename T1 >
double PlaneRigidBody::getGapImpl ( Tensor1< T1, 3 > &  t_coords)

Definition at line 216 of file RigidBodies.hpp.

216  {
217  // test(i) = rotationMat(j, i) * t_off(j);
218  // t_d2(i) = t_coords(i) - rotationMat(j, i) * t_off(j);
219  gAp = tNormal(i) * pointCoords(i);
220  return gAp;
221  };
Tensor1< double, 3 > pointCoords
Definition: RigidBodies.hpp:40

◆ getNormal() [1/3]

Tensor1<double, 3> PlaneRigidBody::getNormal ( Tensor1< double, 3 > &  t_coords,
Tensor1< double, 3 > &  t_disp 
)
virtual

Implements RigidBodyData.

Definition at line 183 of file RigidBodies.hpp.

184  {
185  return getNormalImpl(t_coords, t_disp);
186  }
Tensor1< double, 3 > & getNormalImpl(Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_disp)

◆ getNormal() [2/3]

Tensor1<double, 3> PlaneRigidBody::getNormal ( Tensor1< TPack3, 3 > &  t_coords,
Tensor1< double, 3 > &  t_disp 
)
virtual

Implements RigidBodyData.

Definition at line 179 of file RigidBodies.hpp.

180  {
181  return getNormalImpl(t_coords, t_disp);
182  }

◆ getNormal() [3/3]

Tensor1<double, 3> PlaneRigidBody::getNormal ( Tensor1< TPack3, 3 > &  t_coords,
Tensor1< TPack1, 3 > &  t_disp 
)
virtual

Implements RigidBodyData.

Definition at line 175 of file RigidBodies.hpp.

176  {
177  return getNormalImpl(t_coords, t_disp);
178  }

◆ getNormalImpl()

template<typename T1 , typename T2 >
Tensor1<double, 3>& PlaneRigidBody::getNormalImpl ( Tensor1< T1, 3 > &  t_coords,
Tensor1< T2, 3 > &  t_disp 
)

Definition at line 204 of file RigidBodies.hpp.

205  {
206  auto t_d = getPointCoords(t_coords, t_disp);
207  // FIXME:
208  // tNormal(i) = rotationMat(j, i) * defaultOrientation(j);
209  tNormal(i) = oRientation(i);
210  if (fLip)
211  tNormal(i) *= -1;
212  // rotationMat(i, j)
213  return tNormal;
214  };
Tensor1< double, 3 > & getPointCoords(Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_disp)

◆ getRollerDataForTag()

MoFEMErrorCode PlaneRigidBody::getRollerDataForTag ( )
virtual

Implements RigidBodyData.

Definition at line 169 of file RigidBodies.hpp.

169  {
171  bodyType = PLANE;
173  }
@ PLANE
Definition: RigidBodies.hpp:18
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453

Member Data Documentation

◆ fLip

PetscBool PlaneRigidBody::fLip

Definition at line 165 of file RigidBodies.hpp.


The documentation for this struct was generated from the following file: