v0.14.0
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/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, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::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  :
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  }

Member Function Documentation

◆ doWork()

PetscErrorCode BCs_RVELagrange_Trac::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 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  auto weak_ptr_dof =
201  getNumeredEntFiniteElementPtr()->getRowDofsByPetscGlobalDofIdx(
202  col_data.getIndices()[0]);
203  const FENumeredDofEntity *dof_ptr;
204  if (auto ptr = weak_ptr_dof.lock())
205  dof_ptr = ptr.get();
206  else
207  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Dof not found");
208 
209  int rank = dof_ptr->getNbOfCoeffs();
210 
211  K.resize(nb_row,nb_col);
212  K.clear();
213 
214  N_mat_row.resize(col_data.getN().size1());
215  for(unsigned int gg = 0;gg<col_data.getN().size1();gg++) {
216  double area;
217  if(hoGeometry) {
218  area = norm_2(getNormalsAtGaussPts(gg))*0.5;
219  } else {
220  area = getArea();
221  }
222  double val = getGaussPts()(2,gg)*area;
223 
224  if(col_type==MBVERTEX) {
225  ierr = commonFunctions.shapeMat(rank,gg,col_data,N_mat_row[gg]); CHKERRQ(ierr);
226  }
227  ierr = commonFunctions.shapeMat(rank,gg,col_data,N_mat_col); CHKERRQ(ierr);
228  ierr = commonFunctions.hMat(getNormalsAtGaussPts(gg),rank,N_mat_row(gg),H_mat); CHKERRQ(ierr);
229 
230  if(gg==0){
231  NTN = prod(trans(N_mat_row(gg)),N_mat_col); //we don't need to define size NTN
232  K = val*prod(H_mat,NTN); //we don't need to defiend size of K
233  } else {
234  noalias(NTN) = prod(trans(N_mat_row(gg)),N_mat_col);
235  K += val*prod(H_mat, NTN);
236  }
237  }
238 
239  // Matrix C
240  int nb_rows=row_data.getIndices().size();
241  int nb_cols=col_data.getIndices().size();
242  ierr = MatSetValues(
243  Aij,nb_rows,&row_data.getIndices()[0],nb_cols,&col_data.getIndices()[0],&K(0,0),ADD_VALUES
244  ); CHKERRQ(ierr);
245 
246  // Matrix CT
247  transK = trans(K);
248  ierr = MatSetValues(
249  Aij,nb_cols,&col_data.getIndices()[0],nb_rows,&row_data.getIndices()[0],&transK(0,0),ADD_VALUES
250  ); CHKERRQ(ierr);
251 
252  } catch (const std::exception& ex) {
253  ostringstream ss;
254  ss << "throw in method: " << ex.what() << endl;
255  SETERRQ(PETSC_COMM_SELF,MOFEM_NOT_IMPLEMENTED,ss.str().c_str());
256  }
257  PetscFunctionReturn(0);
258  }

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:
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:1631
BCs_RVELagrange_Trac::OpRVEBCsLhs::H_mat
MatrixDouble H_mat
Definition: BCs_RVELagrange_Trac.hpp:183
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
BCs_RVELagrange_Trac::CommonFunctions::hMat
PetscErrorCode hMat(VectorDouble face_normal, int rank, MatrixDouble &n_mat, MatrixDouble &h_mat)
Definition: BCs_RVELagrange_Trac.hpp:58
BCs_RVELagrange_Trac::OpRVEBCsLhs::dAta
RVEBC_Data & dAta
Definition: BCs_RVELagrange_Trac.hpp:166
BCs_RVELagrange_Trac::OpRVEBCsLhs::commonFunctions
CommonFunctions & commonFunctions
Definition: BCs_RVELagrange_Trac.hpp:167
BCs_RVELagrange_Trac::OpRVEBCsLhs::Aij
Mat Aij
Definition: BCs_RVELagrange_Trac.hpp:169
BCs_RVELagrange_Disp::common_functions
CommonFunctions common_functions
Definition: BCs_RVELagrange_Disp.hpp:434
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
BCs_RVELagrange_Trac::OpRVEBCsLhs::transK
MatrixDouble transK
Definition: BCs_RVELagrange_Trac.hpp:181
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
BCs_RVELagrange_Trac::OpRVEBCsLhs::N_mat_row
ublas::vector< MatrixDouble > N_mat_row
Definition: BCs_RVELagrange_Trac.hpp:182
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
BCs_RVELagrange_Trac::OpRVEBCsLhs::N_mat_col
MatrixDouble N_mat_col
Definition: BCs_RVELagrange_Trac.hpp:183
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
BCs_RVELagrange_Trac::OpRVEBCsLhs::NTN
MatrixDouble NTN
Definition: BCs_RVELagrange_Trac.hpp:181
BCs_RVELagrange_Trac::OpRVEBCsLhs::hoGeometry
bool hoGeometry
Definition: BCs_RVELagrange_Trac.hpp:168
BCs_RVELagrange_Trac::OpRVEBCsLhs::K
MatrixDouble K
Definition: BCs_RVELagrange_Trac.hpp:181
BCs_RVELagrange_Trac::CommonFunctions::shapeMat
PetscErrorCode shapeMat(int rank, unsigned int gg, DataForcesAndSourcesCore::EntData &col_data, MatrixDouble &n)
Definition: BCs_RVELagrange_Trac.hpp:33