v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
BCs_RVELagrange_Trac::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.hpp"

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

Public Member Functions

 OpRVEBCsLhs (const string field_name, const string lagrang_field_name, Mat aij, RVEBC_Data &data, CommonFunctions &common_functions, bool ho_geometry=false)
 
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
 
CommonFunctionscommonFunctions
 
bool hoGeometry
 
Mat Aij
 
MatrixDouble K
 
MatrixDouble transK
 
MatrixDouble NTN
 
ublas::vector< MatrixDouble > N_mat_row
 
MatrixDouble N_mat_col
 
MatrixDouble H_mat
 

Detailed Description

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

Definition at line 164 of file BCs_RVELagrange_Trac.hpp.

Constructor & Destructor Documentation

◆ OpRVEBCsLhs()

BCs_RVELagrange_Trac::OpRVEBCsLhs::OpRVEBCsLhs ( const string  field_name,
const string  lagrang_field_name,
Mat  aij,
RVEBC_Data data,
CommonFunctions common_functions,
bool  ho_geometry = false 
)
inline

Definition at line 171 of file BCs_RVELagrange_Trac.hpp.

173 :
174 FaceElementForcesAndSourcesCore::UserDataOperator(lagrang_field_name, field_name, UserDataOperator::OPROWCOL),
175 dAta(data),commonFunctions(common_functions),hoGeometry(ho_geometry),Aij(aij) {
176 //This will make sure to loop over all entities
177 // (e.g. for order=2 it will make doWork to loop 16 time)
178 sYmm = false;
179 }
constexpr auto field_name
@ OPROWCOL
operator doWork is executed on FE rows &columns

Member Function Documentation

◆ doWork()

PetscErrorCode BCs_RVELagrange_Trac::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 185 of file BCs_RVELagrange_Trac.hpp.

190 {
191 PetscFunctionBegin;
192
193 try {
194 if(row_data.getIndices().size()==0) PetscFunctionReturn(0);
195 if(col_data.getIndices().size()==0) PetscFunctionReturn(0);
196
197 int nb_row = row_data.getIndices().size();
198 int nb_col = col_data.getIndices().size();
199
200
201 const FENumeredDofEntity *dof_ptr;
202 ierr = getNumeredEntFiniteElementPtr()->getColDofsByPetscGlobalDofIdx(col_data.getIndices()[0],&dof_ptr); CHKERRQ(ierr);
203 int rank = dof_ptr->getNbOfCoeffs();
204
205 K.resize(nb_row,nb_col);
206 K.clear();
207
208 N_mat_row.resize(col_data.getN().size1());
209 for(unsigned int gg = 0;gg<col_data.getN().size1();gg++) {
210 double area;
211 if(hoGeometry) {
212 area = norm_2(getNormalsAtGaussPt(gg))*0.5;
213 } else {
214 area = getArea();
215 }
216 double val = getGaussPts()(2,gg)*area;
217
218 if(col_type==MBVERTEX) {
219 ierr = commonFunctions.shapeMat(rank,gg,col_data,N_mat_row[gg]); CHKERRQ(ierr);
220 }
221 ierr = commonFunctions.shapeMat(rank,gg,col_data,N_mat_col); CHKERRQ(ierr);
222 ierr = commonFunctions.hMat(getNormalsAtGaussPt(gg),rank,N_mat_row(gg),H_mat); CHKERRQ(ierr);
223
224 if(gg==0){
225 NTN = prod(trans(N_mat_row(gg)),N_mat_col); //we don't need to define size NTN
226 K = val*prod(H_mat,NTN); //we don't need to defiend size of K
227 } else {
228 noalias(NTN) = prod(trans(N_mat_row(gg)),N_mat_col);
229 K += val*prod(H_mat, NTN);
230 }
231 }
232
233 // Matrix C
234 int nb_rows=row_data.getIndices().size();
235 int nb_cols=col_data.getIndices().size();
237 Aij,nb_rows,&row_data.getIndices()[0],nb_cols,&col_data.getIndices()[0],&K(0,0),ADD_VALUES
238 ); CHKERRQ(ierr);
239
240 // Matrix CT
241 transK = trans(K);
243 Aij,nb_cols,&col_data.getIndices()[0],nb_rows,&row_data.getIndices()[0],&transK(0,0),ADD_VALUES
244 ); CHKERRQ(ierr);
245
246 } catch (const std::exception& ex) {
247 ostringstream ss;
248 ss << "throw in method: " << ex.what() << endl;
249 SETERRQ(PETSC_COMM_SELF,MOFEM_NOT_IMPLEMENTED,ss.str().c_str());
250 }
251 PetscFunctionReturn(0);
252 }
static PetscErrorCode ierr
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
PetscErrorCode shapeMat(int rank, unsigned int gg, DataForcesAndSurcesCore::EntData &col_data, MatrixDouble &n)
PetscErrorCode hMat(VectorDouble face_normal, int rank, MatrixDouble &n_mat, MatrixDouble &h_mat)
ublas::vector< MatrixDouble > N_mat_row

Member Data Documentation

◆ Aij

Mat BCs_RVELagrange_Trac::OpRVEBCsLhs::Aij

Definition at line 169 of file BCs_RVELagrange_Trac.hpp.

◆ commonFunctions

CommonFunctions& BCs_RVELagrange_Trac::OpRVEBCsLhs::commonFunctions

Definition at line 167 of file BCs_RVELagrange_Trac.hpp.

◆ dAta

RVEBC_Data& BCs_RVELagrange_Trac::OpRVEBCsLhs::dAta

Definition at line 166 of file BCs_RVELagrange_Trac.hpp.

◆ H_mat

MatrixDouble BCs_RVELagrange_Trac::OpRVEBCsLhs::H_mat

Definition at line 183 of file BCs_RVELagrange_Trac.hpp.

◆ hoGeometry

bool BCs_RVELagrange_Trac::OpRVEBCsLhs::hoGeometry

Definition at line 168 of file BCs_RVELagrange_Trac.hpp.

◆ K

MatrixDouble BCs_RVELagrange_Trac::OpRVEBCsLhs::K

Definition at line 181 of file BCs_RVELagrange_Trac.hpp.

◆ N_mat_col

MatrixDouble BCs_RVELagrange_Trac::OpRVEBCsLhs::N_mat_col

Definition at line 183 of file BCs_RVELagrange_Trac.hpp.

◆ N_mat_row

ublas::vector<MatrixDouble > BCs_RVELagrange_Trac::OpRVEBCsLhs::N_mat_row

Definition at line 182 of file BCs_RVELagrange_Trac.hpp.

◆ NTN

MatrixDouble BCs_RVELagrange_Trac::OpRVEBCsLhs::NTN

Definition at line 181 of file BCs_RVELagrange_Trac.hpp.

◆ transK

MatrixDouble BCs_RVELagrange_Trac::OpRVEBCsLhs::transK

Definition at line 181 of file BCs_RVELagrange_Trac.hpp.


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