v0.14.0
Public Member Functions | Public Attributes | List of all members
BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU Struct Reference

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

Inheritance diagram for BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU:
[legend]
Collaboration diagram for BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU:
[legend]

Public Member Functions

 OpRVEBCsPeriodicCalCU (const string field_name, RVEBC_Data_Periodic &data, Vec f, CommonDataPeriodic &common_data_periodic, CommonFunctionsPeriodic &common_functions_periodic, bool ho_geometry=false)
 
PetscErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 

Public Attributes

Vec F
 
RVEBC_Data_PeriodicdAta
 
bool hoGeometry
 
CommonDataPeriodiccommonDataPeriodic
 
CommonFunctionsPeriodiccommonFunctionsPeriodic
 
ublas::vector< int > row_ind
 
VectorDouble fe1
 
VectorDouble fe2
 
MatrixDouble N_mat
 

Detailed Description

Definition at line 825 of file BCs_RVELagrange_Periodic.hpp.

Constructor & Destructor Documentation

◆ OpRVEBCsPeriodicCalCU()

BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::OpRVEBCsPeriodicCalCU ( const string  field_name,
RVEBC_Data_Periodic data,
Vec  f,
CommonDataPeriodic common_data_periodic,
CommonFunctionsPeriodic common_functions_periodic,
bool  ho_geometry = false 
)
inline

Definition at line 833 of file BCs_RVELagrange_Periodic.hpp.

840  :
843  ),
844  F(f),
845  dAta(data),
846  hoGeometry(ho_geometry),
847  commonDataPeriodic(common_data_periodic),
848  commonFunctionsPeriodic(common_functions_periodic) {
849  }

Member Function Documentation

◆ doWork()

PetscErrorCode BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::doWork ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)
inline

Definition at line 857 of file BCs_RVELagrange_Periodic.hpp.

857  {
858  PetscFunctionBegin;
859 // cout<<" OpRVEBCsPeriodicCalCU "<<endl;
860  if(data.getIndices().size()==0) PetscFunctionReturn(0);
861  if(type == MBQUAD) PetscFunctionReturn(0);
862  if(type == MBEDGE && side >= 3) PetscFunctionReturn(0); //ignore second triangel edges (cananical number 6 to 8)
863  if(type == MBTRI && side == 4) PetscFunctionReturn(0); //ignore second triangel face (cananical number 4)
864 
865  int rank = data.getFieldDofs()[0]->getNbOfCoeffs();
866 
867  if(type == MBVERTEX) {
868  int nb=data.getIndices().size()/2; //for one face for nodes only (total indieces are for both faces)
869  row_ind.resize(nb); //1st face of prism -- nodes
870  for(int ii=0; ii<nb; ii++){
871  row_ind[ii]=data.getIndices()[ii];
872  }
873  }
874 
875  for(unsigned int gg = 0;gg<data.getN().size1();gg++) {
876  double area;
877  VectorDouble normal_f3 = getNormalsAtGaussPtsF3(gg);
878  area = cblas_dnrm2(3,&normal_f3[0],1)*0.5;
879  double val = getGaussPts()(2,gg)*area;
880 
881  if(type == MBVERTEX) {
882  ierr = commonFunctionsPeriodic.shapeMat(rank, gg, data, N_mat, 2); CHKERRQ(ierr);}
883  else {
884  ierr = commonFunctionsPeriodic.shapeMat(rank, gg, data, N_mat); CHKERRQ(ierr);}
885 
886  if(gg==0){
887  fe1=-val*prod(trans(N_mat), commonDataPeriodic.DispAtGaussPts[0][gg]);
888  fe2= val*prod(trans(N_mat), commonDataPeriodic.DispAtGaussPts[1][gg]);}
889  else{
890  fe1+=-val*prod(trans(N_mat), commonDataPeriodic.DispAtGaussPts[0][gg]);
891  fe2+= val*prod(trans(N_mat), commonDataPeriodic.DispAtGaussPts[1][gg]);}
892  }
893 
894  if(type == MBVERTEX) {
895  ierr = VecSetValues(F,row_ind.size(),&row_ind[0],&fe1[0],ADD_VALUES); CHKERRQ(ierr);
896  ierr = VecSetValues(F,row_ind.size(),&row_ind[0],&fe2[0],ADD_VALUES); CHKERRQ(ierr);}
897  else{
898  ierr = VecSetValues(F,data.getIndices().size(),&data.getIndices()[0],&fe1[0],ADD_VALUES); CHKERRQ(ierr);
899  ierr = VecSetValues(F,data.getIndices().size(),&data.getIndices()[0],&fe2[0],ADD_VALUES); CHKERRQ(ierr);}
900 
901  PetscFunctionReturn(0);
902  }

Member Data Documentation

◆ commonDataPeriodic

CommonDataPeriodic& BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::commonDataPeriodic

Definition at line 830 of file BCs_RVELagrange_Periodic.hpp.

◆ commonFunctionsPeriodic

CommonFunctionsPeriodic& BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::commonFunctionsPeriodic

Definition at line 831 of file BCs_RVELagrange_Periodic.hpp.

◆ dAta

RVEBC_Data_Periodic& BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::dAta

Definition at line 828 of file BCs_RVELagrange_Periodic.hpp.

◆ F

Vec BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::F

Definition at line 827 of file BCs_RVELagrange_Periodic.hpp.

◆ fe1

VectorDouble BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::fe1

Definition at line 853 of file BCs_RVELagrange_Periodic.hpp.

◆ fe2

VectorDouble BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::fe2

Definition at line 853 of file BCs_RVELagrange_Periodic.hpp.

◆ hoGeometry

bool BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::hoGeometry

Definition at line 829 of file BCs_RVELagrange_Periodic.hpp.

◆ N_mat

MatrixDouble BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::N_mat

Definition at line 854 of file BCs_RVELagrange_Periodic.hpp.

◆ row_ind

ublas::vector<int> BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::row_ind

Definition at line 852 of file BCs_RVELagrange_Periodic.hpp.


The documentation for this struct was generated from the following file:
BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::commonFunctionsPeriodic
CommonFunctionsPeriodic & commonFunctionsPeriodic
Definition: BCs_RVELagrange_Periodic.hpp:831
BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::commonDataPeriodic
CommonDataPeriodic & commonDataPeriodic
Definition: BCs_RVELagrange_Periodic.hpp:830
BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::N_mat
MatrixDouble N_mat
Definition: BCs_RVELagrange_Periodic.hpp:854
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1576
BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::hoGeometry
bool hoGeometry
Definition: BCs_RVELagrange_Periodic.hpp:829
BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::dAta
RVEBC_Data_Periodic & dAta
Definition: BCs_RVELagrange_Periodic.hpp:828
BCs_RVELagrange_Periodic::CommonDataPeriodic::DispAtGaussPts
ublas::vector< ublas::vector< VectorDouble > > DispAtGaussPts
Definition: BCs_RVELagrange_Periodic.hpp:152
convert.type
type
Definition: convert.py:64
BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::fe2
VectorDouble fe2
Definition: BCs_RVELagrange_Periodic.hpp:853
BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::F
Vec F
Definition: BCs_RVELagrange_Periodic.hpp:827
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::row_ind
ublas::vector< int > row_ind
Definition: BCs_RVELagrange_Periodic.hpp:852
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
BCs_RVELagrange_Periodic::CommonFunctionsPeriodic::shapeMat
PetscErrorCode shapeMat(int rank, unsigned int gg, DataForcesAndSourcesCore::EntData &col_data, MatrixDouble &N_mat, int div=1)
Definition: BCs_RVELagrange_Periodic.hpp:91
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
BCs_RVELagrange_Periodic::OpRVEBCsPeriodicCalCU::fe1
VectorDouble fe1
Definition: BCs_RVELagrange_Periodic.hpp:853
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567