\biref operator to calculate the LHS for the RVE bounary conditions
More...
#include <users_modules/homogenisation/src/BCs_RVELagrange_Trac.hpp>
|
| 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, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data) |
|
\biref operator to calculate the LHS for the RVE bounary conditions
Definition at line 164 of file BCs_RVELagrange_Trac.hpp.
◆ 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 |
◆ doWork()
Definition at line 185 of file BCs_RVELagrange_Trac.hpp.
194 if(row_data.getIndices().size()==0) PetscFunctionReturn(0);
195 if(col_data.getIndices().size()==0) PetscFunctionReturn(0);
197 int nb_row = row_data.getIndices().size();
198 int nb_col = col_data.getIndices().size();
201 getNumeredEntFiniteElementPtr()->getRowDofsByPetscGlobalDofIdx(
202 col_data.getIndices()[0]);
203 const FENumeredDofEntity *dof_ptr;
204 if (
auto ptr = weak_ptr_dof.lock())
209 int rank = dof_ptr->getNbOfCoeffs();
211 K.resize(nb_row,nb_col);
214 N_mat_row.resize(col_data.getN().size1());
215 for(
unsigned int gg = 0;gg<col_data.getN().size1();gg++) {
218 area = norm_2(getNormalsAtGaussPts(gg))*0.5;
222 double val = getGaussPts()(2,gg)*area;
224 if(col_type==MBVERTEX) {
240 int nb_rows=row_data.getIndices().size();
241 int nb_cols=col_data.getIndices().size();
243 Aij,nb_rows,&row_data.getIndices()[0],nb_cols,&col_data.getIndices()[0],&
K(0,0),ADD_VALUES
249 Aij,nb_cols,&col_data.getIndices()[0],nb_rows,&row_data.getIndices()[0],&
transK(0,0),ADD_VALUES
252 }
catch (
const std::exception& ex) {
254 ss <<
"throw in method: " << ex.what() << endl;
257 PetscFunctionReturn(0);
◆ Aij
Mat BCs_RVELagrange_Trac::OpRVEBCsLhs::Aij |
◆ commonFunctions
◆ dAta
RVEBC_Data& BCs_RVELagrange_Trac::OpRVEBCsLhs::dAta |
◆ H_mat
MatrixDouble BCs_RVELagrange_Trac::OpRVEBCsLhs::H_mat |
◆ hoGeometry
bool BCs_RVELagrange_Trac::OpRVEBCsLhs::hoGeometry |
MatrixDouble BCs_RVELagrange_Trac::OpRVEBCsLhs::K |
◆ N_mat_col
MatrixDouble BCs_RVELagrange_Trac::OpRVEBCsLhs::N_mat_col |
◆ N_mat_row
ublas::vector<MatrixDouble > BCs_RVELagrange_Trac::OpRVEBCsLhs::N_mat_row |
◆ NTN
MatrixDouble BCs_RVELagrange_Trac::OpRVEBCsLhs::NTN |
◆ transK
MatrixDouble BCs_RVELagrange_Trac::OpRVEBCsLhs::transK |
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.