v0.13.0
Public Types | Public Member Functions | Public Attributes | List of all members
RigidBodyData Struct Referenceabstract

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

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

Public Types

using TPack3 = PackPtr< double *, 3 >
 
using TPack1 = PackPtr< double *, 1 >
 

Public Member Functions

 RigidBodyData (VectorDouble c_coords, VectorDouble roller_disp, int id)
 
 RigidBodyData ()=delete
 
virtual ~RigidBodyData ()
 
virtual Tensor1< double, 3 > getNormal (Tensor1< TPack3, 3 > &t_coords, Tensor1< TPack1, 3 > &t_disp)=0
 
virtual Tensor1< double, 3 > getNormal (Tensor1< TPack3, 3 > &t_coords, Tensor1< double, 3 > &t_disp)=0
 
virtual Tensor1< double, 3 > getNormal (Tensor1< double, 3 > &t_coords, Tensor1< double, 3 > &t_disp)=0
 
virtual Tensor2< double, 3, 3 > getDiffNormal (Tensor1< TPack3, 3 > &t_coords, Tensor1< TPack1, 3 > &t_disp, Tensor1< TPack1, 3 > &t_normal)=0
 
virtual double getGap (Tensor1< TPack3, 3 > &t_coords)=0
 
virtual Tensor1< double, 3 > getdGap (Tensor1< TPack3, 3 > &t_coords, Tensor1< TPack1, 3 > &t_normal)=0
 
virtual MoFEMErrorCode getBodyOptions ()=0
 
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)
 
virtual MoFEMErrorCode getRollerDataForTag ()=0
 

Public Attributes

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
 

Detailed Description

Definition at line 29 of file RigidBodies.hpp.

Member Typedef Documentation

◆ TPack1

using RigidBodyData::TPack1 = PackPtr<double *, 1>

Definition at line 66 of file RigidBodies.hpp.

◆ TPack3

using RigidBodyData::TPack3 = PackPtr<double *, 3>

Definition at line 65 of file RigidBodies.hpp.

Constructor & Destructor Documentation

◆ RigidBodyData() [1/2]

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

Definition at line 52 of file RigidBodies.hpp.

53  : originCoords(c_coords), rollerDisp(roller_disp), iD(id),
54  defaultOrientation(0., 0., 1.), oRientation(0., 0., 1.) {
56  BodyDispScaled.clear();
58  std::fill(dataForTags.begin(), dataForTags.end(), 0);
59  bodyType = PLANE;
60  }
@ PLANE
Definition: RigidBodies.hpp:18
Tensor2_Expr< Kronecker_Delta< T >, T, Dim0, Dim1, i, j > kronecker_delta(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
Rank 2.
VectorDouble rollerDisp
Definition: RigidBodies.hpp:35
Tensor1< double, 3 > oRientation
Definition: RigidBodies.hpp:43
VectorDouble BodyDispScaled
Definition: RigidBodies.hpp:36
Index< 'i', 3 > i
Definition: RigidBodies.hpp:31
Tensor2< double, 3, 3 > rotationMat
Definition: RigidBodies.hpp:44
array< double, 9 > dataForTags
Definition: RigidBodies.hpp:50
VectorDouble originCoords
Definition: RigidBodies.hpp:34
Index< 'j', 3 > j
Definition: RigidBodies.hpp:32
const int iD
Definition: RigidBodies.hpp:30
Tensor1< double, 3 > defaultOrientation
Definition: RigidBodies.hpp:42

◆ RigidBodyData() [2/2]

RigidBodyData::RigidBodyData ( )
delete

◆ ~RigidBodyData()

virtual RigidBodyData::~RigidBodyData ( )
virtual

Definition at line 63 of file RigidBodies.hpp.

63 {}

Member Function Documentation

◆ computeRotationMatrix()

MoFEMErrorCode RigidBodyData::computeRotationMatrix ( )

Definition at line 85 of file RigidBodies.hpp.

85  {
86  // https://math.stackexchange.com/questions/180418
88  oRientation.normalize();
91  const double s = v.l2();
92  if (abs(s) < 1e-16)
94  const double c = defaultOrientation(i) * oRientation(i);
95  const double coef = 1. / (1. + c);
96  Tensor2<double, 3, 3> t_v;
97  t_v(i, j) = levi_civita(i, j, k) * v(k);
98  rotationMat(i, j) += t_v(i, j) + (coef * t_v(i, k)) * t_v(k, j);
99  // Tensor1<double, 3> test;
100  // test(i) = rotationMat(j, i) * defaultOrientation(j);
101 
103  };
#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
double v
phase velocity of light in medium (cm/ns)
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
Index< 'k', 3 > k
Definition: RigidBodies.hpp:33

◆ getBodyOffset()

Tensor1<double, 3> RigidBodyData::getBodyOffset ( )

Definition at line 105 of file RigidBodies.hpp.

105  {
107  return Tensor1<double, 3>(pos(0), pos(1), pos(2));
108  };
UBlasVector< double > VectorDouble
Definition: Types.hpp:79

◆ getBodyOptions()

virtual MoFEMErrorCode RigidBodyData::getBodyOptions ( )
pure virtual

◆ getdGap()

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

◆ getDiffNormal()

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

◆ getGap()

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

◆ getNormal() [1/3]

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

◆ getNormal() [2/3]

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

◆ getNormal() [3/3]

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

◆ getPointCoords()

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

Definition at line 111 of file RigidBodies.hpp.

112  {
113  auto t_off = getBodyOffset();
114  pointCoords(i) = t_coords(i) - t_off(i) + t_disp(i);
115  return pointCoords;
116  };
Tensor1< double, 3 > pointCoords
Definition: RigidBodies.hpp:40
Tensor1< double, 3 > getBodyOffset()

◆ getRollerDataForTag()

virtual MoFEMErrorCode RigidBodyData::getRollerDataForTag ( )
pure virtual

◆ saveBasicDataOnTag()

MoFEMErrorCode RigidBodyData::saveBasicDataOnTag ( moab::Interface &  moab_debug,
EntityHandle vertex 
)

Definition at line 118 of file RigidBodies.hpp.

119  {
121 
122  Tag tag0;
123  int def_int = 0;
124  std::array<double, 9> def_real;
125  std::fill(def_real.begin(), def_real.end(), 0);
126 
128 
129  CHKERR moab_debug.tag_get_handle("_ROLLER_ID", 1, MB_TYPE_INTEGER, tag0,
130  MB_TAG_CREAT | MB_TAG_SPARSE, &def_int);
131  CHKERR moab_debug.tag_set_data(tag0, &vertex, 1, &this->iD);
132 
133  CHKERR moab_debug.tag_get_handle("_BODY_TYPE", 1, MB_TYPE_INTEGER, tag0,
134  MB_TAG_CREAT | MB_TAG_SPARSE, &def_int);
135  CHKERR moab_debug.tag_set_data(tag0, &vertex, 1, &this->bodyType);
136 
137  CHKERR moab_debug.tag_get_handle("_ORIENTATION", 3, MB_TYPE_DOUBLE, tag0,
138  MB_TAG_CREAT | MB_TAG_SPARSE, &def_real);
139  CHKERR moab_debug.tag_set_data(tag0, &vertex, 1, &this->oRientation(0));
140 
141  CHKERR moab_debug.tag_get_handle("_ROLLER_DATA", 9, MB_TYPE_DOUBLE, tag0,
142  MB_TAG_CREAT | MB_TAG_SPARSE, &def_real);
143  CHKERR moab_debug.tag_set_data(tag0, &vertex, 1, this->dataForTags.data());
144 
145 
146  // tags
147  // int roller_id
148  // int body_type
149  // vec3 orientation
150  // tensor9 roller_data
151  // 0 radius1 = 5
152  // 1 inner_radius = 1
153  // 2 height = 1
154  // 3 angle = 45
155  // 4 angle_a = 45
156  // 5 angle_b = 45
157  //
158 
160  }
#define CHKERR
Inline error check.
Definition: definitions.h:548
virtual MoFEMErrorCode getRollerDataForTag()=0

Member Data Documentation

◆ BodyDispScaled

VectorDouble RigidBodyData::BodyDispScaled

Definition at line 36 of file RigidBodies.hpp.

◆ bodyType

int RigidBodyData::bodyType

Definition at line 49 of file RigidBodies.hpp.

◆ closestPoint

Tensor1<double, 3> RigidBodyData::closestPoint

Definition at line 41 of file RigidBodies.hpp.

◆ dataForTags

array<double, 9> RigidBodyData::dataForTags

Definition at line 50 of file RigidBodies.hpp.

◆ defaultOrientation

Tensor1<double, 3> RigidBodyData::defaultOrientation

Definition at line 42 of file RigidBodies.hpp.

◆ dGap

Tensor1<double, 3> RigidBodyData::dGap

Definition at line 39 of file RigidBodies.hpp.

◆ diffNormal

Tensor2<double, 3, 3> RigidBodyData::diffNormal

Definition at line 45 of file RigidBodies.hpp.

◆ gAp

double RigidBodyData::gAp

Definition at line 37 of file RigidBodies.hpp.

◆ i

Index<'i', 3> RigidBodyData::i

Definition at line 31 of file RigidBodies.hpp.

◆ iD

const int RigidBodyData::iD

Definition at line 30 of file RigidBodies.hpp.

◆ j

Index<'j', 3> RigidBodyData::j

Definition at line 32 of file RigidBodies.hpp.

◆ k

Index<'k', 3> RigidBodyData::k

Definition at line 33 of file RigidBodies.hpp.

◆ oRientation

Tensor1<double, 3> RigidBodyData::oRientation

Definition at line 43 of file RigidBodies.hpp.

◆ originCoords

VectorDouble RigidBodyData::originCoords

Definition at line 34 of file RigidBodies.hpp.

◆ pointCoords

Tensor1<double, 3> RigidBodyData::pointCoords

Definition at line 40 of file RigidBodies.hpp.

◆ positionDataParamName

string RigidBodyData::positionDataParamName

Definition at line 47 of file RigidBodies.hpp.

◆ rollerDisp

VectorDouble RigidBodyData::rollerDisp

Definition at line 35 of file RigidBodies.hpp.

◆ rotationMat

Tensor2<double, 3, 3> RigidBodyData::rotationMat

Definition at line 44 of file RigidBodies.hpp.

◆ tNormal

Tensor1<double, 3> RigidBodyData::tNormal

Definition at line 38 of file RigidBodies.hpp.


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