v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
BCs_RVELagrange_Trac_Rigid_Trans::OpRVEBCsLhs Struct Reference

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

#include "users_modules/mofem_um_homogenisation/src/BCs_RVELagrange_Trac_Rigid_Trans.hpp"

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

Public Member Functions

 OpRVEBCsLhs (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, DataForcesAndSurcesCore::EntData &row_data, DataForcesAndSurcesCore::EntData &col_data)
 

Public Attributes

RVEBC_DatadAta
 
Mat Aij
 
MatrixDouble Mat_face
 
MatrixDouble Mat_face_Tran
 

Detailed Description

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

Definition at line 24 of file BCs_RVELagrange_Trac_Rigid_Trans.hpp.

Constructor & Destructor Documentation

◆ OpRVEBCsLhs()

BCs_RVELagrange_Trac_Rigid_Trans::OpRVEBCsLhs::OpRVEBCsLhs ( const string  field_name,
const string  lagrang_field_name,
Mat  aij,
RVEBC_Data data 
)
inline

Definition at line 28 of file BCs_RVELagrange_Trac_Rigid_Trans.hpp.

28 :
29 FaceElementForcesAndSourcesCore::UserDataOperator(
30 lagrang_field_name,field_name,UserDataOperator::OPROWCOL
31 ),
32 dAta(data),Aij(aij) {
33 //This will make sure to loop over all intities (e.g. for order=2 it will make doWork to loop 16 time)
34 sYmm = false;
35 }
constexpr auto field_name
@ OPROWCOL
operator doWork is executed on FE rows &columns

Member Function Documentation

◆ doWork()

PetscErrorCode BCs_RVELagrange_Trac_Rigid_Trans::OpRVEBCsLhs::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
DataForcesAndSurcesCore::EntData &  row_data,
DataForcesAndSurcesCore::EntData &  col_data 
)
inline

Definition at line 40 of file BCs_RVELagrange_Trac_Rigid_Trans.hpp.

45 {
46 PetscFunctionBegin;
47
48 try {
49 if(row_data.getIndices().size()==0) PetscFunctionReturn(0);
50 if(col_data.getIndices().size()==0) PetscFunctionReturn(0);
51 if(col_type==MBVERTEX) {
52
53
54 const FENumeredDofEntity *dof_ptr;
55 ierr = getNumeredEntFiniteElementPtr()->getColDofsByPetscGlobalDofIdx(col_data.getIndices()[0],&dof_ptr); CHKERRQ(ierr);
56 int rank = dof_ptr->getNbOfCoeffs();
57
58 Mat_face.resize(rank,3*rank);
59 Mat_face.clear();
60 Mat_face_Tran.resize(3*rank,rank);
61 Mat_face_Tran.clear();
62
63 switch(rank) {
64 case 3:
65 for(int nn=0; nn<3; nn++){
66 Mat_face(0,3*nn+0) = 1.0;
67 Mat_face(1,3*nn+1) = 1.0;
68 Mat_face(2,3*nn+2) = 1.0;
69 }
70 break;
71 case 1:
72 Mat_face(0,0) = 1.0;
73 Mat_face(0,1) = 1.0;
74 Mat_face(0,2) = 1.0;
75 break;
76 default:
77 SETERRQ(PETSC_COMM_SELF,1,"not implemented");
78 }
79
80 // Matrix C1
81 int nb_rows=row_data.getIndices().size();
82 int nb_cols=col_data.getIndices().size();
84 Aij,nb_rows,&row_data.getIndices()[0],nb_cols,&col_data.getIndices()[0],&Mat_face(0,0),ADD_VALUES
85 ); CHKERRQ(ierr);
86
87 // Matrix C1T
88 noalias(Mat_face_Tran) = trans(Mat_face);
90 Aij,nb_cols,&col_data.getIndices()[0],nb_rows,&row_data.getIndices()[0],&Mat_face_Tran(0,0),ADD_VALUES
91 ); CHKERRQ(ierr);
92
93 }
94 } catch (const std::exception& ex) {
95 ostringstream ss;
96 ss << "throw in method: " << ex.what() << endl;
97 SETERRQ(PETSC_COMM_SELF,1,ss.str().c_str());
98 }
99 PetscFunctionReturn(0);
100 }
static PetscErrorCode ierr
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.

Member Data Documentation

◆ Aij

Mat BCs_RVELagrange_Trac_Rigid_Trans::OpRVEBCsLhs::Aij

Definition at line 27 of file BCs_RVELagrange_Trac_Rigid_Trans.hpp.

◆ dAta

RVEBC_Data& BCs_RVELagrange_Trac_Rigid_Trans::OpRVEBCsLhs::dAta

Definition at line 26 of file BCs_RVELagrange_Trac_Rigid_Trans.hpp.

◆ Mat_face

MatrixDouble BCs_RVELagrange_Trac_Rigid_Trans::OpRVEBCsLhs::Mat_face

Definition at line 37 of file BCs_RVELagrange_Trac_Rigid_Trans.hpp.

◆ Mat_face_Tran

MatrixDouble BCs_RVELagrange_Trac_Rigid_Trans::OpRVEBCsLhs::Mat_face_Tran

Definition at line 38 of file BCs_RVELagrange_Trac_Rigid_Trans.hpp.


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