v0.14.0
Classes | Public Member Functions | Public Attributes | List of all members
BCs_RVELagrange_Disp Struct Reference

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

Inheritance diagram for BCs_RVELagrange_Disp:
[legend]
Collaboration diagram for BCs_RVELagrange_Disp:
[legend]

Classes

struct  CommonFunctions
 
struct  MyTriFE
 
struct  OpDmatRhs
 
struct  OpRVEBCsLhs
 \biref operator to calculate the LHS for the RVE bounary conditions More...
 
struct  OpRVEBCsRhs
 \biref operator to calculate the RHS of the constrain for the RVE boundary conditions More...
 
struct  OpRVEBCsRhsForceExternal_CU
 \biref operator to calculate and assemble CU More...
 
struct  OpRVEBCsRhsHomoC
 \biref operator to calculate the RHS for the calculation of the homgoensied stiffness matrix More...
 
struct  OpRVEHomoStress
 \biref operator to calculate the RVE homogenised stress More...
 
struct  RVEBC_Data
 

Public Member Functions

MyTriFEgetLoopFeRVEBCLhs ()
 
MyTriFEgetLoopFeRVEBCRhs ()
 
MyTriFEgetLoopFeRVEBCRhsResidual ()
 
MyTriFEgetLoopFeRVEBCStress ()
 
MyTriFEgetLoopFeRVEBCRhsHomoC ()
 
 BCs_RVELagrange_Disp (MoFEM::Interface &m_field)
 
PetscErrorCode addLagrangiangElement (const string element_name, const string field_name, const string lagrang_field_name, const string mesh_nodals_positions)
 
PetscErrorCode setRVEBCsOperatorsNonlinear (string field_name, string lagrang_field_name, string mesh_nodals_positions, Mat aij, vector< Vec > &fvec, Vec f, VectorDouble given_strain)
 
PetscErrorCode setRVEBCsOperators (string field_name, string lagrang_field_name, string mesh_nodals_positions, Mat aij, vector< Vec > &fvec)
 
PetscErrorCode setRVEBCsHomoStressOperators (string field_name, string lagrang_field_name, string mesh_nodals_positions, Vec Stress_Homo)
 

Public Attributes

boost::ptr_vector< MethodForForceScalingmethodsOp
 
MoFEM::InterfacemField
 
MyTriFE feRVEBCRhs
 
MyTriFE feRVEBCLhs
 
MyTriFE feRVEBCStress
 
MyTriFE feRVEBCRhsResidual
 
MyTriFE feRVEBCRhsHomoC
 
map< int, RVEBC_DatasetOfRVEBC
 maps side set id with appropriate FluxData More...
 
CommonFunctions common_functions
 

Detailed Description

Definition at line 18 of file BCs_RVELagrange_Disp.hpp.

Constructor & Destructor Documentation

◆ BCs_RVELagrange_Disp()

BCs_RVELagrange_Disp::BCs_RVELagrange_Disp ( MoFEM::Interface m_field)
inline

Definition at line 44 of file BCs_RVELagrange_Disp.hpp.

44  :
45  mField(m_field),
46  feRVEBCRhs(m_field),
47  feRVEBCLhs(m_field),
48  feRVEBCRhsResidual(m_field),
49  feRVEBCStress(m_field),
50  feRVEBCRhsHomoC(m_field){
51  }

Member Function Documentation

◆ addLagrangiangElement()

PetscErrorCode BCs_RVELagrange_Disp::addLagrangiangElement ( const string  element_name,
const string  field_name,
const string  lagrang_field_name,
const string  mesh_nodals_positions 
)
inline

Definition at line 58 of file BCs_RVELagrange_Disp.hpp.

63  {
64  PetscFunctionBegin;
65 
66 
67 
68 
69  ierr = mField.add_finite_element(element_name,MF_ZERO); CHKERRQ(ierr);
70 
71  //============================================================================================================
72  //C row as Lagrange_mul_disp and col as DISPLACEMENT
73  ierr = mField.modify_finite_element_add_field_row(element_name,lagrang_field_name); CHKERRQ(ierr);
75  //CT col as Lagrange_mul_disp and row as DISPLACEMENT
76  ierr = mField.modify_finite_element_add_field_col(element_name,lagrang_field_name); CHKERRQ(ierr);
78  //data
79  ierr = mField.modify_finite_element_add_field_data(element_name,lagrang_field_name); CHKERRQ(ierr);
81  //============================================================================================================
82 
83  if(mField.check_field(mesh_nodals_positions)) { //for high order geometry
84  ierr = mField.modify_finite_element_add_field_data(element_name,mesh_nodals_positions); CHKERRQ(ierr);
85  }
86 
87  //this is alternative method for setting boundary conditions, to bypass bu in cubit file reader.
88  //not elegant, but good enough
90  if(it->getName().compare(0,12,"AllBoundSurf") == 0 || it->getMeshsetId() == 103) {
91  // cout<<"Hi from addLagrangiangElement "<<endl;
92  rval = mField.get_moab().get_entities_by_type(it->meshset,MBTRI,setOfRVEBC[it->getMeshsetId()].tRis,true); CHKERRQ_MOAB(rval);
93  ierr = mField.add_ents_to_finite_element_by_type(setOfRVEBC[it->getMeshsetId()].tRis,MBTRI,element_name); CHKERRQ(ierr);
94  }
95  }
96 
97  PetscFunctionReturn(0);
98  }

◆ getLoopFeRVEBCLhs()

MyTriFE& BCs_RVELagrange_Disp::getLoopFeRVEBCLhs ( )
inline

Definition at line 35 of file BCs_RVELagrange_Disp.hpp.

35 { return feRVEBCLhs; }

◆ getLoopFeRVEBCRhs()

MyTriFE& BCs_RVELagrange_Disp::getLoopFeRVEBCRhs ( )
inline

Definition at line 36 of file BCs_RVELagrange_Disp.hpp.

36 { return feRVEBCRhs; }

◆ getLoopFeRVEBCRhsHomoC()

MyTriFE& BCs_RVELagrange_Disp::getLoopFeRVEBCRhsHomoC ( )
inline

Definition at line 41 of file BCs_RVELagrange_Disp.hpp.

41 { return feRVEBCRhsHomoC; }

◆ getLoopFeRVEBCRhsResidual()

MyTriFE& BCs_RVELagrange_Disp::getLoopFeRVEBCRhsResidual ( )
inline

Definition at line 38 of file BCs_RVELagrange_Disp.hpp.

38 { return feRVEBCRhsResidual; }

◆ getLoopFeRVEBCStress()

MyTriFE& BCs_RVELagrange_Disp::getLoopFeRVEBCStress ( )
inline

Definition at line 39 of file BCs_RVELagrange_Disp.hpp.

39 { return feRVEBCStress; }

◆ setRVEBCsHomoStressOperators()

PetscErrorCode BCs_RVELagrange_Disp::setRVEBCsHomoStressOperators ( string  field_name,
string  lagrang_field_name,
string  mesh_nodals_positions,
Vec  Stress_Homo 
)
inline

Definition at line 671 of file BCs_RVELagrange_Disp.hpp.

676  {
677  PetscFunctionBegin;
678 
679  bool ho_geometry = false;
680  if(mField.check_field(mesh_nodals_positions)) {
681  ho_geometry = true;
682  }
683 
684  map<int,RVEBC_Data>::iterator sit = setOfRVEBC.begin();
685  for(;sit!=setOfRVEBC.end();sit++) {
686  feRVEBCStress.getOpPtrVector().push_back(new OpRVEHomoStress(field_name, lagrang_field_name, Stress_Homo, sit->second, ho_geometry));
687  }
688 
689  PetscFunctionReturn(0);
690  }

◆ setRVEBCsOperators()

PetscErrorCode BCs_RVELagrange_Disp::setRVEBCsOperators ( string  field_name,
string  lagrang_field_name,
string  mesh_nodals_positions,
Mat  aij,
vector< Vec > &  fvec 
)
inline

Definition at line 580 of file BCs_RVELagrange_Disp.hpp.

585  {
586  PetscFunctionBegin;
587 
588  bool ho_geometry = false;
589  if(mField.check_field(mesh_nodals_positions)) {
590  ho_geometry = true;
591  }
592  map<int,RVEBC_Data>::iterator sit = setOfRVEBC.begin();
593  for(;sit!=setOfRVEBC.end();sit++) {
594 
595  //Lhs (calculate LHS)
596  feRVEBCLhs.getOpPtrVector().push_back(
597  new OpRVEBCsLhs(field_name,lagrang_field_name, aij, sit->second,ho_geometry)
598  );
599 
600  //Rhs (F[6] used to calculate Homgenised stiffness matrix)
601  feRVEBCRhs.getOpPtrVector().push_back(
602  new OpRVEBCsRhsHomoC(field_name,lagrang_field_name,fvec,sit->second,ho_geometry)
603  );
604  }
605  PetscFunctionReturn(0);
606  }

◆ setRVEBCsOperatorsNonlinear()

PetscErrorCode BCs_RVELagrange_Disp::setRVEBCsOperatorsNonlinear ( string  field_name,
string  lagrang_field_name,
string  mesh_nodals_positions,
Mat  aij,
vector< Vec > &  fvec,
Vec  f,
VectorDouble  given_strain 
)
inline

Definition at line 535 of file BCs_RVELagrange_Disp.hpp.

540  {
541  PetscFunctionBegin;
542 
543  bool ho_geometry = false;
544  if(mField.check_field(mesh_nodals_positions)) {
545  ho_geometry = true;
546  }
547  map<int,RVEBC_Data>::iterator sit = setOfRVEBC.begin();
548  for(;sit!=setOfRVEBC.end();sit++) {
549 
550  //Lhs (calculate LHS)
551  feRVEBCLhs.getOpPtrVector().push_back(
552  new OpRVEBCsLhs(field_name,lagrang_field_name, aij, sit->second,ho_geometry)
553  );
554 
555 
556  //Rhs (calculate RHS)
557  feRVEBCRhs.getOpPtrVector().push_back(
558  new OpRVEBCsRhs(field_name,lagrang_field_name,f,given_strain, methodsOp, sit->second,ho_geometry)
559  );
560 
561 
562  //Rhs (F[6] used to calculate Homgenised stiffness matrix)
563  feRVEBCRhsHomoC.getOpPtrVector().push_back(
564  new OpRVEBCsRhsHomoC(field_name,lagrang_field_name,fvec,sit->second,ho_geometry)
565  );
566 
567  //for internal forces (i.e. RHS)
568  //calculate and assemble CU and CT Lam
569  feRVEBCRhsResidual.getOpPtrVector().push_back(
570  new OpRVEBCsRhsForceExternal_CU(field_name,lagrang_field_name,f,common_functions,sit->second,ho_geometry)
571  );
572 
573 
574  }
575 
576  PetscFunctionReturn(0);
577  }

Member Data Documentation

◆ common_functions

CommonFunctions BCs_RVELagrange_Disp::common_functions

Definition at line 434 of file BCs_RVELagrange_Disp.hpp.

◆ feRVEBCLhs

MyTriFE BCs_RVELagrange_Disp::feRVEBCLhs

Definition at line 30 of file BCs_RVELagrange_Disp.hpp.

◆ feRVEBCRhs

MyTriFE BCs_RVELagrange_Disp::feRVEBCRhs

Definition at line 29 of file BCs_RVELagrange_Disp.hpp.

◆ feRVEBCRhsHomoC

MyTriFE BCs_RVELagrange_Disp::feRVEBCRhsHomoC

Definition at line 33 of file BCs_RVELagrange_Disp.hpp.

◆ feRVEBCRhsResidual

MyTriFE BCs_RVELagrange_Disp::feRVEBCRhsResidual

Definition at line 32 of file BCs_RVELagrange_Disp.hpp.

◆ feRVEBCStress

MyTriFE BCs_RVELagrange_Disp::feRVEBCStress

Definition at line 31 of file BCs_RVELagrange_Disp.hpp.

◆ methodsOp

boost::ptr_vector<MethodForForceScaling> BCs_RVELagrange_Disp::methodsOp

Definition at line 25 of file BCs_RVELagrange_Disp.hpp.

◆ mField

MoFEM::Interface& BCs_RVELagrange_Disp::mField

Definition at line 28 of file BCs_RVELagrange_Disp.hpp.

◆ setOfRVEBC

map<int,RVEBC_Data> BCs_RVELagrange_Disp::setOfRVEBC

maps side set id with appropriate FluxData

Definition at line 56 of file BCs_RVELagrange_Disp.hpp.


The documentation for this struct was generated from the following file:
BCs_RVELagrange_Disp::feRVEBCRhs
MyTriFE feRVEBCRhs
Definition: BCs_RVELagrange_Disp.hpp:29
SIDESET
@ SIDESET
Definition: definitions.h:147
CHKERRQ_MOAB
#define CHKERRQ_MOAB(a)
check error code of MoAB function
Definition: definitions.h:454
MoFEM::CoreInterface::modify_finite_element_add_field_row
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
BCs_RVELagrange_Disp::feRVEBCRhsHomoC
MyTriFE feRVEBCRhsHomoC
Definition: BCs_RVELagrange_Disp.hpp:33
BCs_RVELagrange_Disp::feRVEBCStress
MyTriFE feRVEBCStress
Definition: BCs_RVELagrange_Disp.hpp:31
MoFEM::CoreInterface::add_ents_to_finite_element_by_type
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM::CoreInterface::modify_finite_element_add_field_col
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
BCs_RVELagrange_Disp::methodsOp
boost::ptr_vector< MethodForForceScaling > methodsOp
Definition: BCs_RVELagrange_Disp.hpp:25
MoFEM::CoreInterface::check_field
virtual bool check_field(const std::string &name) const =0
check if field is in database
BCs_RVELagrange_Disp::common_functions
CommonFunctions common_functions
Definition: BCs_RVELagrange_Disp.hpp:434
BCs_RVELagrange_Disp::setOfRVEBC
map< int, RVEBC_Data > setOfRVEBC
maps side set id with appropriate FluxData
Definition: BCs_RVELagrange_Disp.hpp:56
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MF_ZERO
@ MF_ZERO
Definition: definitions.h:98
BCs_RVELagrange_Disp::feRVEBCRhsResidual
MyTriFE feRVEBCRhsResidual
Definition: BCs_RVELagrange_Disp.hpp:32
_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
Definition: MeshsetsManager.hpp:71
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
BCs_RVELagrange_Disp::mField
MoFEM::Interface & mField
Definition: BCs_RVELagrange_Disp.hpp:28
BCs_RVELagrange_Disp::feRVEBCLhs
MyTriFE feRVEBCLhs
Definition: BCs_RVELagrange_Disp.hpp:30