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

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

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

Public Member Functions

 CylinderRigidBody (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 , 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 >
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)
 
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

double rAdius
 
double hEight
 
- 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 384 of file RigidBodies.hpp.

Constructor & Destructor Documentation

◆ CylinderRigidBody()

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

Definition at line 387 of file RigidBodies.hpp.

388 : RigidBodyData(c_coords, roller_disp, id), rAdius(1), hEight(INFINITY) {}
RigidBodyData()=delete

Member Function Documentation

◆ getBodyOptions()

MoFEMErrorCode CylinderRigidBody::getBodyOptions ( )
inlinevirtual

Implements RigidBodyData.

Definition at line 572 of file RigidBodies.hpp.

572 {
574 PetscBool rflg;
575 PetscInt nb_dirs = 3;
576 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "", "");
577 string param = "-direction" + to_string(iD);
578 CHKERR PetscOptionsGetRealArray(NULL, NULL, param.c_str(), &oRientation(0),
579 &nb_dirs, &rflg);
580 param = "-radius" + to_string(iD);
581 CHKERR PetscOptionsScalar(param.c_str(), "set roller radius", "", rAdius,
582 &rAdius, PETSC_NULL);
583 param = "-height" + to_string(iD);
584 CHKERR PetscOptionsScalar(param.c_str(), "set cylinder height", "", hEight,
585 &hEight, PETSC_NULL);
586 if (rflg)
588 ierr = PetscOptionsEnd();
589 CHKERRG(ierr);
591 }
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
#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
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
Tensor1< double, 3 > oRientation
Definition: RigidBodies.hpp:44
MoFEMErrorCode computeRotationMatrix()
Definition: RigidBodies.hpp:91
const int iD
Definition: RigidBodies.hpp:30

◆ getdGap()

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

Implements RigidBodyData.

Definition at line 421 of file RigidBodies.hpp.

422 {
423 return getdGapImpl(t_coords, t_normal);
424 }
Tensor1< double, 3 > getdGapImpl(Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_normal)

◆ getdGapImpl()

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

Definition at line 563 of file RigidBodies.hpp.

564 {
565 // dGap(i) =
566 // pointCoords(j) * diffNormal(i, j) +
567 // tNormal(j) * (kronecker_delta(i, j) - 2. * rAdius * diffNormal(i,
568 // j));
569 return dGap;
570 };
Tensor1< double, 3 > dGap
Definition: RigidBodies.hpp:40

◆ getDiffNormal()

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

Implements RigidBodyData.

Definition at line 411 of file RigidBodies.hpp.

413 {
414 return getDiffNormalImpl(t_coords, t_disp, t_normal);
415 }
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 > CylinderRigidBody::getDiffNormalImpl ( Tensor1< T1, 3 > &  t_coords,
Tensor1< T2, 3 > &  t_disp,
Tensor1< T3, 3 > &  t_normal 
)
inline

Definition at line 535 of file RigidBodies.hpp.

537 {
538
539 Tensor2<double, 3, 3> diff_normal;
540 Tensor1<double, 3> norm_diff;
541 const double clean_gap = gAp;
542 const double eps = 1e-8;
543 for (int ii = 0; ii != 3; ++ii) {
544 Tensor1<double, 3> pert_disp{t_disp(0), t_disp(1), t_disp(2)};
545 pert_disp(ii) += eps;
546 auto pert_normal = getNormal(t_coords, pert_disp);
547 norm_diff(i) = (pert_normal(i) - t_normal(i)) / eps;
548 dGap(ii) = (gAp - clean_gap) / eps;
549
550 for (int jj = 0; jj != 3; ++jj) {
551 diff_normal(jj, ii) = norm_diff(jj);
552 }
553 }
554 return diff_normal;
555 };
static const double eps
Tensor1< double, 3 > getNormal(Tensor1< TPack3, 3 > &t_coords, Tensor1< TPack1, 3 > &t_disp)
Index< 'i', 3 > i
Definition: RigidBodies.hpp:31

◆ getGap()

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

Implements RigidBodyData.

Definition at line 417 of file RigidBodies.hpp.

417 {
418 return getGapImpl(t_coords);
419 }
double getGapImpl(Tensor1< T1, 3 > &t_coords)

◆ getGapImpl()

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

Definition at line 557 of file RigidBodies.hpp.

557 {
558 // gAp = tNormal(i) * (pointCoords(i) - tNormal(i) * rAdius);
559 return gAp;
560 };

◆ getNormal() [1/3]

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

Implements RigidBodyData.

Definition at line 406 of file RigidBodies.hpp.

407 {
408 return getNormalImpl(t_coords, t_disp);
409 }
Tensor1< double, 3 > & getNormalImpl(Tensor1< T1, 3 > &t_coords, Tensor1< T2, 3 > &t_disp)

◆ getNormal() [2/3]

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

Implements RigidBodyData.

Definition at line 402 of file RigidBodies.hpp.

403 {
404 return getNormalImpl(t_coords, t_disp);
405 }

◆ getNormal() [3/3]

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

Implements RigidBodyData.

Definition at line 398 of file RigidBodies.hpp.

399 {
400 return getNormalImpl(t_coords, t_disp);
401 }

◆ getNormalImpl()

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

Definition at line 427 of file RigidBodies.hpp.

428 {
429 auto t_d = getPointCoords(t_coords, t_disp);
431 pointCoords(i) = rotationMat(i, j) * t_d(j);
432 Tensor1<double, 3> c_center(0., 0., 0.);
433 Tensor1<double, 3> q_minus_c(0., 0., 0.);
434 Tensor1<double, 3> K(0., 0., 0.);
435 auto Q = pointCoords;
436 auto P_Q = pointCoords;
437 auto Q_K = pointCoords;
438
439 // https://www.geometrictools.com/Documentation/DistanceToCircle3.pdf
440 const double half_height = hEight / 2.;
441 const double rad2 = rAdius * rAdius;
442 tNormal(i) = 0.;
443
444 if (pointCoords(2) > half_height) {
445 tNormal(2) = 1;
446 if (pow(pointCoords(1), 2) + pow(pointCoords(0), 2) < rad2) {
447 gAp = tNormal(i) * pointCoords(i) - half_height;
448 } else {
449 c_center(2) = half_height;
450 Q(2) = c_center(2); // Q
451 q_minus_c(i) = Q(i) - c_center(i);
452 q_minus_c.normalize();
453 K(i) = c_center(i) + rAdius * q_minus_c(i);
454 P_Q(i) = pointCoords(i) - Q(i);
455 Q_K(i) = Q(i) - K(i);
456 gAp = sqrt(P_Q.l2() * P_Q.l2() + Q_K.l2() * Q_K.l2());
457 }
458 } else if (pointCoords(2) < -half_height) {
459 tNormal(2) = -1;
460 if (pow(pointCoords(1), 2) + pow(pointCoords(0), 2) < rad2) {
461 gAp = tNormal(i) * pointCoords(i) - half_height;
462 } else {
463 c_center(2) = -half_height;
464 Q(2) = c_center(2); // Q
465 q_minus_c(i) = Q(i) - c_center(i);
466 q_minus_c.normalize();
467 K(i) = c_center(i) + rAdius * q_minus_c(i);
468 P_Q(i) = pointCoords(i) - Q(i);
469 Q_K(i) = Q(i) - K(i);
470 gAp = sqrt(P_Q.l2() * P_Q.l2() + Q_K.l2() * Q_K.l2());
471 }
472 } else {
473 pointCoords(2) = 0;
476 gAp = tNormal(i) * (pointCoords(i) - tNormal(i) * rAdius);
477 }
478
479 // auto length = [&](auto v) { return sqrt(v(0) * v(0) + v(1) * v(1)); };
480
481 // auto sd_cylinder = [&](Tensor1<double, 3> p) {
482 // typedef Tensor1<double, 2> vec2;
483 // Index<'i', 2> i;
484 // vec2 d(length(p) - rAdius, abs(p(2)) - half_height);
485 // vec2 mx_d(max(d(0), 0.0), max(d(1), 0.0));
486 // return min(max(d(0), d(1)), 0.0) + length(mx_d);
487 // };
488
489 // const double eps = 1e-6;
490 // auto cyl_distance = sd_cylinder(pointCoords);
491 // gAp = cyl_distance;
492
493 // auto get_d = [&](int a) {
494 // auto pt = pointCoords;
495 // pt(a) += eps;
496 // return (sd_cylinder(pt) - cyl_distance) / eps;
497 // };
498
499 // for (int dd = 0; dd != 3; ++dd)
500 // tNormal(dd) = get_d(dd);
501
502 // for debugging
503 // auto my_rand_mn = [](double M, double N) {
504 // return M + (rand() / (RAND_MAX / (N - M)));
505 // };
506 // auto t_off = getBodyOffset();
507 // std::ofstream afile("cylinder_pts.csv", std::ios::ate);
508 // if (afile.is_open()) {
509 // afile << "X,Y,Z\n";
510 // Tensor1<double, 3> p(0, 0, 0);
511 // for (int n = 0; n != 1e6; n++) {
512 // double x = my_rand_mn(-300, 300);
513 // double y = my_rand_mn(-300, 300);
514 // double z = my_rand_mn(-300, 300);
515 // Tensor1<double, 3> coords(x, y, z);
516 // if (sd_cylinder(coords) < 0) {
517
518 // afile << x + t_off(0) << "," << y + t_off(1) << "," << z + t_off(2)
519 // << "\n";
520 // }
521 // }
522 // afile.close();
523 // }
524 // cout << "kill process " << endl;
525
527
528 auto normal_copy = tNormal;
529 tNormal(i) = rotationMat(j, i) * normal_copy(j);
530
531 return tNormal;
532 };
static Index< 'K', 3 > K
Tensor1< T, Tensor_Dim > normalize()
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)
Tensor2< double, 3, 3 > rotationMat
Definition: RigidBodies.hpp:45
Index< 'j', 3 > j
Definition: RigidBodies.hpp:32

◆ getRollerDataForTag()

MoFEMErrorCode CylinderRigidBody::getRollerDataForTag ( )
inlinevirtual

Implements RigidBodyData.

Definition at line 390 of file RigidBodies.hpp.

390 {
393 dataForTags[0] = rAdius;
394 dataForTags[2] = hEight;
396 }
@ CYLINDER
Definition: RigidBodies.hpp:20
#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
array< double, 9 > dataForTags
Definition: RigidBodies.hpp:55

Member Data Documentation

◆ hEight

double CylinderRigidBody::hEight

Definition at line 386 of file RigidBodies.hpp.

◆ rAdius

double CylinderRigidBody::rAdius

Definition at line 385 of file RigidBodies.hpp.


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