\biref operator to calculate the LHS for the RVE boundary conditions
More...
#include <users_modules/homogenisation/src/BCs_RVELagrange_Trac_Rigid_Rot.hpp>
\biref operator to calculate the LHS for the RVE boundary conditions
Definition at line 27 of file BCs_RVELagrange_Trac_Rigid_Rot.hpp.
◆ 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 |
◆ doWork()
Definition at line 47 of file BCs_RVELagrange_Trac_Rigid_Rot.hpp.
56 if(row_data.getIndices().size()==0) PetscFunctionReturn(0);
57 if(col_data.getIndices().size()==0) PetscFunctionReturn(0);
58 if(col_type==MBVERTEX){
68 double coords_face[9];
69 EntityHandle face_tri = getNumeredEntFiniteElementPtr()->getEnt();
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];
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
91 Aij,nb_cols,&col_data.getIndices()[0],nb_rows,&row_data.getIndices()[0],&
Mat_face_Tran(0,0),ADD_VALUES
96 }
catch (
const std::exception& ex) {
98 ss <<
"throw in method: " << ex.what() << endl;
99 SETERRQ(PETSC_COMM_SELF,1,ss.str().c_str());
101 PetscFunctionReturn(0);
◆ Aij
Mat BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::Aij |
◆ dAta
RVEBC_Data& BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::dAta |
◆ Mat_face
MatrixDouble BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::Mat_face |
◆ Mat_face_Tran
MatrixDouble BCs_RVELagrange_Trac_Rigid_Rot::OpRVEBCsLhs::Mat_face_Tran |
◆ mField
The documentation for this struct was generated from the following file:
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.