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

\biref operator to calculate the RHS of the constrain for the RVE boundary conditions More...

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

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

Public Member Functions

 OpDmatRhs (const string field_name, const string lagrang_field_name, RVEBC_Data_Periodic &data, CommonDataPeriodic &common_data_periodic, CommonFunctionsPeriodic &common_functions_periodic, bool ho_geometry=false)
 
PetscErrorCode calculateDmat (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 

Public Attributes

RVEBC_Data_PeriodicdAta
 
bool hoGeometry
 
CommonDataPeriodiccommonDataPeriodic
 
CommonFunctionsPeriodiccommonFunctionsPeriodic
 
ublas::vector< MatrixDouble > X_mat
 
MatrixDouble N_mat
 

Detailed Description

\biref operator to calculate the RHS of the constrain for the RVE boundary conditions

Definition at line 493 of file BCs_RVELagrange_Periodic.hpp.

Constructor & Destructor Documentation

◆ OpDmatRhs()

BCs_RVELagrange_Periodic::OpDmatRhs::OpDmatRhs ( const string  field_name,
const string  lagrang_field_name,
RVEBC_Data_Periodic data,
CommonDataPeriodic common_data_periodic,
CommonFunctionsPeriodic common_functions_periodic,
bool  ho_geometry = false 
)
inline

Definition at line 500 of file BCs_RVELagrange_Periodic.hpp.

507  :
509  lagrang_field_name,UserDataOperator::OPROW
510  ),
511  dAta(data),
512  hoGeometry(ho_geometry),
513  commonDataPeriodic(common_data_periodic),
514  commonFunctionsPeriodic(common_functions_periodic) {
515  }

Member Function Documentation

◆ calculateDmat()

PetscErrorCode BCs_RVELagrange_Periodic::OpDmatRhs::calculateDmat ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)
inline

Definition at line 520 of file BCs_RVELagrange_Periodic.hpp.

520  {
521  PetscFunctionBegin;
522 // if(data.getIndices().size()==0) PetscFunctionReturn(0);
523 // if(type == MBEDGE && side >= 6) PetscFunctionReturn(0); //ignore second triangel edges (cananical number 6 to 8)
524 // if(type == MBTRI && side == 4) PetscFunctionReturn(0); //ignore second triangel face (cananical number 4)
525 
526  int rank = data.getFieldDofs()[0]->getNbOfCoeffs();
527 
528  double x,y,z,x1,y1,z1;
529  for(unsigned int gg = 0;gg<data.getN().size1();gg++) {
530  double area;
531  VectorDouble normal_f3 = getNormalsAtGaussPtsF3(gg);
532  area = cblas_dnrm2(3,&normal_f3[0],1)*0.5;
533  double val = getGaussPts()(2,gg)*area;
534 
535  // Gauss point coordinates of the prism for both faces with HO geometry
536  x = getHOCoordsAtGaussPtsF3()(gg,0);
537  y = getHOCoordsAtGaussPtsF3()(gg,1);
538  z = getHOCoordsAtGaussPtsF3()(gg,2);
539 
540  x1 = getHOCoordsAtGaussPtsF4()(gg,0);
541  y1 = getHOCoordsAtGaussPtsF4()(gg,1);
542  z1 = getHOCoordsAtGaussPtsF4()(gg,2);
543 
544 // cout<<"x "<<x << " y "<< y << " z "<< z <<endl;
545 // cout<<"x1 "<<x1 << " y1 "<< y1 << " z1 "<< z1 <<endl;
546 
547  commonDataPeriodic.D_mat.resize(2);
548  X_mat.resize(2);
549 
550  switch(rank) {
551  case 3: //mech problem
552  X_mat[0].resize(rank,6,false); X_mat[0].clear();
553  X_mat[0](0,0)=2.0*x; X_mat[0](0,3)=y; X_mat[0](0,5)=z;
554  X_mat[0](1,1)=2.0*y; X_mat[0](1,3)=x; X_mat[0](1,4)=z;
555  X_mat[0](2,2)=2.0*z; X_mat[0](2,4)=y; X_mat[0](2,5)=x;
556  X_mat[0]=0.5*X_mat[0];
557 
558  X_mat[1].resize(rank,6,false); X_mat[1].clear();
559  X_mat[1](0,0)=2.0*x1; X_mat[1](0,3)=y1; X_mat[1](0,5)=z1;
560  X_mat[1](1,1)=2.0*y1; X_mat[1](1,3)=x1; X_mat[1](1,4)=z1;
561  X_mat[1](2,2)=2.0*z1; X_mat[1](2,4)=y1; X_mat[1](2,5)=x1;
562  X_mat[1]=0.5*X_mat[1];
563  break;
564  case 1: //moisture transport or thermal problem
565  X_mat[0].resize(rank,3,false); X_mat[0].clear();
566  X_mat[0](0,0)=x; X_mat[0](0,1)=y; X_mat[0](0,2)=z;
567 
568  X_mat[1].resize(rank,3,false); X_mat[1].clear();
569  X_mat[1](0,0)=x1; X_mat[1](0,1)=y1; X_mat[1](0,2)=z1;
570  break;
571  default:
572  SETERRQ(PETSC_COMM_SELF,1,"not implemented");
573  }
574 
575  //
576  if(type == MBVERTEX) {
577  ierr = commonFunctionsPeriodic.shapeMat(rank, gg, data, N_mat, 2); CHKERRQ(ierr);
578  } else {
579  ierr = commonFunctionsPeriodic.shapeMat(rank, gg, data, N_mat); CHKERRQ(ierr);
580  }
581 
582  if(gg==0) {
583  commonDataPeriodic.D_mat[0]=-1*val*prod(trans(N_mat),X_mat[0]);
584  commonDataPeriodic.D_mat[1]= val*prod(trans(N_mat),X_mat[1]);
585  } else {
586  commonDataPeriodic.D_mat[0]+=-1*val*prod(trans(N_mat),X_mat[0]);
587  commonDataPeriodic.D_mat[1]+= val*prod(trans(N_mat),X_mat[1]);
588  }
589  }
590 
591 // cout<<endl<<endl;
592 
593  PetscFunctionReturn(0);
594  }

Member Data Documentation

◆ commonDataPeriodic

CommonDataPeriodic& BCs_RVELagrange_Periodic::OpDmatRhs::commonDataPeriodic

Definition at line 497 of file BCs_RVELagrange_Periodic.hpp.

◆ commonFunctionsPeriodic

CommonFunctionsPeriodic& BCs_RVELagrange_Periodic::OpDmatRhs::commonFunctionsPeriodic

Definition at line 498 of file BCs_RVELagrange_Periodic.hpp.

◆ dAta

RVEBC_Data_Periodic& BCs_RVELagrange_Periodic::OpDmatRhs::dAta

Definition at line 495 of file BCs_RVELagrange_Periodic.hpp.

◆ hoGeometry

bool BCs_RVELagrange_Periodic::OpDmatRhs::hoGeometry

Definition at line 496 of file BCs_RVELagrange_Periodic.hpp.

◆ N_mat

MatrixDouble BCs_RVELagrange_Periodic::OpDmatRhs::N_mat

Definition at line 518 of file BCs_RVELagrange_Periodic.hpp.

◆ X_mat

ublas::vector<MatrixDouble > BCs_RVELagrange_Periodic::OpDmatRhs::X_mat

Definition at line 517 of file BCs_RVELagrange_Periodic.hpp.


The documentation for this struct was generated from the following file:
BCs_RVELagrange_Periodic::OpDmatRhs::commonDataPeriodic
CommonDataPeriodic & commonDataPeriodic
Definition: BCs_RVELagrange_Periodic.hpp:497
BCs_RVELagrange_Periodic::OpDmatRhs::N_mat
MatrixDouble N_mat
Definition: BCs_RVELagrange_Periodic.hpp:518
convert.type
type
Definition: convert.py:64
BCs_RVELagrange_Periodic::OpDmatRhs::hoGeometry
bool hoGeometry
Definition: BCs_RVELagrange_Periodic.hpp:496
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::CommonDataPeriodic::D_mat
ublas::vector< MatrixDouble > D_mat
Definition: BCs_RVELagrange_Periodic.hpp:140
BCs_RVELagrange_Periodic::OpDmatRhs::dAta
RVEBC_Data_Periodic & dAta
Definition: BCs_RVELagrange_Periodic.hpp:495
BCs_RVELagrange_Periodic::OpDmatRhs::X_mat
ublas::vector< MatrixDouble > X_mat
Definition: BCs_RVELagrange_Periodic.hpp:517
BCs_RVELagrange_Periodic::OpDmatRhs::commonFunctionsPeriodic
CommonFunctionsPeriodic & commonFunctionsPeriodic
Definition: BCs_RVELagrange_Periodic.hpp:498
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567