static char help[] =
"mesh smoothing\n\n";
int main(
int argc,
char *argv[]) {
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",
"",
CHKERR PetscOptionsScalar(
"-quality_reduction_tol",
"",
"Tolerance of quality reduction",
tol, &
tol,
PETSC_NULL);
CHKERR PetscOptionsScalar(
"-gamma_factor",
"",
PETSC_NULL);
ierr = PetscOptionsEnd();
moab::Core moab_core;
moab::Interface &moab = moab_core;
CHKERR DMRegister_MoFEM(
"DMMOFEM");
{
1);
1);
{
true);
"LAMBDA_EDGE");
->synchroniseFieldEntities("LAMBDA_EDGE");
"EDGE_SLIDING", "MESH_NODE_POSITIONS");
"EDGE_SLIDING", "MESH_NODE_POSITIONS");
"EDGE_SLIDING", "MESH_NODE_POSITIONS");
"LAMBDA_EDGE");
"LAMBDA_EDGE");
"LAMBDA_EDGE");
"EDGE_SLIDING");
}
}
}
}
struct ElementsAndOperators {
Vec minQualityVec;
double *minQualityPtr;
&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 {
SMOOTHING_TAG = 1,
SURFACE_CONSTRAINS_TAG,
EDGE_CONSTRAINS_TAG
};
boost::shared_ptr<Smoother> smootherFe;
boost::shared_ptr<Smoother::MyVolumeFE>
feSmootherRhs;
boost::shared_ptr<Smoother::MyVolumeFE>
feSmootherLhs;
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) {}
if (type != MBVERTEX)
double q = Tools::volumeLengthQuality(&*data.
getFieldData().begin());
*minQualityPtr = fmin(*minQualityPtr, q);
}
};
smootherFe = boost::shared_ptr<Smoother>(
new Smoother(mField));
volumeLengthAdouble = boost::shared_ptr<VolumeLengthQuality<adouble> >(
volumeLengthDouble = boost::shared_ptr<VolumeLengthQuality<double> >(
smootherFe->setOfBlocks[0].tEts.merge(tets);
smootherFe->setOfBlocks[0].materialDoublePtr = volumeLengthDouble;
smootherFe->setOfBlocks[0].materialAdoublePtr = volumeLengthAdouble;
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;
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));
"MESH_NODE_POSITIONS", smootherFe->setOfBlocks[0],
smootherFe->commonData, smootherFe->smootherData));
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));
"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));
}
fixMaterialEnts = boost::shared_ptr<DirichletFixFieldAtEntitiesBc>(
fixed_vertex));
fixMaterialEnts->fieldNames.push_back("LAMBDA_SURFACE");
fixMaterialEnts->fieldNames.push_back("LAMBDA_EDGE");
}
skinOrientation = boost::shared_ptr<
surfaceConstrain = boost::shared_ptr<SurfaceSlidingConstrains>(
skinOrientation,
surfaceConstrain->setOperators(SURFACE_CONSTRAINS_TAG, "LAMBDA_SURFACE",
"MESH_NODE_POSITIONS");
true);
CHKERR skin.find_skin(0, tets,
false, skin_faces);
edgeConstrain = boost::shared_ptr<EdgeSlidingConstrains>(
CHKERR edgeConstrain->setOperators(EDGE_CONSTRAINS_TAG, edges,
skin_faces, "LAMBDA_EDGE",
"MESH_NODE_POSITIONS");
}
}
boost::shared_ptr<ForcesAndSourcesCore> null;
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);
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);
}
*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 elements_and_operators.addFEtoDM(dm);
struct Solve {
CHKERR DMCreateGlobalVector(dm, &
F);
CHKERR DMoFEMMeshToLocalVector(dm,
D, INSERT_VALUES, SCATTER_FORWARD);
SNES solver;
CHKERR SNESCreate(PETSC_COMM_WORLD, &solver);
CHKERR SNESSetFromOptions(solver);
CHKERR DMoFEMMeshToGlobalVector(dm,
D, INSERT_VALUES, SCATTER_REVERSE);
}
CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
"MESH_NODE_POSITIONS", it)) {
if (it->get()->getEntType() != MBVERTEX)
continue;
for(int dd = 0;dd!=3;++dd)
coords[dd] = it->get()->getEntFieldData()[dd];
}
}
CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
"MESH_NODE_POSITIONS", it)) {
if (it->get()->getEntType() != MBVERTEX)
continue;
for(int dd = 0;dd!=3;++dd)
it->get()->getEntFieldData()[dd] = coords[dd];
}
}
private:
CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
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);
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",
elements_and_operators.volumeLengthDouble->gAmma = gamma;
elements_and_operators.volumeLengthAdouble->gAmma = gamma;
"");
}
}
Implementing surface sliding constrains.
Implementation of Volume-Length-Quality measure with barrier.
static PetscErrorCode ierr
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define _IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(MFIELD, NAME, IT)
loop over all dofs from a moFEM field and particular field
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
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
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
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
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
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.
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.
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
bool checkMeshset(const int ms_id, const CubitBCType cubit_bc_type) const
check for CUBIT Id and CUBIT type
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VectorBoundedArray< double, 3 > VectorDouble3
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
implementation of Data Operators for Forces and Sources
virtual int get_comm_size() const =0
virtual moab::Interface & get_moab()=0
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.
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorDouble & getFieldData() const
get dofs values
Interface for managing meshsets containing materials and boundary conditions.
Simple interface for fast problem set-up.
MoFEMErrorCode buildProblem()
Build problem.
const std::string getBoundaryFEName() const
Get the Boundary FE Name.
MoFEMErrorCode defineFiniteElements()
Define finite elements.
MoFEMErrorCode buildFiniteElements()
Build finite elements.
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.
std::vector< std::string > & getOtherFiniteElements()
Get the Other Finite Elements.
MoFEMErrorCode getOptions()
get options
MoFEMErrorCode getDM(DM *dm)
Get DM.
MoFEMErrorCode buildFields()
Build fields.
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name)
Load mesh file.
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
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.
MoFEMErrorCode defineProblem(const PetscBool is_partitioned=PETSC_TRUE)
define problem
const std::string getDomainFEName() const
Get the Domain FE Name.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Volume finite element base.
Class implemented by user to detect face orientation.
Shape preserving constrains, i.e. nodes sliding on body surface.