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

Operator calculate element stiffens matrix. More...

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

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

Public Member Functions

 OpLhs (const std::string field_name, CommonData &common_data, PhysicalEquation &physical_eqations)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 

Public Attributes

CommonDatacommonData
 
PhysicalEquationphysicalEqations
 
MatrixDouble K
 
MatrixDouble D
 
MatrixDouble ND
 

Detailed Description

Operator calculate element stiffens matrix.

Definition at line 475 of file CohesiveInterfaceElement.hpp.

Constructor & Destructor Documentation

◆ OpLhs()

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

Definition at line 479 of file CohesiveInterfaceElement.hpp.

479  :
480  FlatPrismElementForcesAndSourcesCore::UserDataOperator(field_name,ForcesAndSourcesCore::UserDataOperator::OPROWCOL),
481  commonData(common_data),physicalEqations(physical_eqations) { sYmm = false; }

Member Function Documentation

◆ doWork()

MoFEMErrorCode CohesiveElement::CohesiveInterfaceElement::OpLhs::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)
inline

Definition at line 484 of file CohesiveInterfaceElement.hpp.

489  {
491 
492  try {
493  int nb_row = row_data.getIndices().size();
494  if(nb_row == 0) MoFEMFunctionReturnHot(0);
495  int nb_col = col_data.getIndices().size();
496  if(nb_col == 0) MoFEMFunctionReturnHot(0);
497  if(physicalEqations.pRisms.find(getNumeredEntFiniteElementPtr()->getEnt())
498  == physicalEqations.pRisms.end()) {
500  }
501  //std::cerr << row_side << " " << row_type << " " << row_data.getN() << std::endl;
502  //std::cerr << col_side << " " << col_type << " " << col_data.getN() << std::endl;
503  ND.resize(nb_row,3);
504  K.resize(nb_row,nb_col);
505  K.clear();
506  int nb_gauss_pts = row_data.getN().size1();
507  for(int gg = 0;gg<nb_gauss_pts;gg++) {
509  double w = getGaussPts()(2,gg)*cblas_dnrm2(3,&getNormalsAtGaussPtsF3()(gg,0),1)*0.5;
510  ND.clear();
511  for(int nn = 0; nn<nb_row/3;nn++) {
512  for(int dd = 0;dd<3;dd++) {
513  for(int DD = 0;DD<3;DD++) {
514  ND(3*nn+dd,DD) += row_data.getN(gg)[nn]*D(dd,DD);
515  }
516  }
517  }
518  for(int nn = 0; nn<nb_row/3; nn++) {
519  for(int dd = 0;dd<3;dd++) {
520  for(int NN = 0; NN<nb_col/3; NN++) {
521  for(int DD = 0; DD<3;DD++) {
522  K(3*nn+dd,3*NN+DD) += w*ND(3*nn+dd,DD)*col_data.getN(gg)[NN];
523  }
524  }
525  }
526  }
527  }
528  ierr = MatSetValues(getFEMethod()->snes_B,
529  nb_row,&row_data.getIndices()[0],
530  nb_col,&col_data.getIndices()[0],
531  &K(0,0),ADD_VALUES
532  ); CHKERRG(ierr);
533  } catch (const std::exception& ex) {
534  std::ostringstream ss;
535  ss << "throw in method: " << ex.what() << std::endl;
536  SETERRQ(PETSC_COMM_SELF,1,ss.str().c_str());
537  }
539  }

Member Data Documentation

◆ commonData

CommonData& CohesiveElement::CohesiveInterfaceElement::OpLhs::commonData

Definition at line 477 of file CohesiveInterfaceElement.hpp.

◆ D

MatrixDouble CohesiveElement::CohesiveInterfaceElement::OpLhs::D

Definition at line 483 of file CohesiveInterfaceElement.hpp.

◆ K

MatrixDouble CohesiveElement::CohesiveInterfaceElement::OpLhs::K

Definition at line 483 of file CohesiveInterfaceElement.hpp.

◆ ND

MatrixDouble CohesiveElement::CohesiveInterfaceElement::OpLhs::ND

Definition at line 483 of file CohesiveInterfaceElement.hpp.

◆ physicalEqations

PhysicalEquation& CohesiveElement::CohesiveInterfaceElement::OpLhs::physicalEqations

Definition at line 478 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:460
MoFEM::MatSetValues
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
Definition: EntitiesFieldData.hpp:1644
CohesiveElement::CohesiveInterfaceElement::OpLhs::physicalEqations
PhysicalEquation & physicalEqations
Definition: CohesiveInterfaceElement.hpp:478
CohesiveElement::CohesiveInterfaceElement::OpLhs::K
MatrixDouble K
Definition: CohesiveInterfaceElement.hpp:483
CohesiveElement::CohesiveInterfaceElement::OpLhs::D
MatrixDouble D
Definition: CohesiveInterfaceElement.hpp:483
CohesiveElement::CohesiveInterfaceElement::OpLhs::commonData
CommonData & commonData
Definition: CohesiveInterfaceElement.hpp:477
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
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
CohesiveElement::CohesiveInterfaceElement::PhysicalEquation::pRisms
Range pRisms
Definition: CohesiveInterfaceElement.hpp:62
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
sdf_wavy_2d.w
int w
Definition: sdf_wavy_2d.py:6
CohesiveElement::CohesiveInterfaceElement::OpLhs::ND
MatrixDouble ND
Definition: CohesiveInterfaceElement.hpp:483
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
CohesiveElement::CohesiveInterfaceElement::PhysicalEquation::calculateTangentStiffeness
virtual MoFEMErrorCode calculateTangentStiffeness(MatrixDouble &tangent_matrix, int gg, CommonData &common_data, const FEMethod *fe_method)
Calculate tangent stiffness.
Definition: CohesiveInterfaceElement.hpp:231
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496