v0.15.0
Loading...
Searching...
No Matches
Smoother::OpLhsSmoother Struct Reference

#include "tools/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, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 
 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, EntitiesFieldData::EntData &row_data, EntitiesFieldData::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 (EntitiesFieldData::EntData &col_data, int gg)
 Directive of Piola Kirchhoff stress over spatial DOFs.
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes
 
const EntityHandlegetConn ()
 get element connectivity
 
double getVolume () const
 element volume (linear geometry)
 
doublegetVolume ()
 element volume (linear geometry)
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian
 
VectorDoublegetCoords ()
 nodal coordinates
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, 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.
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle.
 
int getFEDim () const
 Get dimension of finite element.
 
EntityType getFEType () const
 Get dimension of finite element.
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer.
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity.
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element.
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices.
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices.
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object.
 
int getOpType () const
 Get operator types.
 
void setOpType (const OpType type)
 Set operator type.
 
void addOpType (const OpType type)
 Add operator type.
 
int getNinTheLoop () const
 get number of finite element in the loop
 
int getLoopSize () const
 get size of elements in the loop
 
std::string getFEName () const
 Get name of the element.
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
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
 
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 getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element
 
auto getFTensor0IntegrationWeight ()
 Get integration weights.
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3)
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry.
 
double getMeasure () const
 get measure of element
 
doublegetMeasure ()
 get measure of element
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, boost::shared_ptr< Range > fe_range=nullptr, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopRange (const string &fe_name, ForcesAndSourcesCore *range_fe, boost::shared_ptr< Range > fe_range, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 Iterate over range of elements.
 
- 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, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side.
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side.
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not.
 
void setSymm ()
 set if operator is executed taking in account symmetry
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem
 

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
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure.
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity.
 
booldoVertices
 \deprectaed If false skip vertices
 
booldoEdges
 \deprectaed If false skip edges
 
booldoQuads
 \deprectaed
 
booldoTris
 \deprectaed
 
booldoTets
 \deprectaed
 
booldoPrisms
 \deprectaed
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType
 
using DoWorkRhsHookFunType
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Examples
mesh_smoothing.cpp.

Definition at line 212 of file Smoother.hpp.

Constructor & Destructor Documentation

◆ OpLhsSmoother() [1/2]

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 )
inline

Definition at line 218 of file Smoother.hpp.

227 data, common_data),
228 smootherData(smoother_data),
229 fieldCrackAreaTangentConstrains(crack_area_tangent_constrains) {}
constexpr auto field_name
const std::string fieldCrackAreaTangentConstrains
Definition Smoother.hpp:216
SmootherBlockData & smootherData
Definition Smoother.hpp:215

◆ OpLhsSmoother() [2/2]

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 )
inline

Definition at line 218 of file Smoother.hpp.

227 data, common_data),
228 smootherData(smoother_data),
229 fieldCrackAreaTangentConstrains(crack_area_tangent_constrains) {}

Member Function Documentation

◆ aSemble() [1/2]

MoFEMErrorCode Smoother::OpLhsSmoother::aSemble ( int row_side,
int col_side,
EntityType row_type,
EntityType col_type,
EntitiesFieldData::EntData & row_data,
EntitiesFieldData::EntData & col_data )
inlinevirtual

Reimplemented from NonlinearElasticElement::OpLhsPiolaKirchhoff_dx.

Definition at line 233 of file Smoother.hpp.

236 {
238
239 int nb_row = row_data.getIndices().size();
240 int nb_col = col_data.getIndices().size();
241 int *row_indices_ptr = &row_data.getIndices()[0];
242 int *col_indices_ptr = &col_data.getIndices()[0];
243
244 if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
245 rowIndices.resize(nb_row, false);
246 noalias(rowIndices) = row_data.getIndices();
248 row_indices_ptr = &rowIndices[0];
249 }
250 rowFrontIndices.resize(nb_row, false);
251 noalias(rowFrontIndices) = row_data.getIndices();
252 VectorDofs &dofs = row_data.getFieldDofs();
253 VectorDofs::iterator dit = dofs.begin();
254 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
255 if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) !=
257 rowIndices[ii] = -1;
258 } else {
259 rowFrontIndices[ii] = -1;
260 }
261 }
262 }
263
264 CHKERR MatSetValues(getFEMethod()->snes_B, nb_row, row_indices_ptr,
265 nb_col, col_indices_ptr, &k(0, 0), ADD_VALUES);
266
268
269 const auto bit_number_for_crack_area_tangent_constrain =
271 const auto bit_number_for_mesh_position =
272 getFEMethod()->getFieldBitNumber("MESH_NODE_POSITIONS");
273
274 // get tangent vector array
275 double *f_tangent_front_mesh_array;
277 &f_tangent_front_mesh_array);
278
279 auto row_dofs = getFEMethod()->getRowDofsPtr();
280
281 // iterate nodes on tet
282 for (int nn = 0; nn < 4; nn++) {
283
284 // get indices with Lagrange multiplier at node nn
285 auto dit = row_dofs->get<Unique_mi_tag>().lower_bound(
286 DofEntity::getLoFieldEntityUId(
287 bit_number_for_crack_area_tangent_constrain, getConn()[nn]));
288 auto hi_dit = row_dofs->get<Unique_mi_tag>().upper_bound(
289 DofEntity::getHiFieldEntityUId(
290 bit_number_for_crack_area_tangent_constrain, getConn()[nn]));
291
292 // continue if Lagrange are on element
293 if (std::distance(dit, hi_dit) > 0) {
294
295 // get mesh node positions at node nn
296 auto diit = row_dofs->get<Unique_mi_tag>().lower_bound(
297 DofEntity::getLoFieldEntityUId(bit_number_for_mesh_position,
298 getConn()[nn]));
299
300 auto hi_diit = row_dofs->get<Unique_mi_tag>().upper_bound(
301 DofEntity::getHiFieldEntityUId(bit_number_for_mesh_position,
302 getConn()[nn]));
303
304 // iterate over dofs on node nn
305 for (; diit != hi_diit; diit++) {
306 // iterate overt dofs in element column
307 for (int ddd = 0; ddd < nb_col; ddd++) {
308 // check consistency, node has to be at crack front
309 if (rowFrontIndices[3 * nn + diit->get()->getDofCoeffIdx()] !=
310 diit->get()->getPetscGlobalDofIdx()) {
311 SETERRQ(
312 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
313 "data inconsistency %d != %d",
314 rowFrontIndices[3 * nn + diit->get()->getDofCoeffIdx()],
315 diit->get()->getPetscGlobalDofIdx());
316 }
317 // dof is not on this partition
318 if (diit->get()->getPetscLocalDofIdx() == -1) {
319 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
320 "data inconsistency");
321 }
322 double g =
323 f_tangent_front_mesh_array[diit->get()
324 ->getPetscLocalDofIdx()] *
325 k(3 * nn + diit->get()->getDofCoeffIdx(), ddd);
326 int lambda_idx = dit->get()->getPetscGlobalDofIdx();
327 CHKERR MatSetValues(getFEMethod()->snes_B, 1, &lambda_idx, 1,
328 &col_indices_ptr[ddd], &g, ADD_VALUES);
329 }
330 }
331 }
332 }
333 CHKERR VecRestoreArray(smootherData.tangentFrontF,
334 &f_tangent_front_mesh_array);
335 }
336
338 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
constexpr double g
unsigned int getFieldBitNumber(std::string field_name) const
auto getRowDofsPtr() const
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
ublas::vector< int > rowFrontIndices
Definition Smoother.hpp:231

◆ aSemble() [2/2]

MoFEMErrorCode Smoother::OpLhsSmoother::aSemble ( int row_side,
int col_side,
EntityType row_type,
EntityType col_type,
EntitiesFieldData::EntData & row_data,
EntitiesFieldData::EntData & col_data )
inlinevirtual

Reimplemented from NonlinearElasticElement::OpLhsPiolaKirchhoff_dx.

Definition at line 233 of file Smoother.hpp.

236 {
238
239 int nb_row = row_data.getIndices().size();
240 int nb_col = col_data.getIndices().size();
241 int *row_indices_ptr = &row_data.getIndices()[0];
242 int *col_indices_ptr = &col_data.getIndices()[0];
243
244 if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
245 rowIndices.resize(nb_row, false);
246 noalias(rowIndices) = row_data.getIndices();
248 row_indices_ptr = &rowIndices[0];
249 }
250 rowFrontIndices.resize(nb_row, false);
251 noalias(rowFrontIndices) = row_data.getIndices();
252 VectorDofs &dofs = row_data.getFieldDofs();
253 VectorDofs::iterator dit = dofs.begin();
254 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
255 if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) !=
257 rowIndices[ii] = -1;
258 } else {
259 rowFrontIndices[ii] = -1;
260 }
261 }
262 }
263
264 CHKERR MatSetValues(getFEMethod()->snes_B, nb_row, row_indices_ptr,
265 nb_col, col_indices_ptr, &k(0, 0), ADD_VALUES);
266
268
269 const auto bit_number_for_crack_area_tangent_constrain =
271 const auto bit_number_for_mesh_position =
272 getFEMethod()->getFieldBitNumber("MESH_NODE_POSITIONS");
273
274 // get tangent vector array
275 double *f_tangent_front_mesh_array;
277 &f_tangent_front_mesh_array);
278
279 auto row_dofs = getFEMethod()->getRowDofsPtr();
280
281 // iterate nodes on tet
282 for (int nn = 0; nn < 4; nn++) {
283
284 // get indices with Lagrange multiplier at node nn
285 auto dit = row_dofs->get<Unique_mi_tag>().lower_bound(
286 DofEntity::getLoFieldEntityUId(
287 bit_number_for_crack_area_tangent_constrain, getConn()[nn]));
288 auto hi_dit = row_dofs->get<Unique_mi_tag>().upper_bound(
289 DofEntity::getHiFieldEntityUId(
290 bit_number_for_crack_area_tangent_constrain, getConn()[nn]));
291
292 // continue if Lagrange are on element
293 if (std::distance(dit, hi_dit) > 0) {
294
295 // get mesh node positions at node nn
296 auto diit = row_dofs->get<Unique_mi_tag>().lower_bound(
297 DofEntity::getLoFieldEntityUId(bit_number_for_mesh_position,
298 getConn()[nn]));
299
300 auto hi_diit = row_dofs->get<Unique_mi_tag>().upper_bound(
301 DofEntity::getHiFieldEntityUId(bit_number_for_mesh_position,
302 getConn()[nn]));
303
304 // iterate over dofs on node nn
305 for (; diit != hi_diit; diit++) {
306 // iterate overt dofs in element column
307 for (int ddd = 0; ddd < nb_col; ddd++) {
308 // check consistency, node has to be at crack front
309 if (rowFrontIndices[3 * nn + diit->get()->getDofCoeffIdx()] !=
310 diit->get()->getPetscGlobalDofIdx()) {
311 SETERRQ(
312 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
313 "data inconsistency %d != %d",
314 rowFrontIndices[3 * nn + diit->get()->getDofCoeffIdx()],
315 diit->get()->getPetscGlobalDofIdx());
316 }
317 // dof is not on this partition
318 if (diit->get()->getPetscLocalDofIdx() == -1) {
319 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
320 "data inconsistency");
321 }
322 double g =
323 f_tangent_front_mesh_array[diit->get()
324 ->getPetscLocalDofIdx()] *
325 k(3 * nn + diit->get()->getDofCoeffIdx(), ddd);
326 int lambda_idx = dit->get()->getPetscGlobalDofIdx();
327 CHKERR MatSetValues(getFEMethod()->snes_B, 1, &lambda_idx, 1,
328 &col_indices_ptr[ddd], &g, ADD_VALUES);
329 }
330 }
331 }
332 }
333 CHKERR VecRestoreArray(smootherData.tangentFrontF,
334 &f_tangent_front_mesh_array);
335 }
336
338 }

Member Data Documentation

◆ fieldCrackAreaTangentConstrains

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

Definition at line 216 of file Smoother.hpp.

◆ rowFrontIndices

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

Definition at line 231 of file Smoother.hpp.

◆ smootherData

SmootherBlockData & Smoother::OpLhsSmoother::smootherData

Definition at line 215 of file Smoother.hpp.


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