19#ifndef __BCS_RVELAGRANGE_TRAC_HPP
20#define __BCS_RVELAGRANGE_TRAC_HPP
34 int rank,
unsigned int gg,DataForcesAndSurcesCore::EntData &col_data,MatrixDouble &
n
38 int shape_size = col_data.getN().size2();
40 n.resize(rank,shape_size*rank);
44 for(
int ii=0; ii<shape_size; ii++){
45 double val = col_data.getN()(gg,ii);
46 for(
int jj=0; jj<rank; jj++){
52 PetscFunctionReturn(0);
59 VectorDouble face_normal,
71 if(face_normal[0]>0) {
75 }
else if(face_normal[0]<0) {
79 }
else if(face_normal[1]>0) {
83 }
else if(face_normal[1]<0) {
87 }
else if(face_normal[2]>0) {
91 }
else if(face_normal[2]<0) {
96 SETERRQ(PETSC_COMM_SELF,MOFEM_IMPOSIBLE_CASE,
"Can not be ?!");
99 int num_col = n_mat.size2();
100 h_mat.resize(6,num_col);
103 for(
int bb = 0; bb<num_col/3; bb++) {
105 for(
int rr = 0; rr<6; rr++) {
106 for(
int cc = 0; cc<3; cc++) {
119 if(face_normal[0]>0) {
122 if(face_normal[0]<0) {
125 if(face_normal[1]>0) {
128 if(face_normal[1]<0) {
131 if(face_normal[2]>0) {
134 if(face_normal[2]<0) {
138 int num_col = n_mat.size2();
139 h_mat.resize(3,num_col);
142 for(
int bb = 0; bb<num_col; bb++) {
144 for(
int rr = 0; rr<3; rr++) {
145 for(
int cc = 0; cc<1; cc++) {
157 PetscFunctionReturn(0);
164 struct OpRVEBCsLhs:
public FaceElementForcesAndSourcesCore::UserDataOperator {
186 int row_side,
int col_side,
187 EntityType row_type,EntityType col_type,
188 DataForcesAndSurcesCore::EntData &row_data,
189 DataForcesAndSurcesCore::EntData &col_data
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 const FENumeredDofEntity *dof_ptr;
202 ierr = getNumeredEntFiniteElementPtr()->getColDofsByPetscGlobalDofIdx(col_data.getIndices()[0],&dof_ptr); CHKERRQ(
ierr);
203 int rank = dof_ptr->getNbOfCoeffs();
205 K.resize(nb_row,nb_col);
208 N_mat_row.resize(col_data.getN().size1());
209 for(
unsigned int gg = 0;gg<col_data.getN().size1();gg++) {
212 area = norm_2(getNormalsAtGaussPt(gg))*0.5;
216 double val = getGaussPts()(2,gg)*area;
218 if(col_type==MBVERTEX) {
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
243 Aij,nb_cols,&col_data.getIndices()[0],nb_rows,&row_data.getIndices()[0],&
transK(0,0),ADD_VALUES
246 }
catch (
const std::exception& ex) {
248 ss <<
"throw in method: " << ex.what() << endl;
251 PetscFunctionReturn(0);
276 PetscErrorCode
doWork(
int side,EntityType type,DataForcesAndSurcesCore::EntData &data) {
279 if(type!=MBVERTEX) PetscFunctionReturn(0);
280 if(data.getIndices().size()==0) PetscFunctionReturn(0);
283 const FENumeredDofEntity *dof_ptr;
284 ierr = getNumeredEntFiniteElementPtr()->getColDofsByPetscGlobalDofIdx(data.getIndices()[0],&dof_ptr); CHKERRQ(
ierr);
285 int rank = dof_ptr->getNbOfCoeffs();
288 for(
unsigned int gg = 0;gg<data.getN().size1();gg++) {
291 area = norm_2(getNormalsAtGaussPt(gg))*0.5;
295 double val = getGaussPts()(2,gg)*area;
297 x = getHoCoordsAtGaussPts()(gg,0);
298 y = getHoCoordsAtGaussPts()(gg,1);
299 z = getHoCoordsAtGaussPts()(gg,2);
303 X_mat.resize(rank,6,
false);
311 X_mat.resize(rank,3,
false);
332 PetscFunctionReturn(0);
351 PetscErrorCode
doWork(
int side,EntityType type,DataForcesAndSurcesCore::EntData &data) {
354 if(data.getIndices().size()==0) PetscFunctionReturn(0);
359 for(
int ii = 0;ii<
F.size();ii++) {
362 f.resize(D_mat.size1(),
false);
364 ierr = VecSetValues(
F[ii],data.getIndices().size(),&data.getIndices()[0],&
f[0],ADD_VALUES); CHKERRQ(
ierr);
367 PetscFunctionReturn(0);
372 string field_name,
string lagrang_field_name,
string mesh_nodals_positions,Mat aij,vector<Vec> &f
375 bool ho_geometry =
false;
379 map<int,RVEBC_Data>::iterator sit =
setOfRVEBC.begin();
393 PetscFunctionReturn(0);
412 PetscErrorCode
doWork(
int side,EntityType type,DataForcesAndSurcesCore::EntData &data) {
414 if(data.getIndices().size()==0) PetscFunctionReturn(0);
415 if(
dAta.
tRis.find(getNumeredEntFiniteElementPtr()->getEnt())==
dAta.
tRis.end()) PetscFunctionReturn(0);
419 const int indices_6[6] = {0, 1, 2, 3, 4, 5};
420 const int indices_3[3] = {0, 1, 2};
431 PetscFunctionReturn(0);
437 string field_name,
string lagrang_field_name,
string mesh_nodals_positions,Vec stress_homo
440 bool ho_geometry =
false;
444 map<int,RVEBC_Data>::iterator sit =
setOfRVEBC.begin();
453 PetscFunctionReturn(0);
static PetscErrorCode ierr
virtual bool check_field(const std::string &name) const =0
check if field is in database
const double n
refractive index of diffusive medium
constexpr auto field_name
map< int, RVEBC_Data > setOfRVEBC
maps side set id with appropriate FluxData
CommonFunctions common_functions
MoFEM::Interface & mField
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)
\biref operator to calculate the LHS for the RVE bounary conditions
CommonFunctions & commonFunctions
PetscErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSurcesCore::EntData &row_data, DataForcesAndSurcesCore::EntData &col_data)
ublas::vector< MatrixDouble > N_mat_row
OpRVEBCsLhs(const string field_name, const string lagrang_field_name, Mat aij, RVEBC_Data &data, CommonFunctions &common_functions, bool ho_geometry=false)
\biref operator to calculate the RHS of the constrain for the RVE boundary conditions
PetscErrorCode doWork(int side, EntityType type, DataForcesAndSurcesCore::EntData &data)
OpRVEBCsRhs_Assemble(const string lagrang_field_name, vector< Vec > &f, RVEBC_Data &data, CommonData &common_data)
VectorDouble applied_strain
\biref operator to calculate the RHS of the constrain for the RVE boundary conditions
OpRVEBCsRhs_Cal(const string field_name, RVEBC_Data &data, CommonData &common_data, CommonFunctions &common_functions, bool ho_geometry=false)
PetscErrorCode doWork(int side, EntityType type, DataForcesAndSurcesCore::EntData &data)
CommonFunctions & commonFunctions
VectorDouble stressVector
PetscErrorCode doWork(int side, EntityType type, DataForcesAndSurcesCore::EntData &data)
OpRVEHomoStress_Assemble(const string lagrang_field_name, Vec stress_homo, RVEBC_Data &data, CommonData &common_data)
BCs_RVELagrange_Trac(MoFEM::Interface &m_field)
PetscErrorCode setRVEBCsOperators(string field_name, string lagrang_field_name, string mesh_nodals_positions, Mat aij, vector< Vec > &f)
PetscErrorCode setRVEBCsHomoStressOperators(string field_name, string lagrang_field_name, string mesh_nodals_positions, Vec stress_homo)
CommonFunctions commonFunctions
Deprecated interface functions.