v0.14.0
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | Protected Member Functions | List of all members
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator Struct Reference

default operator for Contact Prism element More...

#include <src/finite_elements/ContactPrismElementForcesAndSourcesCore.hpp>

Inheritance diagram for MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator:
[legend]
Collaboration diagram for MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator:
[legend]

Public Types

enum  FaceType {
  FACEMASTER = 1 << 0 , FACESLAVE = 1 << 1 , FACEMASTERMASTER = 1 << 2 , FACEMASTERSLAVE = 1 << 3 ,
  FACESLAVEMASTER = 1 << 4 , FACESLAVESLAVE = 1 << 5 , FACELAST = 1 << 6
}
 
- 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 = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > >
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 

Public Member Functions

 UserDataOperator (const FieldSpace space)
 
 UserDataOperator (const std::string &field_name, const char type)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const char face_type)
 
 UserDataOperator (const std::string &field_name, const char type, const char face_type)
 
int getFaceType () const
 Get operator types. More...
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 
double getAreaMaster ()
 get face aRea Master More...
 
double getAreaSlave ()
 get face aRea Slave More...
 
VectorAdaptor getNormalMaster ()
 get face normal vector to Master face More...
 
VectorAdaptor getTangentMasterOne ()
 get first face tangent vector to Master face More...
 
VectorAdaptor getTangentMasterTwo ()
 get second face tangent vector to Master face More...
 
VectorAdaptor getNormalSlave ()
 get face normal vector to Slave face More...
 
VectorAdaptor getTangentSlaveOne ()
 get first face tangent vector to Slave face More...
 
VectorAdaptor getTangentSlaveTwo ()
 get second face tangent vector to Slave face More...
 
MatrixDoublegetGaussPtsMaster ()
 get Gauss point at Master face More...
 
MatrixDoublegetGaussPtsSlave ()
 get Gauss point at Slave face More...
 
auto getFTensor0IntegrationWeightSlave ()
 Get integration weights for slave side. More...
 
auto getFTensor0IntegrationWeightMaster ()
 Get integration weights for master side. More...
 
VectorDouble getCoordsMaster ()
 get triangle coordinates More...
 
VectorDouble getCoordsSlave ()
 get triangle coordinates More...
 
MatrixDoublegetCoordsAtGaussPtsMaster ()
 get coordinates at Gauss pts on full prism. More...
 
MatrixDoublegetCoordsAtGaussPtsSlave ()
 get coordinates at Gauss pts on full prism. More...
 
const ContactPrismElementForcesAndSourcesCoregetContactPrismElementForcesAndSourcesCore ()
 return pointer to triangle finite element object More...
 
MoFEMErrorCode loopSideVolumes (const string fe_name, VolumeElementForcesAndSourcesCoreOnContactPrismSide &fe_method, const int side_type, const EntityHandle ent_for_side)
 
- 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. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
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 More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, 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. More...
 
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. More...
 
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. More...
 
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. More...
 
- 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. More...
 
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. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Public Attributes

char faceType
 
- 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. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 

Protected Member Functions

ForcesAndSourcesCoregetSidePtrFE () const
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 

Additional Inherited Members

- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

default operator for Contact Prism element

Examples
continuity_check_on_contact_prism_side_ele.cpp.

Definition at line 208 of file ContactPrismElementForcesAndSourcesCore.hpp.

Member Enumeration Documentation

◆ FaceType

enum MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::FaceType
   \brief Controls loop over faces and face combination on element

    FACEMASTER is used when column or row data needs to be accessed

located at master face FACESLAVE is used when column or row data needs to be accessed located at slave face FACEMASTERMASTER is used for accessing simultaneously row and col data located at master face. FACEMASTERSLAVE is used for accessing simultaneously row data that is located on master face and col data located at slave face. FACESLAVEMASTER is used for accessing simultaneously row data that is located on slave face and col data located at master face. FACESLAVESLAVE is used for accessing simultaneously row and col data located at slave face.

Enumerator
FACEMASTER 
FACESLAVE 
FACEMASTERMASTER 
FACEMASTERSLAVE 
FACESLAVEMASTER 
FACESLAVESLAVE 
FACELAST 

Definition at line 251 of file ContactPrismElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ UserDataOperator() [1/5]

MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const FieldSpace  space)
inline

Definition at line 211 of file ContactPrismElementForcesAndSourcesCore.hpp.

212 : ForcesAndSourcesCore::UserDataOperator(space) {}

◆ UserDataOperator() [2/5]

MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const std::string &  field_name,
const char  type 
)
inline

Definition at line 214 of file ContactPrismElementForcesAndSourcesCore.hpp.

215 : ForcesAndSourcesCore::UserDataOperator(field_name, type) {}
constexpr auto field_name

◆ UserDataOperator() [3/5]

MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const std::string &  row_field_name,
const std::string &  col_field_name,
const char  type 
)
inline

Definition at line 217 of file ContactPrismElementForcesAndSourcesCore.hpp.

219 : ForcesAndSourcesCore::UserDataOperator(row_field_name, col_field_name,
220 type) {}

◆ UserDataOperator() [4/5]

MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const std::string &  row_field_name,
const std::string &  col_field_name,
const char  type,
const char  face_type 
)
inline

Definition at line 222 of file ContactPrismElementForcesAndSourcesCore.hpp.

225 : ForcesAndSourcesCore::UserDataOperator(row_field_name, col_field_name,
226 type),
227 faceType(face_type) {}

◆ UserDataOperator() [5/5]

MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const std::string &  field_name,
const char  type,
const char  face_type 
)
inline

Definition at line 229 of file ContactPrismElementForcesAndSourcesCore.hpp.

231 : ForcesAndSourcesCore::UserDataOperator(field_name, type),
232 faceType(face_type) {}

Member Function Documentation

◆ getAreaMaster()

double ContactPrismElementForcesAndSourcesCore::UserDataOperator::getAreaMaster ( )
inline

get face aRea Master

Definition at line 411 of file ContactPrismElementForcesAndSourcesCore.hpp.

411 {
412 return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[0];
413}
std::array< double, 2 > aRea
Array storing master and slave faces areas.

◆ getAreaSlave()

double ContactPrismElementForcesAndSourcesCore::UserDataOperator::getAreaSlave ( )
inline

get face aRea Slave

Definition at line 416 of file ContactPrismElementForcesAndSourcesCore.hpp.

416 {
417 return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[1];
418}

◆ getContactPrismElementForcesAndSourcesCore()

const ContactPrismElementForcesAndSourcesCore * ContactPrismElementForcesAndSourcesCore::UserDataOperator::getContactPrismElementForcesAndSourcesCore ( )
inline

return pointer to triangle finite element object

Definition at line 513 of file ContactPrismElementForcesAndSourcesCore.hpp.

514 {
515 return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE);
516}

◆ getCoordsAtGaussPtsMaster()

MatrixDouble & ContactPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPtsMaster ( )
inline

get coordinates at Gauss pts on full prism.

Matrix has size (nb integration points on master)x(3), i.e. coordinates on face Master

Definition at line 500 of file ContactPrismElementForcesAndSourcesCore.hpp.

◆ getCoordsAtGaussPtsSlave()

MatrixDouble & ContactPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPtsSlave ( )
inline

get coordinates at Gauss pts on full prism.

Matrix has size (nb integration points on slave)x(3), i.e. coordinates on face Slave

Definition at line 506 of file ContactPrismElementForcesAndSourcesCore.hpp.

◆ getCoordsMaster()

VectorDouble ContactPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsMaster ( )
inline

get triangle coordinates

Vector has 9 elements, i.e. coordinates on Master face

Definition at line 487 of file ContactPrismElementForcesAndSourcesCore.hpp.

487 {
488 double *data = &(
490 return VectorAdaptor(9, ublas::shallow_array_adaptor<double>(9, data));
491}
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115

◆ getCoordsSlave()

VectorDouble ContactPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsSlave ( )
inline

get triangle coordinates

Vector has 9 elements, i.e. coordinates on Slave face

Definition at line 494 of file ContactPrismElementForcesAndSourcesCore.hpp.

494 {
495 double *data = &(
497 return VectorAdaptor(9, ublas::shallow_array_adaptor<double>(9, data));
498}

◆ getFaceType()

int ContactPrismElementForcesAndSourcesCore::UserDataOperator::getFaceType ( ) const
inline

Get operator types.

Returns
Return operator type

Definition at line 406 of file ContactPrismElementForcesAndSourcesCore.hpp.

406 {
407 return faceType;
408}

◆ getFTensor0IntegrationWeightMaster()

auto ContactPrismElementForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeightMaster ( )
inline

Get integration weights for master side.

for(int gg = 0; gg!=getGaussPts.size2(); ++gg) {
// integrate something
++t_w;
}
auto getFTensor0IntegrationWeight()
Get integration weights.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Returns
FTensor::Tensor0<FTensor::PackPtr<double *, 1>>

Definition at line 480 of file ContactPrismElementForcesAndSourcesCore.hpp.

◆ getFTensor0IntegrationWeightSlave()

auto ContactPrismElementForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeightSlave ( )
inline

Get integration weights for slave side.

for(int gg = 0; gg!=getGaussPts.size2(); ++gg) {
// integrate something
++t_w;
}
Returns
FTensor::Tensor0<FTensor::PackPtr<double *, 1>>

Definition at line 474 of file ContactPrismElementForcesAndSourcesCore.hpp.

◆ getGaussPtsMaster()

MatrixDouble & ContactPrismElementForcesAndSourcesCore::UserDataOperator::getGaussPtsMaster ( )
inline

get Gauss point at Master face

Definition at line 463 of file ContactPrismElementForcesAndSourcesCore.hpp.

◆ getGaussPtsSlave()

MatrixDouble & ContactPrismElementForcesAndSourcesCore::UserDataOperator::getGaussPtsSlave ( )
inline

get Gauss point at Slave face

Definition at line 469 of file ContactPrismElementForcesAndSourcesCore.hpp.

◆ getNormalMaster()

VectorAdaptor ContactPrismElementForcesAndSourcesCore::UserDataOperator::getNormalMaster ( )
inline

get face normal vector to Master face

Definition at line 421 of file ContactPrismElementForcesAndSourcesCore.hpp.

421 {
422 double *data = &(
424 return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
425}
VectorDouble normal
vector storing vector normal to master or slave element

◆ getNormalSlave()

VectorAdaptor ContactPrismElementForcesAndSourcesCore::UserDataOperator::getNormalSlave ( )
inline

get face normal vector to Slave face

Definition at line 442 of file ContactPrismElementForcesAndSourcesCore.hpp.

442 {
443 double *data = &(
445 return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
446}

◆ getNumeredEntFiniteElementPtr()

boost::shared_ptr< const NumeredEntFiniteElement > ContactPrismElementForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr ( ) const
inline

Definition at line 399 of file ContactPrismElementForcesAndSourcesCore.hpp.

400 {
403};
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr

◆ getSidePtrFE()

ForcesAndSourcesCore * MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getSidePtrFE ( ) const
inlineprotected

◆ getTangentMasterOne()

VectorAdaptor ContactPrismElementForcesAndSourcesCore::UserDataOperator::getTangentMasterOne ( )
inline

get first face tangent vector to Master face

Definition at line 427 of file ContactPrismElementForcesAndSourcesCore.hpp.

428 {
429 double *data = &(static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
430 ->tangentMasterOne[0]);
431 return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
432}

◆ getTangentMasterTwo()

VectorAdaptor ContactPrismElementForcesAndSourcesCore::UserDataOperator::getTangentMasterTwo ( )
inline

get second face tangent vector to Master face

Definition at line 434 of file ContactPrismElementForcesAndSourcesCore.hpp.

435 {
436 double *data = &(static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
437 ->tangentMasterTwo[0]);
438 return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
439}

◆ getTangentSlaveOne()

VectorAdaptor ContactPrismElementForcesAndSourcesCore::UserDataOperator::getTangentSlaveOne ( )
inline

get first face tangent vector to Slave face

Definition at line 448 of file ContactPrismElementForcesAndSourcesCore.hpp.

449 {
450 double *data = &(static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
451 ->tangentSlaveOne[0]);
452 return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
453}

◆ getTangentSlaveTwo()

VectorAdaptor ContactPrismElementForcesAndSourcesCore::UserDataOperator::getTangentSlaveTwo ( )
inline

get second face tangent vector to Slave face

Definition at line 455 of file ContactPrismElementForcesAndSourcesCore.hpp.

456 {
457 double *data = &(static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
458 ->tangentSlaveTwo[0]);
459 return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
460}

◆ loopSideVolumes()

MoFEMErrorCode ContactPrismElementForcesAndSourcesCore::UserDataOperator::loopSideVolumes ( const string  fe_name,
VolumeElementForcesAndSourcesCoreOnContactPrismSide fe_method,
const int  side_type,
const EntityHandle  ent_for_side 
)

User call this function to loop over elements on the side of face. This function calls MoFEM::VolumeElementForcesAndSourcesCoreOnContactPrismSide with is operator to do calculations.

Parameters
fe_nameName of the element
methodFinite element object
side_typestates the side from which side element will work (0 for master 1 for slave)
Returns
error code

Definition at line 1189 of file ContactPrismElementForcesAndSourcesCore.cpp.

1192 {
1193 return loopSide(fe_name, &fe_method, side_type, ent_for_side);
1194}
MoFEMErrorCode loopSide(const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, 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 elemen...

Member Data Documentation

◆ faceType

char MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::faceType

Definition at line 261 of file ContactPrismElementForcesAndSourcesCore.hpp.


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