v0.14.0
Public Member Functions | Public Attributes | List of all members
BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs Struct Reference

\biref operator to calculate the LHS for the RVE boundary conditions More...

#include <users_modules/homogenisation/src/BCs_RVELagrange_Trac_Rigid_Rot.hpp>

Inheritance diagram for BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs:
[legend]
Collaboration diagram for BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs:
[legend]

Public Member Functions

 OpRVEBCsLhs (MoFEM::Interface &_mField, const string field_name, const string lagrang_field_name, Mat aij, RVEBC_Data &data)
 
PetscErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 

Public Attributes

MoFEM::InterfacemField
 
Mat Aij
 
RVEBC_DatadAta
 
MatrixDouble Mat_face
 
MatrixDouble Mat_face_Tran
 

Detailed Description

\biref operator to calculate the LHS for the RVE boundary conditions

Definition at line 27 of file BCs_RVELagrange_Trac_Rigid_Rot.hpp.

Constructor & Destructor Documentation

◆ OpRVEBCsLhs()

BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::OpRVEBCsLhs ( MoFEM::Interface _mField,
const string  field_name,
const string  lagrang_field_name,
Mat  aij,
RVEBC_Data data 
)
inline

Definition at line 33 of file BCs_RVELagrange_Trac_Rigid_Rot.hpp.

35  :
37  lagrang_field_name, field_name, UserDataOperator::OPROWCOL),
38  mField(_mField),Aij(aij), dAta(data){
39  //This will make sure to loop over all entities
40  //(e.g. for order=2 it will make doWork to loop 16 time)
41  sYmm = false;
42  }

Member Function Documentation

◆ doWork()

PetscErrorCode BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
DataForcesAndSourcesCore::EntData row_data,
DataForcesAndSourcesCore::EntData col_data 
)
inline

Definition at line 47 of file BCs_RVELagrange_Trac_Rigid_Rot.hpp.

52  {
53  PetscFunctionBegin;
54 
55  try {
56  if(row_data.getIndices().size()==0) PetscFunctionReturn(0);
57  if(col_data.getIndices().size()==0) PetscFunctionReturn(0);
58  if(col_type==MBVERTEX){
59 
60 
61  Mat_face.resize(3,9);
62  Mat_face.clear();
63  Mat_face_Tran.resize(9,3);
64  Mat_face_Tran.clear();
65 
66  int num_nodes;
67  const EntityHandle* conn_face;
68  double coords_face[9];
69  EntityHandle face_tri = getNumeredEntFiniteElementPtr()->getEnt(); //handle of finite element
70 
71 
72  rval = mField.get_moab().get_connectivity(face_tri,conn_face,num_nodes,true); CHKERRQ_MOAB(rval);
73  rval = mField.get_moab().get_coords(conn_face,num_nodes,coords_face); CHKERRQ_MOAB(rval);
74 
75  for(int nn=0; nn<3; nn++){
76  Mat_face(0,3*nn+1)=-coords_face[3*nn+2]; Mat_face(0,3*nn+2)= coords_face[3*nn+1];
77  Mat_face(1,3*nn+0)= coords_face[3*nn+2]; Mat_face(1,3*nn+2)=-coords_face[3*nn+0];
78  Mat_face(2,3*nn+0)=-coords_face[3*nn+1]; Mat_face(2,3*nn+1)= coords_face[3*nn+0];
79  }
80 
81  // Matrix C1
82  int nb_rows=row_data.getIndices().size();
83  int nb_cols=col_data.getIndices().size();
85  Aij,nb_rows,&row_data.getIndices()[0],nb_cols,&col_data.getIndices()[0],&Mat_face(0,0),ADD_VALUES
86  ); CHKERRQ(ierr);
87 
88  // Matrix C1T
89  noalias(Mat_face_Tran) = trans(Mat_face);
91  Aij,nb_cols,&col_data.getIndices()[0],nb_rows,&row_data.getIndices()[0],&Mat_face_Tran(0,0),ADD_VALUES
92  ); CHKERRQ(ierr);
93 
94  }
95 
96  } catch (const std::exception& ex) {
97  ostringstream ss;
98  ss << "throw in method: " << ex.what() << endl;
99  SETERRQ(PETSC_COMM_SELF,1,ss.str().c_str());
100  }
101  PetscFunctionReturn(0);
102  }

Member Data Documentation

◆ Aij

Mat BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::Aij

Definition at line 30 of file BCs_RVELagrange_Trac_Rigid_Rot.hpp.

◆ dAta

RVEBC_Data& BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::dAta

Definition at line 31 of file BCs_RVELagrange_Trac_Rigid_Rot.hpp.

◆ Mat_face

MatrixDouble BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::Mat_face

Definition at line 44 of file BCs_RVELagrange_Trac_Rigid_Rot.hpp.

◆ Mat_face_Tran

MatrixDouble BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::Mat_face_Tran

Definition at line 45 of file BCs_RVELagrange_Trac_Rigid_Rot.hpp.

◆ mField

MoFEM::Interface& BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::mField

Definition at line 29 of file BCs_RVELagrange_Trac_Rigid_Rot.hpp.


The documentation for this struct was generated from the following file:
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
EntityHandle
BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::Mat_face
MatrixDouble Mat_face
Definition: BCs_RVELagrange_Trac_Rigid_Rot.hpp:44
CHKERRQ_MOAB
#define CHKERRQ_MOAB(a)
check error code of MoAB function
Definition: definitions.h:467
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::dAta
RVEBC_Data & dAta
Definition: BCs_RVELagrange_Trac_Rigid_Rot.hpp:31
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::Aij
Mat Aij
Definition: BCs_RVELagrange_Trac_Rigid_Rot.hpp:30
BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::Mat_face_Tran
MatrixDouble Mat_face_Tran
Definition: BCs_RVELagrange_Trac_Rigid_Rot.hpp:45
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::mField
MoFEM::Interface & mField
Definition: BCs_RVELagrange_Trac_Rigid_Rot.hpp:29
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76