v0.9.0
Public Member Functions | Public Attributes | List of all members
DirichletFixFieldAtEntitiesBc Struct Reference

Fix dofs on entities. More...

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

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

Public Member Functions

 DirichletFixFieldAtEntitiesBc (MoFEM::Interface &m_field, const std::string &field_name, Mat aij, Vec x, Vec f, Range &ents)
 
 DirichletFixFieldAtEntitiesBc (MoFEM::Interface &m_field, const std::string &field_name, Range &ents)
 
MoFEMErrorCode iNitalize ()
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop 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)
 
 DirichletDisplacementBc (MoFEM::Interface &m_field, const std::string &field_name)
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
virtual MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()
 
MoFEMErrorCode setKspCtx (const KSPContext &ctx)
 set operator type More...
 
MoFEMErrorCode setKsp (KSP ksp)
 set solver More...
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of 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
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()
 
MoFEMErrorCode setSnesCtx (const SNESContext &ctx)
 Set SNES context. More...
 
MoFEMErrorCode setSnes (SNES snes)
 Set SNES instance. More...
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. More...
 

Public Attributes

Range eNts
 
std::vector< std::string > fieldNames
 
- 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< intdofsIndices
 
std::vector< doubledofsValues
 
std::vector< doubledofsXValues
 
boost::ptr_vector< MethodForForceScalingmethodsOp
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 Pointer to finite element rows dofs view. More...
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 Pointer to finite element columns dofs view. More...
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 Pointer to finite element data dofs. More...
 
boost::shared_ptr< const FieldEntity_vector_viewrowFieldEntsPtr
 Pointer to finite element field entities row view. More...
 
boost::shared_ptr< const FieldEntity_vector_viewcolFieldEntsPtr
 Pointer to finite element field entities column view. More...
 
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_viewdataFieldEntsPtr
 Pointer to finite element field entities data view. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process 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
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec ksp_f
 the right hand side vector More...
 
Mat ksp_A
 matrix More...
 
Mat ksp_B
 preconditioner matrix More...
 
- 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
TSContext ts_ctx
 
TS ts
 time solver 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...
 
PetscInt ts_step
 time step More...
 
PetscReal ts_a
 shift for U_tt (see PETSc Time Solver) More...
 
PetscReal ts_v
 shift for U_t shift for U_t More...
 
PetscReal ts_t
 time 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::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
}
 
- Protected Member Functions inherited from MoFEM::UnknownInterface
boost::typeindex::type_index getClassIdx (const MOFEMuuid &uid) const
 Get type name for interface Id. More...
 
MOFEMuuid getUId (const boost::typeindex::type_index &class_idx) const
 Get interface Id for class name. More...
 

Detailed Description

Fix dofs on entities.

Examples
mesh_smoothing.cpp, and minimal_surface_area.cpp.

Definition at line 120 of file DirichletBC.hpp.

Constructor & Destructor Documentation

◆ DirichletFixFieldAtEntitiesBc() [1/2]

DirichletFixFieldAtEntitiesBc::DirichletFixFieldAtEntitiesBc ( MoFEM::Interface m_field,
const std::string &  field_name,
Mat  aij,
Vec  x,
Vec  f,
Range &  ents 
)

Definition at line 124 of file DirichletBC.hpp.

127  : DirichletDisplacementBc(m_field, field_name, aij, x, f), eNts(ents) {
128  fieldNames.push_back(fieldName);
129  }
const std::string fieldName
field name to set Dirichlet BC
Definition: DirichletBC.hpp:46
std::vector< std::string > fieldNames
DirichletDisplacementBc(MoFEM::Interface &m_field, const std::string &field_name, Mat Aij, Vec X, Vec F)
Definition: DirichletBC.cpp:50

◆ DirichletFixFieldAtEntitiesBc() [2/2]

DirichletFixFieldAtEntitiesBc::DirichletFixFieldAtEntitiesBc ( MoFEM::Interface m_field,
const std::string &  field_name,
Range &  ents 
)

Definition at line 131 of file DirichletBC.hpp.

133  : DirichletDisplacementBc(m_field, field_name), eNts(ents) {
134  fieldNames.push_back(fieldName);
135  }
const std::string fieldName
field name to set Dirichlet BC
Definition: DirichletBC.hpp:46
std::vector< std::string > fieldNames
DirichletDisplacementBc(MoFEM::Interface &m_field, const std::string &field_name, Mat Aij, Vec X, Vec F)
Definition: DirichletBC.cpp:50

Member Function Documentation

◆ iNitalize()

MoFEMErrorCode DirichletFixFieldAtEntitiesBc::iNitalize ( )
virtual

Reimplemented from DirichletDisplacementBc.

Definition at line 512 of file DirichletBC.cpp.

512  {
514  if (mapZeroRows.empty()) {
515 
516  for (const auto &field_name : fieldNames) {
517 
518  auto for_each_dof = [&](auto &dof) {
520  mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
522  };
523 
525  for_each_dof);
526  }
527 
528  dofsIndices.resize(mapZeroRows.size());
529  dofsValues.resize(mapZeroRows.size());
530  int ii = 0;
531  for (auto mit = mapZeroRows.begin(); mit != mapZeroRows.end();
532  ++mit, ++ii) {
533  dofsIndices[ii] = mit->first;
534  dofsValues[ii] = mit->second;
535  }
536 
537  }
539 }
std::vector< std::string > fieldNames
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
const Problem * problemPtr
raw pointer to problem
std::vector< int > dofsIndices
Definition: DirichletBC.hpp:55
std::vector< double > dofsValues
Definition: DirichletBC.hpp:56
#define CHKERR
Inline error check.
Definition: definitions.h:596
static MoFEMErrorCode set_numered_dofs_on_ents(const Problem *problem_ptr, const string &field_name, Range &ents, boost::function< MoFEMErrorCode(const boost::shared_ptr< MoFEM::NumeredDofEntity > &dof)> for_each_dof)
Definition: DirichletBC.cpp:23
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
std::map< DofIdx, FieldData > mapZeroRows
Definition: DirichletBC.hpp:54

◆ postProcess()

MoFEMErrorCode DirichletFixFieldAtEntitiesBc::postProcess ( )
virtual

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 564 of file DirichletBC.cpp.

564  {
566  if (snes_ctx == CTX_SNESNONE && ts_ctx == CTX_TSNONE) {
567  if (snes_B) {
568  CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
569  CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
570  CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
571  dofsIndices.empty() ? PETSC_NULL
572  : &*dofsIndices.begin(),
573  dIag, PETSC_NULL, PETSC_NULL);
574  }
575  if (snes_f) {
576  CHKERR VecAssemblyBegin(snes_f);
577  CHKERR VecAssemblyEnd(snes_f);
578  int ii = 0;
579  for (std::vector<int>::iterator vit = dofsIndices.begin();
580  vit != dofsIndices.end(); vit++, ii++) {
581  CHKERR VecSetValue(snes_f, *vit, dofsValues[ii], INSERT_VALUES);
582  }
583  CHKERR VecAssemblyBegin(snes_f);
584  CHKERR VecAssemblyEnd(snes_f);
585  }
586  }
587 
588  switch (snes_ctx) {
589  case CTX_SNESNONE: {
590  } break;
591  case CTX_SNESSETFUNCTION: {
592  CHKERR VecAssemblyBegin(snes_f);
593  CHKERR VecAssemblyEnd(snes_f);
594  int ii = 0;
595  for (std::vector<int>::iterator vit = dofsIndices.begin();
596  vit != dofsIndices.end(); vit++, ii++) {
597  CHKERR VecSetValue(snes_f, *vit, dofsValues[ii], INSERT_VALUES);
598  }
599  CHKERR VecAssemblyBegin(snes_f);
600  CHKERR VecAssemblyEnd(snes_f);
601  } break;
602  case CTX_SNESSETJACOBIAN: {
603  CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
604  CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
605  CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
606  dofsIndices.empty() ? PETSC_NULL
607  : &*dofsIndices.begin(),
608  dIag, PETSC_NULL, PETSC_NULL);
609  } break;
610  default:
611  SETERRQ(PETSC_COMM_SELF, 1, "unknown snes stage");
612  }
613 
615 }
Vec snes_f
residual
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
std::vector< int > dofsIndices
Definition: DirichletBC.hpp:55
double dIag
diagonal value set on zeroed column and rows
Definition: DirichletBC.hpp:47
std::vector< double > dofsValues
Definition: DirichletBC.hpp:56
#define CHKERR
Inline error check.
Definition: definitions.h:596
Mat snes_B
preconditioner of jacobian matrix
SNESContext snes_ctx
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
TSContext ts_ctx

◆ preProcess()

MoFEMErrorCode DirichletFixFieldAtEntitiesBc::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 541 of file DirichletBC.cpp.

541  {
543 
544  switch (ts_ctx) {
545  case CTX_TSSETIFUNCTION: {
547  snes_x = ts_u;
548  snes_f = ts_F;
549  break;
550  }
551  case CTX_TSSETIJACOBIAN: {
553  snes_B = ts_B;
554  break;
555  }
556  default:
557  break;
558  }
559 
560  CHKERR iNitalize();
562 }
Vec snes_f
residual
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
Vec ts_F
residual vector
Vec snes_x
state vector
#define CHKERR
Inline error check.
Definition: definitions.h:596
Vec ts_u
state vector
Mat snes_B
preconditioner of jacobian matrix
SNESContext snes_ctx
Mat ts_B
Preconditioner for ts_A.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
TSContext ts_ctx

Member Data Documentation

◆ eNts

Range DirichletFixFieldAtEntitiesBc::eNts

Definition at line 122 of file DirichletBC.hpp.

◆ fieldNames

std::vector<std::string> DirichletFixFieldAtEntitiesBc::fieldNames

Definition at line 123 of file DirichletBC.hpp.


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