v0.14.0
Public Member Functions | Public Attributes | List of all members
CohesiveElement::CohesiveInterfaceElement::OpRhs Struct Reference

Operator calculate right hand side vector. More...

#include <users_modules/basic_finite_elements/cohesive_interface/src/CohesiveInterfaceElement.hpp>

Inheritance diagram for CohesiveElement::CohesiveInterfaceElement::OpRhs:
[legend]
Collaboration diagram for CohesiveElement::CohesiveInterfaceElement::OpRhs:
[legend]

Public Member Functions

 OpRhs (const std::string field_name, CommonData &common_data, PhysicalEquation &physical_eqations)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 

Public Attributes

CommonDatacommonData
 
PhysicalEquationphysicalEqations
 
VectorDouble traction
 
VectorDouble Nf
 

Detailed Description

Operator calculate right hand side vector.

Definition at line 431 of file CohesiveInterfaceElement.hpp.

Constructor & Destructor Documentation

◆ OpRhs()

CohesiveElement::CohesiveInterfaceElement::OpRhs::OpRhs ( const std::string  field_name,
CommonData common_data,
PhysicalEquation physical_eqations 
)
inline

Definition at line 435 of file CohesiveInterfaceElement.hpp.

435  :
436  FlatPrismElementForcesAndSourcesCore::UserDataOperator(field_name,ForcesAndSourcesCore::UserDataOperator::OPROW),
437  commonData(common_data),physicalEqations(physical_eqations) {}

Member Function Documentation

◆ doWork()

MoFEMErrorCode CohesiveElement::CohesiveInterfaceElement::OpRhs::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)
inline

Definition at line 440 of file CohesiveInterfaceElement.hpp.

440  {
442 
443  try {
444  int nb_dofs = data.getIndices().size();
445  if(nb_dofs == 0) MoFEMFunctionReturnHot(0);
446  if(physicalEqations.pRisms.find(getNumeredEntFiniteElementPtr()->getEnt()) == physicalEqations.pRisms.end()) {
448  }
449  Nf.resize(nb_dofs);
450  Nf.clear();
451  int nb_gauss_pts = data.getN().size1();
452  for(int gg = 0;gg<nb_gauss_pts;gg++) {
454  double w = getGaussPts()(2,gg)*cblas_dnrm2(3,&getNormalsAtGaussPtsF3()(gg,0),1)*0.5;
455  for(int nn = 0;nn<nb_dofs/3;nn++) {
456  for(int dd = 0;dd<3;dd++) {
457  Nf[3*nn+dd] += w*data.getN(gg)[nn]*traction[dd];
458  }
459  }
460  }
461  ierr = VecSetValues(getFEMethod()->snes_f,
462  data.getIndices().size(),&data.getIndices()[0],&Nf[0],ADD_VALUES); CHKERRG(ierr);
463  } catch (const std::exception& ex) {
464  std::ostringstream ss;
465  ss << "throw in method: " << ex.what() << std::endl;
466  SETERRQ(PETSC_COMM_SELF,1,ss.str().c_str());
467  }
469  }

Member Data Documentation

◆ commonData

CommonData& CohesiveElement::CohesiveInterfaceElement::OpRhs::commonData

Definition at line 433 of file CohesiveInterfaceElement.hpp.

◆ Nf

VectorDouble CohesiveElement::CohesiveInterfaceElement::OpRhs::Nf

Definition at line 439 of file CohesiveInterfaceElement.hpp.

◆ physicalEqations

PhysicalEquation& CohesiveElement::CohesiveInterfaceElement::OpRhs::physicalEqations

Definition at line 434 of file CohesiveInterfaceElement.hpp.

◆ traction

VectorDouble CohesiveElement::CohesiveInterfaceElement::OpRhs::traction

Definition at line 439 of file CohesiveInterfaceElement.hpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
PlasticOps::w
double w(double eqiv, double dot_tau, double f, double sigma_y, double sigma_Y)
Definition: PlasticOpsGeneric.hpp:276
CohesiveElement::CohesiveInterfaceElement::OpRhs::physicalEqations
PhysicalEquation & physicalEqations
Definition: CohesiveInterfaceElement.hpp:434
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1576
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::dd
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
CohesiveElement::CohesiveInterfaceElement::OpRhs::commonData
CommonData & commonData
Definition: CohesiveInterfaceElement.hpp:433
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
CohesiveElement::CohesiveInterfaceElement::PhysicalEquation::calculateTraction
virtual MoFEMErrorCode calculateTraction(VectorDouble &traction, int gg, CommonData &common_data, const FEMethod *fe_method)
Calculate tractions.
Definition: CohesiveInterfaceElement.hpp:203
CohesiveElement::CohesiveInterfaceElement::PhysicalEquation::pRisms
Range pRisms
Definition: CohesiveInterfaceElement.hpp:62
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
CohesiveElement::CohesiveInterfaceElement::OpRhs::Nf
VectorDouble Nf
Definition: CohesiveInterfaceElement.hpp:439
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
CohesiveElement::CohesiveInterfaceElement::OpRhs::traction
VectorDouble traction
Definition: CohesiveInterfaceElement.hpp:439
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483