v0.14.0
mesh_smoothing.cpp
/** \file mesh_smoothing.cpp
* \brief Improve mesh quality using Volume-length quality measure with barrier
* \example mesh_smoothing.cpp
*
*/
#include <MoFEM.hpp>
using namespace MoFEM;
#include <Smoother.hpp>
#include <SurfaceSlidingConstrains.hpp>
#include <VolumeLengthQuality.hpp>
using namespace MoFEM;
static char help[] = "mesh smoothing\n\n";
PetscBool flg_myfile = PETSC_TRUE;
char mesh_file_name[255];
PetscBool output_vtk = PETSC_TRUE;
double tol = 0.1;
double gamma_factor = 0.8;
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "Mesh cut options", "none");
CHKERR PetscOptionsString("-my_file", "mesh file name", "", "mesh.h5m",
CHKERR PetscOptionsInt("-edges_block_set", "edges side set", "",
CHKERR PetscOptionsInt("-vertex_block_set", "vertex side set", "",
CHKERR PetscOptionsBool("-output_vtk", "if true outout vtk file", "",
output_vtk, &output_vtk, PETSC_NULL);
CHKERR PetscOptionsScalar("-quality_reduction_tol", "",
"Tolerance of quality reduction", tol, &tol,
PETSC_NULL);
CHKERR PetscOptionsScalar("-gamma_factor", "",
"Gamma factor", gamma_factor, &gamma_factor,
PETSC_NULL);
ierr = PetscOptionsEnd();
// Create MoAB database
moab::Core moab_core; // create database
moab::Interface &moab = moab_core; // create interface to database
// Create MoFEM database and link it to MoAB
MoFEM::Core mofem_core(moab); // create database
MoFEM::Interface &m_field = mofem_core; // create interface to database
// Register DM Manager
CHKERR DMRegister_MoFEM("DMMOFEM"); // register MoFEM DM in PETSc
// Get simple interface is simplified version enabling quick and
// easy construction of problem.
Simple *simple_interface;
// Query interface and get pointer to Simple interface
CHKERR m_field.getInterface(simple_interface);
// Build problem with simple interface
{
// Get options for simple interface from command line
CHKERR simple_interface->getOptions();
// Load mesh file to database
CHKERR simple_interface->loadFile();
// Add domain filed "U" in space H1 and Legendre base, Ainsworth recipe is
// used to construct base functions.
CHKERR simple_interface->addDomainField("MESH_NODE_POSITIONS", H1,
// Add Lagrange multiplier field on body boundary
CHKERR simple_interface->addBoundaryField("LAMBDA_SURFACE", H1,
// Set fields order domain and boundary fields.
CHKERR simple_interface->setFieldOrder("MESH_NODE_POSITIONS",
1); // to approximate function
CHKERR simple_interface->setFieldOrder("LAMBDA_SURFACE",
1); // to Lagrange multipliers
simple_interface->getDomainFEName() = "SMOOTHING";
simple_interface->getBoundaryFEName() = "SURFACE_SLIDING";
// Other fields and finite elements added to database directly
{
// Declare approximation fields
CHKERR m_field.add_field("LAMBDA_EDGE", H1, AINSWORTH_LOBATTO_BASE, 2,
MB_TAG_SPARSE, MF_ZERO);
Range edges;
->getEntitiesByDimension(edges_block_set, BLOCKSET, 1, edges,
true);
CHKERR m_field.add_ents_to_field_by_type(edges, MBEDGE,
"LAMBDA_EDGE");
->synchroniseFieldEntities("LAMBDA_EDGE");
CHKERR m_field.set_field_order(0, MBVERTEX, "LAMBDA_EDGE", 1);
CHKERR m_field.add_finite_element("EDGE_SLIDING");
"EDGE_SLIDING", "MESH_NODE_POSITIONS");
"EDGE_SLIDING", "MESH_NODE_POSITIONS");
"EDGE_SLIDING", "MESH_NODE_POSITIONS");
"LAMBDA_EDGE");
"LAMBDA_EDGE");
"LAMBDA_EDGE");
"EDGE_SLIDING");
simple_interface->getOtherFiniteElements().push_back("EDGE_SLIDING");
}
}
CHKERR simple_interface->defineFiniteElements();
CHKERR simple_interface->defineProblem();
CHKERR simple_interface->buildFields();
// Remove vertices form LAMBDA_SURFACE which are on the edges
Range edges;
CHKERR m_field.getInterface<MeshsetsManager>()->getEntitiesByDimension(
edges_block_set, BLOCKSET, 1, edges, true);
Range verts;
CHKERR m_field.get_moab().get_connectivity(edges, verts, true);
CHKERR m_field.remove_ents_from_field("LAMBDA_SURFACE", verts);
}
CHKERR simple_interface->buildFiniteElements();
CHKERR simple_interface->buildProblem();
}
struct ElementsAndOperators {
Vec minQualityVec;
double *minQualityPtr;
ElementsAndOperators(MoFEM::Interface &m_field) : mField(m_field) {
ierr = VecCreateMPI(PETSC_COMM_WORLD, 1, m_field.get_comm_size(),
&minQualityVec);
CHKERRABORT(PETSC_COMM_WORLD, ierr);
ierr = VecGetArray(minQualityVec, &minQualityPtr);
CHKERRABORT(PETSC_COMM_WORLD, ierr);
}
virtual ~ElementsAndOperators() {
ierr = VecRestoreArray(minQualityVec, &minQualityPtr);
CHKERRABORT(PETSC_COMM_WORLD, ierr);
ierr = VecDestroy(&minQualityVec);
CHKERRABORT(PETSC_COMM_WORLD, ierr);
}
double getMinQuality() const { return *minQualityPtr; }
enum Tags {
SURFACE_CONSTRAINS_TAG,
EDGE_CONSTRAINS_TAG
};
boost::shared_ptr<Smoother> smootherFe;
boost::shared_ptr<Smoother::MyVolumeFE>
feSmootherRhs; ///< Integrate smoothing operators
boost::shared_ptr<Smoother::MyVolumeFE>
feSmootherLhs; ///< Integrate smoothing operators
boost::shared_ptr<VolumeLengthQuality<double> > volumeLengthDouble;
boost::shared_ptr<VolumeLengthQuality<adouble> > volumeLengthAdouble;
boost::shared_ptr<SurfaceSlidingConstrains> surfaceConstrain;
boost::shared_ptr<SurfaceSlidingConstrains::DriverElementOrientation>
skinOrientation;
boost::shared_ptr<EdgeSlidingConstrains> edgeConstrain;
boost::shared_ptr<DirichletFixFieldAtEntitiesBc> fixMaterialEnts;
boost::shared_ptr<MoFEM::VolumeElementForcesAndSourcesCore> minQualityFe;
double *minQualityPtr;
MinQualityOp(double *min_quality_ptr)
"MESH_NODE_POSITIONS", UserDataOperator::OPROW),
minQualityPtr(min_quality_ptr) {}
MoFEMErrorCode doWork(int side, EntityType type,
if (type != MBVERTEX)
double q = Tools::volumeLengthQuality(&*data.getFieldData().begin());
*minQualityPtr = fmin(*minQualityPtr, q);
}
};
MoFEMErrorCode createSmoothingFE() {
smootherFe = boost::shared_ptr<Smoother>(new Smoother(mField));
volumeLengthAdouble = boost::shared_ptr<VolumeLengthQuality<adouble> >(
volumeLengthDouble = boost::shared_ptr<VolumeLengthQuality<double> >(
Range tets;
CHKERR mField.get_moab().get_entities_by_type(0, MBTET, tets);
smootherFe->setOfBlocks[0].tEts.merge(tets);
smootherFe->setOfBlocks[0].materialDoublePtr = volumeLengthDouble;
smootherFe->setOfBlocks[0].materialAdoublePtr = volumeLengthAdouble;
// set element data
smootherFe->commonData.spatialPositions = "MESH_NODE_POSITIONS";
smootherFe->commonData.meshPositions = "NONE";
smootherFe->feRhs.meshPositionsFieldName = "NONE";
smootherFe->feLhs.meshPositionsFieldName = "NONE";
smootherFe->feRhs.addToRule = 0;
smootherFe->feLhs.addToRule = 0;
feSmootherRhs = smootherFe->feRhsPtr;
feSmootherLhs = smootherFe->feLhsPtr;
// Smoother right hand side
feSmootherRhs->getOpPtrVector().push_back(
"MESH_NODE_POSITIONS", smootherFe->commonData));
feSmootherRhs->getOpPtrVector().push_back(
"MESH_NODE_POSITIONS", smootherFe->setOfBlocks.at(0),
smootherFe->commonData, SMOOTHING_TAG, false));
feSmootherRhs->getOpPtrVector().push_back(new Smoother::OpRhsSmoother(
"MESH_NODE_POSITIONS", smootherFe->setOfBlocks[0],
smootherFe->commonData, smootherFe->smootherData));
// Smoother left hand side
feSmootherLhs->getOpPtrVector().push_back(
"MESH_NODE_POSITIONS", smootherFe->commonData));
feSmootherLhs->getOpPtrVector().push_back(
"MESH_NODE_POSITIONS", smootherFe->setOfBlocks.at(0),
smootherFe->commonData, SMOOTHING_TAG, true));
feSmootherLhs->getOpPtrVector().push_back(new Smoother::OpLhsSmoother(
"MESH_NODE_POSITIONS", "MESH_NODE_POSITIONS",
smootherFe->setOfBlocks.at(0), smootherFe->commonData,
smootherFe->smootherData, "LAMBDA_CRACKFRONT_AREA_TANGENT"));
minQualityFe =
boost::shared_ptr<MoFEM::VolumeElementForcesAndSourcesCore>(
minQualityFe->getOpPtrVector().push_back(
new MinQualityOp(minQualityPtr));
Range fixed_vertex;
CHKERR mField.getInterface<MeshsetsManager>()->getEntitiesByDimension(
vertex_block_set, BLOCKSET, 0, fixed_vertex, true);
}
fixMaterialEnts = boost::shared_ptr<DirichletFixFieldAtEntitiesBc>(
new DirichletFixFieldAtEntitiesBc(mField, "MESH_NODE_POSITIONS",
fixed_vertex));
fixMaterialEnts->fieldNames.push_back("LAMBDA_SURFACE");
fixMaterialEnts->fieldNames.push_back("LAMBDA_EDGE");
}
MoFEMErrorCode createConstrians() {
skinOrientation = boost::shared_ptr<
surfaceConstrain = boost::shared_ptr<SurfaceSlidingConstrains>(
skinOrientation,
new SurfaceSlidingConstrains(mField, *skinOrientation));
surfaceConstrain->setOperators(SURFACE_CONSTRAINS_TAG, "LAMBDA_SURFACE",
"MESH_NODE_POSITIONS");
Range edges;
->getEntitiesByDimension(edges_block_set, BLOCKSET, 1, edges,
true);
Range tets;
CHKERR mField.get_moab().get_entities_by_type(0, MBTET, tets);
Skinner skin(&mField.get_moab());
Range skin_faces; // skin faces from 3d ents
CHKERR skin.find_skin(0, tets, false, skin_faces);
edgeConstrain = boost::shared_ptr<EdgeSlidingConstrains>(
new EdgeSlidingConstrains(mField));
CHKERR edgeConstrain->setOperators(EDGE_CONSTRAINS_TAG, edges,
skin_faces, "LAMBDA_EDGE",
"MESH_NODE_POSITIONS");
// CHKERR EdgeSlidingConstrains::CalculateEdgeBase::saveEdges(
// mField.get_moab(), "out_edges.vtk", edges);
}
}
MoFEMErrorCode addFEtoDM(DM dm) {
boost::shared_ptr<ForcesAndSourcesCore> null;
CHKERR DMMoFEMSNESSetFunction(dm, DM_NO_ELEMENT, null, fixMaterialEnts,
null);
CHKERR DMMoFEMSNESSetFunction(dm, "SMOOTHING", feSmootherRhs, null,
null);
CHKERR DMMoFEMSNESSetFunction(dm, "SURFACE_SLIDING",
surfaceConstrain->feRhsPtr, null, null);
CHKERR DMMoFEMSNESSetFunction(dm, "EDGE_SLIDING",
edgeConstrain->feRhsPtr, null, null);
fixMaterialEnts);
CHKERR DMMoFEMSNESSetJacobian(dm, DM_NO_ELEMENT, null, fixMaterialEnts,
null);
CHKERR DMMoFEMSNESSetJacobian(dm, "SMOOTHING", feSmootherLhs, null,
null);
CHKERR DMMoFEMSNESSetJacobian(dm, "SURFACE_SLIDING",
surfaceConstrain->feLhsPtr, null, null);
CHKERR DMMoFEMSNESSetJacobian(dm, "EDGE_SLIDING",
edgeConstrain->feLhsPtr, null, null);
fixMaterialEnts);
// MoFEM::SnesCtx *snes_ctx;
// DMMoFEMGetSnesCtx(dm,&snes_ctx);
// snes_ctx->vErify = true;
}
MoFEMErrorCode calcuteMinQuality(DM dm) {
*minQualityPtr = 1;
CHKERR DMoFEMLoopFiniteElements(dm, "SMOOTHING", minQualityFe.get());
CHKERR VecMin(minQualityVec, PETSC_NULL, minQualityPtr);
}
};
ElementsAndOperators elements_and_operators(m_field);
CHKERR elements_and_operators.createSmoothingFE();
CHKERR elements_and_operators.createConstrians();
DM dm;
CHKERR simple_interface->getDM(&dm);
CHKERR elements_and_operators.addFEtoDM(dm);
struct Solve {
MoFEMErrorCode operator()(DM dm) const {
// Create the right hand side vector and vector of unknowns
Vec F, D;
CHKERR DMCreateGlobalVector(dm, &F);
// Create unknown vector by creating duplicate copy of F vector. only
// structure is duplicated no values.
CHKERR VecDuplicate(F, &D);
CHKERR zeroLambdaFields(dm);
CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_FORWARD);
// Create solver and link it to DM
SNES solver;
CHKERR SNESCreate(PETSC_COMM_WORLD, &solver);
CHKERR SNESSetFromOptions(solver);
CHKERR SNESSetDM(solver, dm);
// Set-up solver, is type of solver and pre-conditioners
CHKERR SNESSetUp(solver);
// At solution process, KSP solver using DM creates matrices, Calculate
// values of the left hand side and the right hand side vector. then
// solves system of equations. Results are stored in vector D.
CHKERR SNESSolve(solver, F, D);
// Scatter solution on the mesh. Stores unknown vector on field on the
// mesh.
CHKERR DMoFEMMeshToGlobalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
// Clean data. Solver and vector are not needed any more.
CHKERR SNESDestroy(&solver);
CHKERR VecDestroy(&D);
CHKERR VecDestroy(&F);
}
MoFEMErrorCode setCoordsFromField(DM dm) const {
MoFEM::Interface *m_field_ptr;
CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
"MESH_NODE_POSITIONS", it)) {
if (it->get()->getEntType() != MBVERTEX)
continue;
VectorDouble3 coords(3);
for(int dd = 0;dd!=3;++dd)
coords[dd] = it->get()->getEntFieldData()[dd];
EntityHandle ent = it->get()->getEnt();
CHKERR m_field_ptr->get_moab().set_coords(&ent, 1, &*coords.begin());
}
}
MoFEMErrorCode setFieldFromCoords(DM dm) const {
MoFEM::Interface *m_field_ptr;
CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
"MESH_NODE_POSITIONS", it)) {
if (it->get()->getEntType() != MBVERTEX)
continue;
EntityHandle ent = it->get()->getEnt();
VectorDouble3 coords(3);
CHKERR m_field_ptr->get_moab().get_coords(&ent, 1, &*coords.begin());
for(int dd = 0;dd!=3;++dd)
it->get()->getEntFieldData()[dd] = coords[dd];
}
}
private:
MoFEMErrorCode zeroLambdaFields(DM dm) const {
MoFEM::Interface *m_field_ptr;
CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
CHKERR m_field_ptr->getInterface<FieldBlas>()->setField(
0, MBVERTEX, "LAMBDA_SURFACE");
}
};
Solve solve;
CHKERR solve.setFieldFromCoords(dm);
CHKERR elements_and_operators.calcuteMinQuality(dm);
double min_quality = elements_and_operators.getMinQuality();
PetscPrintf(PETSC_COMM_WORLD, "Min quality = %4.3f\n", min_quality);
double gamma = min_quality > 0 ? gamma_factor * min_quality
: min_quality / gamma_factor;
elements_and_operators.volumeLengthDouble->gAmma = gamma;
elements_and_operators.volumeLengthAdouble->gAmma = gamma;
double min_quality_p, eps;
do {
min_quality_p = min_quality;
CHKERR solve.setCoordsFromField(dm);
CHKERR elements_and_operators.calcuteMinQuality(dm);
min_quality = elements_and_operators.getMinQuality();
eps = (min_quality - min_quality_p) / min_quality;
PetscPrintf(PETSC_COMM_WORLD, "Min quality = %4.3f eps = %4.3f\n",
min_quality, eps);
double gamma = min_quality > 0 ? gamma_factor * min_quality
: min_quality / gamma_factor;
elements_and_operators.volumeLengthDouble->gAmma = gamma;
elements_and_operators.volumeLengthAdouble->gAmma = gamma;
} while (eps > tol);
// if (m_field.getInterface<MeshsetsManager>()->checkMeshset(edges_block_set,
// BLOCKSET)) {
// Range edges;
// CHKERR m_field.getInterface<MeshsetsManager>()->getEntitiesByDimension(
// edges_block_set, BLOCKSET, 1, edges, true);
// Range tets;
// CHKERR moab.get_entities_by_type(0,MBTET,tets);
// Skinner skin(&moab);
// Range skin_faces; // skin faces from 3d ents
// CHKERR skin.find_skin(0, tets, false, skin_faces);
// CHKERR EdgeSlidingConstrains::CalculateEdgeBase::setTags(moab, edges,
// skin_faces);
// CHKERR EdgeSlidingConstrains::CalculateEdgeBase::saveEdges(
// moab, "out_edges.vtk", edges);
// }
CHKERR m_field.getInterface<BitRefManager>()->writeBitLevelByType(
BitRefLevel().set(0), BitRefLevel().set(), MBTET, "out.vtk", "VTK",
"");
}
}
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
Smoother::OpLhsSmoother
Definition: Smoother.hpp:213
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::Types::VectorDouble3
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:92
EntityHandle
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
MoFEM::Simple::buildProblem
MoFEMErrorCode buildProblem()
Build problem.
Definition: Simple.cpp:669
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
FractureMechanics::SMOOTHING_TAG
@ SMOOTHING_TAG
Definition: CrackPropagation.hpp:45
MoFEM::CoreInterface::remove_ents_from_field
virtual MoFEMErrorCode remove_ents_from_field(const std::string name, const EntityHandle meshset, const EntityType type, int verb=DEFAULT_VERBOSITY)=0
remove entities from field
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1254
MoFEM::Simple::loadFile
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
flg_myfile
PetscBool flg_myfile
Definition: mesh_smoothing.cpp:21
MoFEM.hpp
MoFEM::DMoFEMMeshToLocalVector
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:523
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
gamma_factor
double gamma_factor
Definition: mesh_smoothing.cpp:27
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEM::CoreInterface::add_ents_to_field_by_type
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
DirichletFixFieldAtEntitiesBc
Fix dofs on entities.
Definition: DirichletBC.hpp:258
MoFEM::Simple::getOptions
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
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::Simple::getDM
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:747
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
vertex_block_set
int vertex_block_set
Definition: mesh_smoothing.cpp:24
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
output_vtk
PetscBool output_vtk
Definition: mesh_smoothing.cpp:25
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
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
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
convert.type
type
Definition: convert.py:64
main
int main(int argc, char *argv[])
Definition: mesh_smoothing.cpp:29
DM_NO_ELEMENT
#define DM_NO_ELEMENT
Definition: DMMoFEM.hpp:10
MoFEM::DMoFEMMeshToGlobalVector
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Definition: DMMoFEM.cpp:535
_IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_
#define _IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(MFIELD, NAME, IT)
loop over all dofs from a moFEM field and particular field
Definition: Interface.hpp:1842
MoFEM::Simple::addDomainField
MoFEMErrorCode addDomainField(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
Definition: Simple.cpp:264
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::CoreInterface::get_comm_size
virtual int get_comm_size() const =0
MoFEM::Simple::getDomainFEName
const std::string getDomainFEName() const
Get the Domain FE Name.
Definition: Simple.hpp:368
MoFEM::solve
static PetscErrorCode solve(Mat mat, Vec x, Vec y)
Definition: Schur.cpp:1223
MoFEM::DMMoFEMSNESSetJacobian
PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES Jacobian evaluation function
Definition: DMMoFEM.cpp:759
AINSWORTH_LOBATTO_BASE
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:62
VolumeLengthQuality
Volume Length Quality.
Definition: VolumeLengthQuality.hpp:32
SurfaceSlidingConstrains
Shape preserving constrains, i.e. nodes sliding on body surface.
Definition: SurfaceSlidingConstrains.hpp:323
MoFEM::ForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: ForcesAndSourcesCore.hpp:482
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_field)=0
set finite element field data
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
help
static char help[]
Definition: mesh_smoothing.cpp:19
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
MoFEM::Simple::setFieldOrder
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:545
MoFEM::Simple::buildFiniteElements
MoFEMErrorCode buildFiniteElements()
Build finite elements.
Definition: Simple.cpp:643
Range
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
Smoother
Definition: Smoother.hpp:10
MF_ZERO
@ MF_ZERO
Definition: definitions.h:111
edges_block_set
int edges_block_set
Definition: mesh_smoothing.cpp:23
NonlinearElasticElement::OpGetCommonDataAtGaussPts
Definition: NonLinearElasticElement.hpp:362
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::CommInterface
Managing BitRefLevels.
Definition: CommInterface.hpp:21
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
MoFEM::DMoFEMGetInterfacePtr
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition: DMMoFEM.cpp:414
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::Simple::getBoundaryFEName
const std::string getBoundaryFEName() const
Get the Boundary FE Name.
Definition: Simple.hpp:375
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
MoFEM::Simple::defineFiniteElements
MoFEMErrorCode defineFiniteElements()
Define finite elements.
Definition: Simple.cpp:452
MoFEM::Simple::buildFields
MoFEMErrorCode buildFields()
Build fields.
Definition: Simple.cpp:554
MoFEM::Tools::volumeLengthQuality
static double volumeLengthQuality(const double *coords)
Calculate tetrahedron volume length quality.
Definition: Tools.cpp:15
MoFEM::Simple::addBoundaryField
MoFEMErrorCode addBoundaryField(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on boundary.
Definition: Simple.cpp:354
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
Smoother::OpRhsSmoother
Definition: Smoother.hpp:159
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::Simple::defineProblem
MoFEMErrorCode defineProblem(const PetscBool is_partitioned=PETSC_TRUE)
define problem
Definition: Simple.cpp:524
BARRIER_AND_QUALITY
@ BARRIER_AND_QUALITY
Definition: VolumeLengthQuality.hpp:17
MoFEM::CoreInterface::set_field_order
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
MoFEM::MeshsetsManager::checkMeshset
bool checkMeshset(const int ms_id, const CubitBCType cubit_bc_type) const
check for CUBIT Id and CUBIT type
Definition: MeshsetsManager.cpp:360
MoFEM::Simple::getOtherFiniteElements
std::vector< std::string > & getOtherFiniteElements()
Get the Other Finite Elements.
Definition: Simple.hpp:427
MoFEM::DMoFEMLoopFiniteElements
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:586
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::FieldBlas
Basic algebra on fields.
Definition: FieldBlas.hpp:21
Smoother::OpJacobianSmoother
Definition: Smoother.hpp:131
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496
MoFEM::CoreInterface::add_field
virtual MoFEMErrorCode add_field(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
SurfaceSlidingConstrains::DriverElementOrientation
Class implemented by user to detect face orientation.
Definition: SurfaceSlidingConstrains.hpp:334
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
tol
double tol
Definition: mesh_smoothing.cpp:26
F
@ F
Definition: free_surface.cpp:394
EdgeSlidingConstrains
Definition: SurfaceSlidingConstrains.hpp:639
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567
MoFEM::DMMoFEMSNESSetFunction
PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES residual evaluation function
Definition: DMMoFEM.cpp:718