v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
STLRigidBody Struct Reference

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

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

Public Member Functions

 STLRigidBody (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 >
Tensor1< double, 3 > getdGapImpl (Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_normal)
 
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)
 
MoFEMErrorCode getBodyOptions ()
 
- Public Member Functions inherited from RigidBodyData
 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

moab::Core mbInstance
 
moab::Interface * bodyMesh
 
shared_ptr< OrientedBoxTreeTool > orientedBox
 
EntityHandle treeMeshset
 
Tag tagNormal
 
EntityHandle triangle
 
- 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
 
VectorDouble BodyDirectionScaled
 
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
 
boost::shared_ptr< TimeAccelerogrammethodOpForRollerPosition
 
boost::shared_ptr< TimeAccelerogrammethodOpForRollerDirection
 
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 1243 of file RigidBodies.hpp.

Constructor & Destructor Documentation

◆ STLRigidBody()

STLRigidBody::STLRigidBody ( VectorDouble  c_coords,
VectorDouble  roller_disp,
int  id 
)
inline

Definition at line 1251 of file RigidBodies.hpp.

1252 : RigidBodyData(c_coords, roller_disp, id) {
1254 }
RigidBodyData()=delete
moab::Core mbInstance
moab::Interface * bodyMesh

Member Function Documentation

◆ getBodyOptions()

MoFEMErrorCode STLRigidBody::getBodyOptions ( )
inlinevirtual

Implements RigidBodyData.

Definition at line 1334 of file RigidBodies.hpp.

1334 {
1336 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "", "");
1337 PetscBool rflg;
1338 string param = "-body_file_name" + to_string(iD);
1339 char mesh_file_name[255];
1340 CHKERR PetscOptionsString(param.c_str(), "file name", "", "mesh.vtk",
1341 mesh_file_name, 255, &rflg);
1342 if (string(mesh_file_name).substr(string(mesh_file_name).length() - 3) !=
1343 "vtk")
1344 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1345 "body surface should be vtk surface");
1346 bodyMesh->load_file(mesh_file_name, 0, "");
1347 Range tris;
1348 CHKERR bodyMesh->get_entities_by_dimension(0, 2, tris, false);
1349 orientedBox =
1350 make_shared<OrientedBoxTreeTool>(bodyMesh, "ROOTSETSURF", true);
1351 CHKERR orientedBox->build(tris, treeMeshset);
1352 char normal_tag_name[255] = "Normals";
1353 if (bodyMesh->tag_get_handle(normal_tag_name, tagNormal) != MB_SUCCESS) {
1354 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1355 "Normal tag name %s cannot be found. Please "
1356 "provide the correct name. or run ./calculate_normals");
1357 }
1358
1359 ierr = PetscOptionsEnd();
1360 CHKERRG(ierr);
1362 }
static PetscErrorCode ierr
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
char mesh_file_name[255]
const int iD
Definition: RigidBodies.hpp:30
shared_ptr< OrientedBoxTreeTool > orientedBox
EntityHandle treeMeshset

◆ getdGap()

Tensor1< double, 3 > STLRigidBody::getdGap ( Tensor1< TPack3, 3 > &  t_coords,
Tensor1< TPack1, 3 > &  t_normal 
)
inlinevirtual

Implements RigidBodyData.

Definition at line 1285 of file RigidBodies.hpp.

1286 {
1287 return getdGapImpl(t_coords, t_normal);
1288 }
Tensor1< double, 3 > getdGapImpl(Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_normal)

◆ getdGapImpl()

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

Definition at line 1307 of file RigidBodies.hpp.

1308 {
1309 return dGap;
1310 };
Tensor1< double, 3 > dGap
Definition: RigidBodies.hpp:40

◆ getDiffNormal()

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

Implements RigidBodyData.

Definition at line 1275 of file RigidBodies.hpp.

1277 {
1278 return getDiffNormalImpl(t_coords, t_disp, t_normal);
1279 }
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 > STLRigidBody::getDiffNormalImpl ( Tensor1< T1, 3 > &  t_coords,
Tensor1< T2, 3 > &  t_disp,
Tensor1< T3, 3 > &  t_normal 
)
inline

Definition at line 1313 of file RigidBodies.hpp.

1315 {
1316 Tensor2<double, 3, 3> diff_normal;
1317 Tensor1<double, 3> norm_diff;
1318 const double clean_gap = gAp;
1319 const double eps = 1e-8;
1320 for (int ii = 0; ii != 3; ++ii) {
1321 Tensor1<double, 3> pert_disp{t_disp(0), t_disp(1), t_disp(2)};
1322 pert_disp(ii) += eps;
1323 auto pert_normal = getNormal(t_coords, pert_disp);
1324 norm_diff(i) = (pert_normal(i) - t_normal(i)) / eps;
1325 dGap(ii) = (gAp - clean_gap) / eps;
1326
1327 for (int jj = 0; jj != 3; ++jj) {
1328 diff_normal(jj, ii) = norm_diff(jj);
1329 }
1330 }
1331 return diff_normal;
1332 };
static const double eps
Index< 'i', 3 > i
Definition: RigidBodies.hpp:31
Tensor1< double, 3 > getNormal(Tensor1< TPack3, 3 > &t_coords, Tensor1< TPack1, 3 > &t_disp)

◆ getGap()

double STLRigidBody::getGap ( Tensor1< TPack3, 3 > &  t_coords)
inlinevirtual

Implements RigidBodyData.

Definition at line 1281 of file RigidBodies.hpp.

1281 {
1282 return getGapImpl(t_coords);
1283 }
double getGapImpl(Tensor1< T1, 3 > &t_coords)

◆ getGapImpl()

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

Definition at line 1302 of file RigidBodies.hpp.

1302 {
1303 return gAp;
1304 };

◆ getNormal() [1/3]

Tensor1< double, 3 > STLRigidBody::getNormal ( Tensor1< double, 3 > &  t_coords,
Tensor1< double, 3 > &  t_disp 
)
inlinevirtual

Implements RigidBodyData.

Definition at line 1270 of file RigidBodies.hpp.

1271 {
1272 return getNormalImpl(t_coords, t_disp);
1273 }
Tensor1< double, 3 > & getNormalImpl(Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_disp)

◆ getNormal() [2/3]

Tensor1< double, 3 > STLRigidBody::getNormal ( Tensor1< TPack3, 3 > &  t_coords,
Tensor1< double, 3 > &  t_disp 
)
inlinevirtual

Implements RigidBodyData.

Definition at line 1266 of file RigidBodies.hpp.

1267 {
1268 return getNormalImpl(t_coords, t_disp);
1269 }

◆ getNormal() [3/3]

Tensor1< double, 3 > STLRigidBody::getNormal ( Tensor1< TPack3, 3 > &  t_coords,
Tensor1< TPack1, 3 > &  t_disp 
)
inlinevirtual

Implements RigidBodyData.

Definition at line 1262 of file RigidBodies.hpp.

1263 {
1264 return getNormalImpl(t_coords, t_disp);
1265 }

◆ getNormalImpl()

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

Definition at line 1291 of file RigidBodies.hpp.

1292 {
1293 auto t_d = getPointCoords(t_coords, t_disp);
1294
1295 CHKERR orientedBox->closest_to_location(&t_d(0), treeMeshset,
1296 &closestPoint(0), triangle);
1297 CHKERR bodyMesh->tag_get_data(tagNormal, &triangle, 1, &tNormal(0));
1298 gAp = tNormal(i) * (pointCoords(i) - closestPoint(i));
1299 return tNormal;
1300 };
Tensor1< double, 3 > pointCoords
Definition: RigidBodies.hpp:41
Tensor1< double, 3 > tNormal
Definition: RigidBodies.hpp:39
Tensor1< double, 3 > & getPointCoords(Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_disp)
Tensor1< double, 3 > closestPoint
Definition: RigidBodies.hpp:42
EntityHandle triangle

◆ getRollerDataForTag()

MoFEMErrorCode STLRigidBody::getRollerDataForTag ( )
inlinevirtual

Implements RigidBodyData.

Definition at line 1256 of file RigidBodies.hpp.

1256 {
1258 bodyType = STL;
1260 }
@ STL
Definition: RigidBodies.hpp:24
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440

Member Data Documentation

◆ bodyMesh

moab::Interface* STLRigidBody::bodyMesh

Definition at line 1245 of file RigidBodies.hpp.

◆ mbInstance

moab::Core STLRigidBody::mbInstance

Definition at line 1244 of file RigidBodies.hpp.

◆ orientedBox

shared_ptr<OrientedBoxTreeTool> STLRigidBody::orientedBox

Definition at line 1246 of file RigidBodies.hpp.

◆ tagNormal

Tag STLRigidBody::tagNormal

Definition at line 1248 of file RigidBodies.hpp.

◆ treeMeshset

EntityHandle STLRigidBody::treeMeshset

Definition at line 1247 of file RigidBodies.hpp.

◆ triangle

EntityHandle STLRigidBody::triangle

Definition at line 1249 of file RigidBodies.hpp.


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