v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
DirichletDisplacementRemoveDofsBc Struct Reference

Set Dirichlet boundary conditions on displacements by removing dofs. More...

#include <users_modules/basic_finite_elements/src/DirichletBC.hpp>

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

Public Member Functions

 DirichletDisplacementRemoveDofsBc (MoFEM::Interface &m_field, const std::string &field_name, const std::string &problem_name, string blockset_name="DISPLACEMENT", bool is_partitioned=false)
 
MoFEMErrorCode iNitialize ()
 
virtual boost::shared_ptr< EntityMethod > getEntMethodPtr (DataFromBc &data)
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
- Public Member Functions inherited from DirichletDisplacementBc
 DirichletDisplacementBc (MoFEM::Interface &m_field, const std::string &field_name, Mat Aij, Vec X, Vec F, string blockset_name="DISPLACEMENT")
 
 DirichletDisplacementBc (MoFEM::Interface &m_field, const std::string &field_name, string blockset_name="DISPLACEMENT")
 
virtual MoFEMErrorCode iNitialize ()
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
MoFEMErrorCode getBcDataFromSetsAndBlocks (std::vector< DataFromBc > &bc_data)
 Get the Bc Data From Sets And Blocks object Use DISPLACEMENT blockset name (default) with 6 atributes: 1,2,3 are values of displacements x,y,z 4,5,6 are flags for x,y,z (0 or 1) More...
 
MoFEMErrorCode getRotationBcFromBlock (std::vector< DataFromBc > &bc_data)
 Get the Rotation Bc From Block object Use ROTATION blockset name with 7 atributes: 1,2,3 are x,y,z coords of the center of rotation 4,5,6 are are angular velocities in x,y,z. More...
 
MoFEMErrorCode calculateRotationForDof (VectorDouble3 &coords, DataFromBc &bc_data)
 Calculate displacements from rotation for particular dof. More...
 
MoFEMErrorCode calculateRotationForDof (EntityHandle ent, DataFromBc &bc_data)
 
MoFEMErrorCode applyScaleBcData (DataFromBc &bc_data)
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 FEMethod ()=default
 
auto getFEName () const
 get finite element name More...
 
auto getDataDofsPtr () const
 
auto getDataVectorDofsPtr () const
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 
auto & getRowFieldEnts () const
 
auto & getRowFieldEntsPtr () const
 
auto & getColFieldEnts () const
 
auto & getColFieldEntsPtr () const
 
auto getRowDofsPtr () const
 
auto getColDofsPtr () const
 
auto getNumberOfNodes () const
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
- Public Member Functions inherited from MoFEM::BasicMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities. More...
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements. More...
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements. More...
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::PetscData
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
virtual MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
 
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 

Public Attributes

boost::shared_ptr< vector< DataFromBc > > bcDataPtr
 
bool isPartitioned
 
string problemName
 
- Public Attributes inherited from DirichletDisplacementBc
MoFEM::InterfacemField
 
const std::string fieldName
 field name to set Dirichlet BC More...
 
double dIag
 diagonal value set on zeroed column and rows More...
 
std::map< DofIdx, FieldData > mapZeroRows
 
std::vector< int > dofsIndices
 
std::vector< doubledofsValues
 
std::vector< doubledofsXValues
 
const std::string blocksetName
 
boost::ptr_vector< MethodForForceScalingmethodsOp
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities. More...
 
int rAnk
 processor rank More...
 
int sIze
 number of processors in communicator More...
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities More...
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities More...
 
const ProblemproblemPtr
 raw pointer to problem More...
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container More...
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities More...
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs More...
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements More...
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing. More...
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing. More...
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator. More...
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec & ksp_f
 
Mat & ksp_A
 
Mat & ksp_B
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 
Vec f
 
Mat A
 
Mat B
 
Vec x
 
Vec x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 snes solver More...
 
Vec & snes_x
 state vector More...
 
Vec & snes_f
 residual More...
 
Mat & snes_A
 jacobian matrix More...
 
Mat & snes_B
 preconditioner of jacobian matrix More...
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 time solver More...
 
TSContext ts_ctx
 
PetscInt ts_step
 time step number More...
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver) More...
 
PetscReal ts_aa
 shift for U_tt shift for U_tt More...
 
PetscReal ts_t
 time More...
 
PetscReal ts_dt
 time step size More...
 
Vec & ts_u
 state vector More...
 
Vec & ts_u_t
 time derivative of state vector More...
 
Vec & ts_u_tt
 second time derivative of state vector More...
 
Vec & ts_F
 residual vector More...
 
Mat & ts_A
 
Mat & ts_B
 Preconditioner for ts_A. More...
 

Additional Inherited Members

- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 pass information about context of KSP/DM for with finite element is computed More...
 
- Public Types inherited from MoFEM::PetscData
enum  DataContext {
  CTX_SET_NONE = 0 , CTX_SET_F = 1 << 0 , CTX_SET_A = 1 << 1 , CTX_SET_B = 1 << 2 ,
  CTX_SET_X = 1 << 3 , CTX_SET_X_T = 1 << 4 , CTX_SET_X_TT = 1 << 6 , CTX_SET_TIME = 1 << 7
}
 
using Switches = std::bitset< 8 >
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 

Detailed Description

Set Dirichlet boundary conditions on displacements by removing dofs.

Definition at line 283 of file DirichletBC.hpp.

Constructor & Destructor Documentation

◆ DirichletDisplacementRemoveDofsBc()

DirichletDisplacementRemoveDofsBc::DirichletDisplacementRemoveDofsBc ( MoFEM::Interface m_field,
const std::string &  field_name,
const std::string &  problem_name,
string  blockset_name = "DISPLACEMENT",
bool  is_partitioned = false 
)
inline

Definition at line 289 of file DirichletBC.hpp.

294 : DirichletDisplacementBc(m_field, field_name, blockset_name),
295 problemName(problem_name), isPartitioned(is_partitioned) {}
constexpr auto field_name
Set Dirichlet boundary conditions on displacements.
Definition: DirichletBC.hpp:57

Member Function Documentation

◆ getEntMethodPtr()

virtual boost::shared_ptr< EntityMethod > DirichletDisplacementRemoveDofsBc::getEntMethodPtr ( DataFromBc data)
inlinevirtual

Reimplemented in DirichletSpatialRemoveDofsBc.

Definition at line 299 of file DirichletBC.hpp.

299 {
300 return boost::make_shared<BcEntMethodDisp>(this, data);
301 }

◆ iNitialize()

MoFEMErrorCode DirichletDisplacementRemoveDofsBc::iNitialize ( )
virtual

Reimplemented from DirichletDisplacementBc.

Definition at line 647 of file DirichletBC.cpp.

647 {
649
651 auto prb_mng = mField.getInterface<ProblemsManager>();
652 auto field_ptr = mField.get_field_structure(fieldName);
653 const int nb_coefficients = field_ptr->getNbOfCoeffs();
654
655 bcDataPtr = boost::make_shared<vector<DataFromBc>>();
658
659 auto get_dim = [&](const Range &ents) {
660 for (auto d : {3, 2, 1})
661 if (ents.num_of_dimension(d))
662 return d;
663 return 0;
664 };
665
666 auto get_adj_ents = [&](const Range &ents) {
667 Range verts;
668 CHKERR mField.get_moab().get_connectivity(ents, verts, true);
669 const auto dim = get_dim(ents);
670 for (size_t d = 1; d < dim; ++d)
671 CHKERR mField.get_moab().get_adjacencies(ents, d, false, verts,
672 moab::Interface::UNION);
673 verts.merge(ents);
674 CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(verts);
675 return verts;
676 };
677
678 auto remove_dofs_from_dirichlet_data = [&]() {
680 for (auto &bc_it : *bcDataPtr)
681 for (auto &ents : bc_it.bc_ents)
682 for (int i = 0; i != nb_coefficients; i++)
683 if (bc_it.bc_flags[i])
684 CHKERR prb_mng->removeDofsOnEntities(problemName, fieldName,
685 // get_adj_ents(ents), i, i);
686 ents, i, i);
687
689 };
690
691 auto remove_dofs_from_dirichlet_data_non_distributed = [&]() {
693 for (auto &bc_it : *bcDataPtr)
694 for (auto &ents : bc_it.bc_ents)
695 for (int i = 0; i != nb_coefficients; i++)
696 if (bc_it.bc_flags[i])
697 CHKERR prb_mng->removeDofsOnEntitiesNotDistributed(
698 // problemName, fieldName, get_adj_ents(ents), i, i);
699 problemName, fieldName, ents, i, i, QUIET);
700
702 };
703
704 if (isPartitioned)
705 CHKERR remove_dofs_from_dirichlet_data();
706 else
707 CHKERR remove_dofs_from_dirichlet_data_non_distributed();
708
710}
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
@ QUIET
Definition: definitions.h:208
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
const int dim
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
FTensor::Index< 'i', SPACE_DIM > i
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: dTensor0.hpp:27
auto get_dim(const Range &ents)
Definition: BcManager.cpp:10
auto get_adj_ents(moab::Interface &moab, const Range &ents)
Definition: BcManager.cpp:17
MoFEMErrorCode getRotationBcFromBlock(std::vector< DataFromBc > &bc_data)
Get the Rotation Bc From Block object Use ROTATION blockset name with 7 atributes: 1,...
const std::string fieldName
field name to set Dirichlet BC
Definition: DirichletBC.hpp:60
MoFEMErrorCode getBcDataFromSetsAndBlocks(std::vector< DataFromBc > &bc_data)
Get the Bc Data From Sets And Blocks object Use DISPLACEMENT blockset name (default) with 6 atributes...
Definition: DirichletBC.cpp:64
MoFEM::Interface & mField
Definition: DirichletBC.hpp:59
boost::shared_ptr< vector< DataFromBc > > bcDataPtr
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
Problem manager is used to build and partition problems.
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ operator()()

MoFEMErrorCode DirichletDisplacementRemoveDofsBc::operator() ( )
inlinevirtual

function is run for every finite element

It is used to calculate element local matrices and assembly. It can be used for post-processing.

Reimplemented from DirichletDisplacementBc.

Definition at line 304 of file DirichletBC.hpp.

304{ return 0; }

◆ postProcess()

MoFEMErrorCode DirichletDisplacementRemoveDofsBc::postProcess ( )
inlinevirtual

function is run at the end of loop

It is used to assembly matrices and vectors, calculating global variables, f.e. total internal energy, ect.

Iterating over dofs: Example1 iterating over dofs in row by name of the field for(IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP(this,"DISPLACEMENT",it)) { ... }

Reimplemented from DirichletDisplacementBc.

Definition at line 305 of file DirichletBC.hpp.

305{ return 0; }

◆ preProcess()

MoFEMErrorCode DirichletDisplacementRemoveDofsBc::preProcess ( )
virtual

function is run at the beginning of loop

It is used to zeroing matrices and vectors, calculation of shape functions on reference element, preprocessing boundary conditions, etc.

Reimplemented from DirichletDisplacementBc.

Definition at line 712 of file DirichletBC.cpp.

712 {
714 if (!bcDataPtr)
716
717 for (auto &bc_it : *bcDataPtr) {
718 Range all_bc_ents;
720 auto method = getEntMethodPtr(bc_it);
721 for (auto &ents : bc_it.bc_ents)
722 all_bc_ents.merge(ents);
723 CHKERR mField.loop_entities(fieldName, *method, &all_bc_ents);
724 }
725
727}
virtual MoFEMErrorCode loop_entities(const std::string field_name, EntityMethod &method, Range const *const ents=nullptr, int verb=DEFAULT_VERBOSITY)=0
Loop over field entities.
MoFEMErrorCode applyScaleBcData(DataFromBc &bc_data)
virtual boost::shared_ptr< EntityMethod > getEntMethodPtr(DataFromBc &data)

Member Data Documentation

◆ bcDataPtr

boost::shared_ptr<vector<DataFromBc> > DirichletDisplacementRemoveDofsBc::bcDataPtr

Definition at line 285 of file DirichletBC.hpp.

◆ isPartitioned

bool DirichletDisplacementRemoveDofsBc::isPartitioned

Definition at line 286 of file DirichletBC.hpp.

◆ problemName

string DirichletDisplacementRemoveDofsBc::problemName

Definition at line 287 of file DirichletBC.hpp.


The documentation for this struct was generated from the following files: