v0.9.2
Public Member Functions | Public Attributes | List of all members
Smoother::OpLhsSmoother Struct Reference

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

Inheritance diagram for Smoother::OpLhsSmoother:
[legend]
Collaboration diagram for Smoother::OpLhsSmoother:
[legend]

Public Member Functions

 OpLhsSmoother (const std::string vel_field, const std::string field_name, NonlinearElasticElement::BlockData &data, NonlinearElasticElement::CommonData &common_data, SmootherBlockData &smoother_data, const std::string crack_area_tangent_constrains)
 
MoFEMErrorCode aSemble (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 
- Public Member Functions inherited from NonlinearElasticElement::OpLhsPiolaKirchhoff_dx
 OpLhsPiolaKirchhoff_dx (const std::string vel_field, const std::string field_name, BlockData &data, CommonData &common_data)
 
virtual MoFEMErrorCode getJac (DataForcesAndSourcesCore::EntData &col_data, int gg)
 Directive of Piola Kirchhoff stress over spatial DOFs. More...
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
MatrixDoublegetHoCoordsAtGaussPts ()
 coordinate at Gauss points (if hierarchical approximation of element geometry) More...
 
MatrixDoublegetHoGaussPtsJac ()
 
MatrixDoublegetHoGaussPtsInvJac ()
 
VectorDoublegetHoGaussPtsDetJac ()
 
auto getFTenosr0HoMeasure ()
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
auto getFTensor1HoCoordsAtGaussPts ()
 Get coordinates at integration points taking geometry from field. More...
 
auto getFTensor2HoGaussPtsJac ()
 
auto getFTensor2HoGaussPtsInvJac ()
 
VolumeElementForcesAndSourcesCoreBasegetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
MoFEMErrorCode getDivergenceOfHDivBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, VectorDouble &div)
 Get divergence of base functions at integration point. More...
 
MoFEMErrorCode getCurlOfHCurlBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, MatrixDouble &curl)
 Get curl of base functions at integration point. More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPLAST, const bool symm=true)
 
 UserDataOperator (const std::string &field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
boost::weak_ptr< SideNumber > getSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement ()
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
const std::string & getFEName () const
 Get name of the element. More...
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
DEPRECATED Vec getSnesF () const
 
DEPRECATED Vec getSnesX () const
 
DEPRECATED Mat getSnesA () const
 
DEPRECATED Mat getSnesB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTSa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
DEPRECATED MoFEMErrorCode getPorblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Public Attributes

SmootherBlockDatasmootherData
 
const std::string fieldCrackAreaTangentConstrains
 
ublas::vector< int > rowFrontIndices
 
- Public Attributes inherited from NonlinearElasticElement::OpLhsPiolaKirchhoff_dx
BlockDatadAta
 
CommonDatacommonData
 
int tAg
 
bool aLe
 
ublas::vector< int > rowIndices
 
ublas::vector< int > colIndices
 
MatrixDouble k
 
MatrixDouble trans_k
 
MatrixDouble jac
 
MatrixDouble F
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType { OPROW = 1 << 0, OPCOL = 1 << 1, OPROWCOL = 1 << 2, OPLAST = 1 << 3 }
 Controls loop over entities on element. More...
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Examples
mesh_smoothing.cpp.

Definition at line 225 of file Smoother.hpp.

Constructor & Destructor Documentation

◆ OpLhsSmoother()

Smoother::OpLhsSmoother::OpLhsSmoother ( const std::string  vel_field,
const std::string  field_name,
NonlinearElasticElement::BlockData data,
NonlinearElasticElement::CommonData common_data,
SmootherBlockData smoother_data,
const std::string  crack_area_tangent_constrains 
)

Definition at line 231 of file Smoother.hpp.

239  : NonlinearElasticElement::OpLhsPiolaKirchhoff_dx(vel_field, field_name,
240  data, common_data),
241  smootherData(smoother_data),
242  fieldCrackAreaTangentConstrains(crack_area_tangent_constrains) {}
SmootherBlockData & smootherData
Definition: Smoother.hpp:228
const std::string fieldCrackAreaTangentConstrains
Definition: Smoother.hpp:229

Member Function Documentation

◆ aSemble()

MoFEMErrorCode Smoother::OpLhsSmoother::aSemble ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
DataForcesAndSourcesCore::EntData row_data,
DataForcesAndSourcesCore::EntData col_data 
)
virtual

Reimplemented from NonlinearElasticElement::OpLhsPiolaKirchhoff_dx.

Definition at line 246 of file Smoother.hpp.

249  {
251 
252  int nb_row = row_data.getIndices().size();
253  int nb_col = col_data.getIndices().size();
254  int *row_indices_ptr = &row_data.getIndices()[0];
255  int *col_indices_ptr = &col_data.getIndices()[0];
256 
257  if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
258  rowIndices.resize(nb_row, false);
259  noalias(rowIndices) = row_data.getIndices();
260  if (!smootherData.sTabilised) {
261  row_indices_ptr = &rowIndices[0];
262  }
263  rowFrontIndices.resize(nb_row, false);
264  noalias(rowFrontIndices) = row_data.getIndices();
265  VectorDofs &dofs = row_data.getFieldDofs();
266  VectorDofs::iterator dit = dofs.begin();
267  for (int ii = 0; dit != dofs.end(); dit++, ii++) {
268  if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) !=
270  rowIndices[ii] = -1;
271  } else {
272  rowFrontIndices[ii] = -1;
273  }
274  }
275  }
276 
277  CHKERR MatSetValues(getFEMethod()->snes_B, nb_row, row_indices_ptr,
278  nb_col, col_indices_ptr, &k(0, 0), ADD_VALUES);
279 
281 
282  // get tangent vector array
283  double *f_tangent_front_mesh_array;
284  CHKERR VecGetArray(smootherData.tangentFrontF,
285  &f_tangent_front_mesh_array);
286  // iterate nodes on tet
287  for (int nn = 0; nn < 4; nn++) {
288 
289  // get indices with Lagrange multiplier at node nn
290  FENumeredDofEntityByNameAndEnt::iterator dit, hi_dit;
291  dit = getFEMethod()
292  ->rowPtr->get<Composite_Name_And_Ent_mi_tag>()
293  .lower_bound(boost::make_tuple(
295  hi_dit = getFEMethod()
296  ->rowPtr->get<Composite_Name_And_Ent_mi_tag>()
297  .upper_bound(boost::make_tuple(
299 
300  // continue if Lagrange are on element
301  if (std::distance(dit, hi_dit) > 0) {
302 
303  FENumeredDofEntityByNameAndEnt::iterator diit, hi_diit;
304 
305  // get mesh node positions at node nn
306  diit = getFEMethod()
307  ->rowPtr->get<Composite_Name_And_Ent_mi_tag>()
308  .lower_bound(boost::make_tuple("MESH_NODE_POSITIONS",
309  getConn()[nn]));
310  hi_diit = getFEMethod()
311  ->rowPtr->get<Composite_Name_And_Ent_mi_tag>()
312  .upper_bound(boost::make_tuple("MESH_NODE_POSITIONS",
313  getConn()[nn]));
314 
315  // iterate over dofs on node nn
316  for (; diit != hi_diit; diit++) {
317  // iterate overt dofs in element column
318  for (int ddd = 0; ddd < nb_col; ddd++) {
319  // check consistency, node has to be at crack front
320  if (rowFrontIndices[3 * nn + diit->get()->getDofCoeffIdx()] !=
321  diit->get()->getPetscGlobalDofIdx()) {
322  SETERRQ2(
323  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
324  "data inconsistency %d != %d",
325  rowFrontIndices[3 * nn + diit->get()->getDofCoeffIdx()],
326  diit->get()->getPetscGlobalDofIdx());
327  }
328  // dof is not on this partition
329  if (diit->get()->getPetscLocalDofIdx() == -1) {
330  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
331  "data inconsistency");
332  }
333  double g =
334  f_tangent_front_mesh_array[diit->get()
335  ->getPetscLocalDofIdx()] *
336  k(3 * nn + diit->get()->getDofCoeffIdx(), ddd);
337  int lambda_idx = dit->get()->getPetscGlobalDofIdx();
338  CHKERR MatSetValues(getFEMethod()->snes_B, 1, &lambda_idx, 1,
339  &col_indices_ptr[ddd], &g, ADD_VALUES);
340  }
341  }
342  }
343  }
344  CHKERR VecRestoreArray(smootherData.tangentFrontF,
345  &f_tangent_front_mesh_array);
346  }
347 
349  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
ublas::vector< int > rowFrontIndices
Definition: Smoother.hpp:244
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MoFEMErrorCode MatSetValues(Mat M, const DataForcesAndSourcesCore::EntData &row_data, const DataForcesAndSourcesCore::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
#define CHKERR
Inline error check.
Definition: definitions.h:602
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
SmootherBlockData & smootherData
Definition: Smoother.hpp:228
const std::string fieldCrackAreaTangentConstrains
Definition: Smoother.hpp:229
boost::shared_ptr< const FENumeredDofEntity_multiIndex > rowPtr
Pointer to finite element rows dofs view.

Member Data Documentation

◆ fieldCrackAreaTangentConstrains

const std::string Smoother::OpLhsSmoother::fieldCrackAreaTangentConstrains

Definition at line 229 of file Smoother.hpp.

◆ rowFrontIndices

ublas::vector<int> Smoother::OpLhsSmoother::rowFrontIndices

Definition at line 244 of file Smoother.hpp.

◆ smootherData

SmootherBlockData& Smoother::OpLhsSmoother::smootherData

Definition at line 228 of file Smoother.hpp.


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