v0.15.0
Loading...
Searching...
No Matches
SurfaceSlidingConstrains::OpJacobian Struct Reference

#include "tools/src/SurfaceSlidingConstrains.hpp"

Inheritance diagram for SurfaceSlidingConstrains::OpJacobian:
[legend]
Collaboration diagram for SurfaceSlidingConstrains::OpJacobian:
[legend]

Public Member Functions

 OpJacobian (int tag, const std::string field_name, boost::shared_ptr< VectorDouble > &active_variables_ptr, boost::shared_ptr< VectorDouble > &results_ptr, boost::shared_ptr< MatrixDouble > &jacobian_ptr, DriverElementOrientation &orientation, bool evaluate_jacobian, const double &alpha)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 
 OpJacobian (int tag, const std::string field_name, boost::shared_ptr< VectorDouble > &active_variables_ptr, boost::shared_ptr< VectorDouble > &results_ptr, boost::shared_ptr< MatrixDouble > &jacobian_ptr, DriverElementOrientation &orientation, bool evaluate_jacobian, const double &alpha)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 
- Public Member Functions inherited from MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
double getArea ()
 get area of face
 
VectorDoublegetNormal ()
 get triangle normal
 
VectorDoublegetTangent1 ()
 get triangle tangent 1
 
VectorDoublegetTangent2 ()
 get triangle tangent 2
 
auto getFTensor1Normal ()
 get normal as tensor
 
auto getFTensor1Tangent1 ()
 get tangentOne as tensor
 
auto getFTensor1Tangent2 ()
 get tangentTwo as tensor
 
int getNumNodes ()
 get element number of nodes
 
const EntityHandlegetConn ()
 get element connectivity
 
VectorDoublegetCoords ()
 get triangle coordinates
 
auto getFTensor1Coords ()
 get get coords at gauss points
 
MatrixDoublegetNormalsAtGaussPts ()
 if higher order geometry return normals at Gauss pts.
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 if higher order geometry return normals at Gauss pts.
 
MatrixDoublegetTangent1AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts.
 
MatrixDoublegetTangent2AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts.
 
auto getFTensor1NormalsAtGaussPts ()
 get normal at integration points
 
auto getFTensor1Tangent1AtGaussPts ()
 get tangent 1 at integration points
 
auto getFTensor1Tangent2AtGaussPts ()
 get tangent 2 at integration points
 
FaceElementForcesAndSourcesCoregetFaceFE ()
 return pointer to Generic Triangle Finite Element object
 
MoFEMErrorCode loopSideVolumes (const string fe_name, VolumeElementForcesAndSourcesCoreOnSide &fe_method)
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement.
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle.
 
int getFEDim () const
 Get dimension of finite element.
 
EntityType getFEType () const
 Get dimension of finite element.
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer.
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity.
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element.
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices.
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices.
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object.
 
int getOpType () const
 Get operator types.
 
void setOpType (const OpType type)
 Set operator type.
 
void addOpType (const OpType type)
 Add operator type.
 
int getNinTheLoop () const
 get number of finite element in the loop
 
int getLoopSize () const
 get size of elements in the loop
 
std::string getFEName () const
 Get name of the element.
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element
 
auto getFTensor0IntegrationWeight ()
 Get integration weights.
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3)
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry.
 
double getMeasure () const
 get measure of element
 
doublegetMeasure ()
 get measure of element
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, boost::shared_ptr< Range > fe_range=nullptr, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopRange (const string &fe_name, ForcesAndSourcesCore *range_fe, boost::shared_ptr< Range > fe_range, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 Iterate over range of elements.
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side.
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side.
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not.
 
void setSymm ()
 set if operator is executed taking in account symmetry
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem
 

Public Attributes

const int tAg
 
boost::shared_ptr< VectorDouble > activeVariablesPtr
 
boost::shared_ptr< VectorDouble > resultsPtr
 
boost::shared_ptr< MatrixDouble > jacobianPtr
 
DriverElementOrientationoRientation
 
bool evaluateJacobian
 
const doubleaLpha
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure.
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity.
 
booldoVertices
 \deprectaed If false skip vertices
 
booldoEdges
 \deprectaed If false skip edges
 
booldoQuads
 \deprectaed
 
booldoTris
 \deprectaed
 
booldoTets
 \deprectaed
 
booldoPrisms
 \deprectaed
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType
 
using DoWorkRhsHookFunType
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Definition at line 414 of file SurfaceSlidingConstrains.hpp.

Constructor & Destructor Documentation

◆ OpJacobian() [1/2]

SurfaceSlidingConstrains::OpJacobian::OpJacobian ( int tag,
const std::string field_name,
boost::shared_ptr< VectorDouble > & active_variables_ptr,
boost::shared_ptr< VectorDouble > & results_ptr,
boost::shared_ptr< MatrixDouble > & jacobian_ptr,
DriverElementOrientation & orientation,
bool evaluate_jacobian,
const double & alpha )
inline

Definition at line 427 of file SurfaceSlidingConstrains.hpp.

434 field_name, UserDataOperator::OPCOL),
435 tAg(tag), activeVariablesPtr(active_variables_ptr),
436 resultsPtr(results_ptr), jacobianPtr(jacobian_ptr),
437 oRientation(orientation), evaluateJacobian(evaluate_jacobian),
438 aLpha(alpha) {}
constexpr auto field_name
boost::shared_ptr< VectorDouble > resultsPtr
boost::shared_ptr< VectorDouble > activeVariablesPtr
boost::shared_ptr< MatrixDouble > jacobianPtr

◆ OpJacobian() [2/2]

SurfaceSlidingConstrains::OpJacobian::OpJacobian ( int tag,
const std::string field_name,
boost::shared_ptr< VectorDouble > & active_variables_ptr,
boost::shared_ptr< VectorDouble > & results_ptr,
boost::shared_ptr< MatrixDouble > & jacobian_ptr,
DriverElementOrientation & orientation,
bool evaluate_jacobian,
const double & alpha )
inline

Definition at line 426 of file SurfaceSlidingConstrains.hpp.

433 field_name, UserDataOperator::OPCOL),
434 tAg(tag), activeVariablesPtr(active_variables_ptr),
435 resultsPtr(results_ptr), jacobianPtr(jacobian_ptr),
436 oRientation(orientation), evaluateJacobian(evaluate_jacobian),
437 aLpha(alpha) {}

Member Function Documentation

◆ doWork() [1/2]

MoFEMErrorCode SurfaceSlidingConstrains::OpJacobian::doWork ( int side,
EntityType type,
EntitiesFieldData::EntData & data )
inline

Definition at line 440 of file SurfaceSlidingConstrains.hpp.

441 {
443 if (type != MBVERTEX)
445
446 VectorInt &indices = data.getIndices();
447
448 trace_on(tAg);
449
450 ublas::vector<adouble> lambda_dofs(3);
451 for (int dd = 0; dd != 3; ++dd) {
452 lambda_dofs[dd] <<= (*activeVariablesPtr)[dd];
453 }
454 ublas::vector<adouble> position_dofs(9);
455 for (int dd = 0; dd != 9; ++dd) {
456 position_dofs[dd] <<= (*activeVariablesPtr)[3 + dd];
457 }
458
459 FTensor::Index<'i', 3> i;
460 FTensor::Index<'j', 3> j;
461 FTensor::Index<'k', 3> k;
465
467 getFEMethod());
469
470 int nb_gauss_pts = data.getN().size1();
471 int nb_base_functions = data.getN().size2();
472 auto t_base1 = data.getFTensor0N();
473 auto t_base2 = data.getFTensor0N();
474 auto t_diff_base = data.getFTensor1DiffN<2>();
476 FTensor::Tensor1<adouble, 3> t_position_ksi;
477 FTensor::Tensor1<adouble, 3> t_position_eta;
481
482 ublas::vector<adouble> c_vec(3);
483 ublas::vector<adouble> f_vec(9);
484 c_vec.clear();
485 f_vec.clear();
486
487 auto t_coord_ref = getFTensor1CoordsAtGaussPts();
488
489 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
490
491 FTensor::Tensor1<adouble *, 3> t_position_dofs(
492 &position_dofs[0], &position_dofs[1], &position_dofs[2], 3);
493 FTensor::Tensor0<adouble *> t_lambda_dof(&lambda_dofs[0]);
494
495 t_position(i) = 0;
496 t_position_ksi(i) = 0;
497 t_position_eta(i) = 0;
498 lambda = 0;
499
500 for (int bb = 0; bb != nb_base_functions; ++bb) {
501 t_position(i) += t_base1 * t_position_dofs(i);
502 t_position_ksi(i) += t_diff_base(N0) * t_position_dofs(i);
503 t_position_eta(i) += t_diff_base(N1) * t_position_dofs(i);
504 lambda += t_base1 * t_lambda_dof;
505 ++t_base1;
506 ++t_position_dofs;
507 ++t_lambda_dof;
508 ++t_diff_base;
509 }
510
511 t_delta(i) = t_position(i) - t_coord_ref(i);
512 t_normal(k) = FTensor::cross(t_position_ksi(i), t_position_eta(j), k);
513
514 double w = getGaussPts()(2, gg) * 0.5 * aLpha;
515 adouble val;
516 FTensor::Tensor0<adouble *> t_c(&c_vec[0]);
517 FTensor::Tensor1<adouble *, 3> t_f(&f_vec[0], &f_vec[1], &f_vec[2], 3);
518
519 adouble area = sqrt(t_normal(i) * t_normal(i));
520
521 for (int bb = 0; bb != nb_base_functions; ++bb) {
522 if (indices[bb] != -1) {
523 val = w * eo * t_base2;
524 t_c += val * t_normal(i) * t_delta(i);
525 val *= lambda;
526 t_f(i) += val * t_normal(i);
527 }
528 ++t_c;
529 ++t_f;
530 ++t_base2;
531 }
532
533 ++t_coord_ref;
534 }
535
536 for (int rr = 0; rr != 3; ++rr)
537 c_vec[rr] >>= (*resultsPtr)[rr];
538
539 for (int rr = 0; rr != 9; ++rr)
540 f_vec(rr) >>= (*resultsPtr)[3 + rr];
541
542 trace_off();
543
544 if (evaluateJacobian) {
545 double *jac_ptr[3 + 9];
546 for (int rr = 0; rr != 3 + 9; ++rr)
547 jac_ptr[rr] = &(*jacobianPtr)(rr, 0);
548
549 // play recorder for jacobians
550 int r =
551 ::jacobian(tAg, 3 + 9, 3 + 9, &(*activeVariablesPtr)[0], jac_ptr);
552 if (r < 0)
553 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
554 "ADOL-C function evaluation with error");
555 }
556
558 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
FTensor::Index< 'i', SPACE_DIM > i
static double lambda
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition ddTensor0.hpp:33
auto cross(const Tensor1_Expr< A, T, 3, i > &a, const Tensor1_Expr< B, U, 3, j > &b, const Index< k, 3 > &)
Definition cross.hpp:10
int r
Definition sdf.py:8
FaceElementForcesAndSourcesCore * getFaceFE()
return pointer to Generic Triangle Finite Element object
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
virtual MoFEMErrorCode getElementOrientation(MoFEM::Interface &m_field, const FEMethod *fe_method_ptr)

◆ doWork() [2/2]

MoFEMErrorCode SurfaceSlidingConstrains::OpJacobian::doWork ( int side,
EntityType type,
EntitiesFieldData::EntData & data )
inline

Definition at line 439 of file SurfaceSlidingConstrains.hpp.

440 {
442 if (type != MBVERTEX)
444
445 VectorInt &indices = data.getIndices();
446
447 trace_on(tAg);
448
449 ublas::vector<adouble> lambda_dofs(3);
450 for (int dd = 0; dd != 3; ++dd) {
451 lambda_dofs[dd] <<= (*activeVariablesPtr)[dd];
452 }
453 ublas::vector<adouble> position_dofs(9);
454 for (int dd = 0; dd != 9; ++dd) {
455 position_dofs[dd] <<= (*activeVariablesPtr)[3 + dd];
456 }
457
458 FTensor::Index<'i', 3> i;
459 FTensor::Index<'j', 3> j;
460 FTensor::Index<'k', 3> k;
464
466 getFEMethod());
468
469 int nb_gauss_pts = data.getN().size1();
470 int nb_base_functions = data.getN().size2();
471 auto t_base1 = data.getFTensor0N();
472 auto t_base2 = data.getFTensor0N();
473 auto t_diff_base = data.getFTensor1DiffN<2>();
475 FTensor::Tensor1<adouble, 3> t_position_ksi;
476 FTensor::Tensor1<adouble, 3> t_position_eta;
480
481 ublas::vector<adouble> c_vec(3);
482 ublas::vector<adouble> f_vec(9);
483 c_vec.clear();
484 f_vec.clear();
485
486 auto t_coord_ref = getFTensor1CoordsAtGaussPts();
487
488 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
489
490 FTensor::Tensor1<adouble *, 3> t_position_dofs(
491 &position_dofs[0], &position_dofs[1], &position_dofs[2], 3);
492 FTensor::Tensor0<adouble *> t_lambda_dof(&lambda_dofs[0]);
493
494 t_position(i) = 0;
495 t_position_ksi(i) = 0;
496 t_position_eta(i) = 0;
497 lambda = 0;
498
499 for (int bb = 0; bb != nb_base_functions; ++bb) {
500 t_position(i) += t_base1 * t_position_dofs(i);
501 t_position_ksi(i) += t_diff_base(N0) * t_position_dofs(i);
502 t_position_eta(i) += t_diff_base(N1) * t_position_dofs(i);
503 lambda += t_base1 * t_lambda_dof;
504 ++t_base1;
505 ++t_position_dofs;
506 ++t_lambda_dof;
507 ++t_diff_base;
508 }
509
510 t_delta(i) = t_position(i) - t_coord_ref(i);
511 t_normal(k) = FTensor::cross(t_position_ksi(i), t_position_eta(j), k);
512
513 double w = getGaussPts()(2, gg) * 0.5 * aLpha;
514 adouble val;
515 FTensor::Tensor0<adouble *> t_c(&c_vec[0]);
516 FTensor::Tensor1<adouble *, 3> t_f(&f_vec[0], &f_vec[1], &f_vec[2], 3);
517
518 adouble area = sqrt(t_normal(i) * t_normal(i));
519
520 for (int bb = 0; bb != nb_base_functions; ++bb) {
521 if (indices[bb] != -1) {
522 val = w * eo * t_base2;
523 t_c += val * t_normal(i) * t_delta(i);
524 val *= lambda;
525 t_f(i) += val * t_normal(i);
526 }
527 ++t_c;
528 ++t_f;
529 ++t_base2;
530 }
531
532 ++t_coord_ref;
533 }
534
535 for (int rr = 0; rr != 3; ++rr)
536 c_vec[rr] >>= (*resultsPtr)[rr];
537
538 for (int rr = 0; rr != 9; ++rr)
539 f_vec(rr) >>= (*resultsPtr)[3 + rr];
540
541 trace_off();
542
543 if (evaluateJacobian) {
544 double *jac_ptr[3 + 9];
545 for (int rr = 0; rr != 3 + 9; ++rr)
546 jac_ptr[rr] = &(*jacobianPtr)(rr, 0);
547
548 // play recorder for jacobians
549 int r =
550 ::jacobian(tAg, 3 + 9, 3 + 9, &(*activeVariablesPtr)[0], jac_ptr);
551 if (r < 0)
552 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
553 "ADOL-C function evaluation with error");
554 }
555
557 }

Member Data Documentation

◆ activeVariablesPtr

boost::shared_ptr< VectorDouble > SurfaceSlidingConstrains::OpJacobian::activeVariablesPtr

Definition at line 419 of file SurfaceSlidingConstrains.hpp.

◆ aLpha

const double & SurfaceSlidingConstrains::OpJacobian::aLpha

Definition at line 425 of file SurfaceSlidingConstrains.hpp.

◆ evaluateJacobian

bool SurfaceSlidingConstrains::OpJacobian::evaluateJacobian

Definition at line 423 of file SurfaceSlidingConstrains.hpp.

◆ jacobianPtr

boost::shared_ptr< MatrixDouble > SurfaceSlidingConstrains::OpJacobian::jacobianPtr

Definition at line 421 of file SurfaceSlidingConstrains.hpp.

◆ oRientation

DriverElementOrientation & SurfaceSlidingConstrains::OpJacobian::oRientation

Definition at line 422 of file SurfaceSlidingConstrains.hpp.

◆ resultsPtr

boost::shared_ptr< VectorDouble > SurfaceSlidingConstrains::OpJacobian::resultsPtr

Definition at line 420 of file SurfaceSlidingConstrains.hpp.

◆ tAg

const int SurfaceSlidingConstrains::OpJacobian::tAg

Definition at line 418 of file SurfaceSlidingConstrains.hpp.


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