v0.14.0
Public Member Functions | Public Attributes | List of all members
PeriodicNitscheConstrains::PeriodicVolume Struct Reference

#include <users_modules/homogenisation/src/NitschePeriodicMethod.hpp>

Inheritance diagram for PeriodicNitscheConstrains::PeriodicVolume:
[legend]
Collaboration diagram for PeriodicNitscheConstrains::PeriodicVolume:
[legend]

Public Member Functions

 PeriodicVolume (MoFEM::Interface &m_field, CommonData &common_data)
 
int getRule (int order)
 
PetscErrorCode setGaussPts (int order)
 

Public Attributes

EntityHandle forFace
 
CommonDatacommonData
 

Detailed Description

Definition at line 235 of file NitschePeriodicMethod.hpp.

Constructor & Destructor Documentation

◆ PeriodicVolume()

PeriodicNitscheConstrains::PeriodicVolume::PeriodicVolume ( MoFEM::Interface m_field,
CommonData common_data 
)
inline

Definition at line 240 of file NitschePeriodicMethod.hpp.

240  :
241  VolumeElementForcesAndSourcesCore(m_field),
242  forFace(0),
243  commonData(common_data) {}

Member Function Documentation

◆ getRule()

int PeriodicNitscheConstrains::PeriodicVolume::getRule ( int  order)
inline

Definition at line 244 of file NitschePeriodicMethod.hpp.

244 { return -1; }

◆ setGaussPts()

PetscErrorCode PeriodicNitscheConstrains::PeriodicVolume::setGaussPts ( int  order)
inline

Definition at line 245 of file NitschePeriodicMethod.hpp.

245  {
246 
247  PetscFunctionBegin;
248  try {
249  gaussPts.resize(4,0,false);
250  EntityHandle tet = numeredEntFiniteElementPtr->getEnt();
251  int vffgg = 0;
252  for(int ff = 0;ff<4;ff++) {
253  EntityHandle face;
254  rval = mField.get_moab().side_element(tet,2,ff,face); CHKERRQ_MOAB(rval);
255  if(face!=forFace) continue;
256  if(commonData.localCoordsMap.find(face)!=commonData.localCoordsMap.end()) {
257  const double coords_tet[12] = { 0,0,0, 1,0,0, 0,1,0, 0,0,1 };
258  int nb_face_gauss_pts = commonData.localCoordsMap[face].size();
259  gaussPts.resize(4,nb_face_gauss_pts);
260  commonData.inTetTetGaussPtsNumber[tet].resize(nb_face_gauss_pts);
261  for(int ffgg = 0;ffgg!=nb_face_gauss_pts;ffgg++,vffgg++) {
263  double N[3];
264  VectorDouble &local_coords = commonData.localCoordsMap[face][ffgg];
265  ierr = ShapeMBTRI(N,&local_coords[0],&local_coords[1],1); CHKERRQ(ierr);
266  for(int dd = 0;dd<3;dd++) {
267  gaussPts(dd,vffgg) =
268  N[0]*coords_tet[3*dataH1.facesNodes(ff,0)+dd]+
269  N[1]*coords_tet[3*dataH1.facesNodes(ff,1)+dd]+
270  N[2]*coords_tet[3*dataH1.facesNodes(ff,2)+dd];
271  }
272  gaussPts(3,vffgg) = 0;
273  }
274  } else {
275  SETERRQ(PETSC_COMM_SELF,MOFEM_DATA_INCONSISTENCY,"data inconstancy");
276  }
277  }
278  //cerr << gaussPts << endl;
279  } catch (const std::exception& ex) {
280  ostringstream ss;
281  ss << "throw in method: " << ex.what() << endl;
282  SETERRQ(PETSC_COMM_SELF,1,ss.str().c_str());
283  }
284  PetscFunctionReturn(0);
285  }

Member Data Documentation

◆ commonData

CommonData& PeriodicNitscheConstrains::PeriodicVolume::commonData

Definition at line 238 of file NitschePeriodicMethod.hpp.

◆ forFace

EntityHandle PeriodicNitscheConstrains::PeriodicVolume::forFace

Definition at line 237 of file NitschePeriodicMethod.hpp.


The documentation for this struct was generated from the following file:
EntityHandle
CHKERRQ_MOAB
#define CHKERRQ_MOAB(a)
check error code of MoAB function
Definition: definitions.h:454
PeriodicNitscheConstrains::CommonData::inTetTetGaussPtsNumber
map< EntityHandle, vector< int > > inTetTetGaussPtsNumber
Definition: NitschePeriodicMethod.hpp:37
PeriodicNitscheConstrains::PeriodicVolume::forFace
EntityHandle forFace
Definition: NitschePeriodicMethod.hpp:237
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
PeriodicNitscheConstrains::PeriodicVolume::commonData
CommonData & commonData
Definition: NitschePeriodicMethod.hpp:238
N
const int N
Definition: speed_test.cpp:3
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
PeriodicNitscheConstrains::CommonData::localCoordsMap
map< EntityHandle, vector< VectorDouble > > localCoordsMap
Definition: NitschePeriodicMethod.hpp:35
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
PeriodicNitscheConstrains::CommonData::inTetFaceGaussPtsNumber
map< EntityHandle, vector< int > > inTetFaceGaussPtsNumber
Definition: NitschePeriodicMethod.hpp:36
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
ShapeMBTRI
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182