|
| v0.14.0
|
Go to the documentation of this file.
7 #ifndef __CUTMESHINTERFACE_HPP__
8 #define __CUTMESHINTERFACE_HPP__
39 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD,
"cut_",
"MOFEM Cut mesh options",
42 CHKERR PetscOptionsInt(
"-linesearch_steps",
43 "number of bisection steps which line search do to "
44 "find optimal merged nodes position",
47 CHKERR PetscOptionsInt(
"-max_merging_cycles",
48 "number of maximal merging cycles",
"",
51 CHKERR PetscOptionsScalar(
"-project_entities_quality_trashold",
52 "project entities quality trashold",
"",
56 ierr = PetscOptionsEnd();
76 double *shift = NULL,
double *origin = NULL,
78 const std::string save_mesh =
"");
115 const double rel_tol,
116 const double abs_tol, Tag
th =
nullptr,
117 const bool debug =
false);
120 const Range fixed_edges,
121 const double rel_tol,
const double abs_tol,
122 Tag
th =
nullptr,
const bool debug =
false);
145 const double tol_geometry,
const double tol_cut_close,
146 Range *fixed_edges = NULL,
Range *corner_nodes = NULL,
147 const bool update_meshsets =
false,
148 const bool debug =
false);
163 const double tol_trim_close,
164 Range *fixed_edges = NULL,
Range *corner_nodes = NULL,
165 const bool update_meshsets =
false,
166 const bool debug =
false);
183 cutAndTrim(
int &first_bit, Tag
th,
const double tol_geometry,
184 const double tol_cut_close,
const double tol_trim_close,
185 Range *fixed_edges = NULL,
Range *corner_nodes = NULL,
186 const bool update_meshsets =
false,
const bool debug =
false);
207 Tag
th,
const double tol_geometry,
208 const double tol_cut_close,
209 const double tol_trim_close,
211 const bool update_meshsets =
false,
212 const bool debug =
false);
231 const bool debug =
false);
244 const bool debug =
false);
258 const bool remove_adj_prims_edges,
260 const std::string edges_file_name =
string());
271 const bool debug =
false);
289 const int front_levels,
290 Range *fixed_edges =
nullptr,
int verb =
QUIET,
291 const bool debug =
false);
302 const bool debug =
false);
316 const double geometry_tol = 0,
317 const double close_tol = 0,
318 const int verb =
QUIET,
319 const bool debug =
false);
333 const bool debug =
false);
353 const double tol = 1e-4,
354 const bool debug =
false);
362 const bool debug =
false);
379 const double tol = 1e-4, Tag
th = NULL,
380 const bool debug =
false);
406 const Range &ents, Tag
th = NULL);
430 const Range &corner_nodes,
Range &merged_nodes,
432 const bool update_meshsets =
false,
434 const bool debug =
false);
447 const Range &corner_nodes, Tag
th = NULL,
448 const bool update_meshsets =
false,
const bool debug =
false);
504 CHKERR moab.create_meshset(MESHSET_SET, meshset);
506 CHKERR moab.write_file(name.c_str(),
"VTK",
"", &meshset, 1);
521 typedef multi_index_container<
526 member<LengthMapData, double, &LengthMapData::lEngth>>,
529 member<LengthMapData, EntityHandle, &LengthMapData::eDge>>,
532 member<LengthMapData, double, &LengthMapData::qUality>>
581 #endif //__CUTMESHINTERFACE_HPP__
const Range & projectZeroDistanceEnts() const
MoFEMErrorCode trimOnly(const BitRefLevel trim_bit, Tag th, const double tol_trim_close, Range *fixed_edges=NULL, Range *corner_nodes=NULL, const bool update_meshsets=false, const bool debug=false)
Trim mesh only.
const Range & getNewCutSurfaces() const
VectorBoundedArray< double, 3 > VectorDouble3
MoFEMErrorCode setFront(const Range surface)
MoFEMErrorCode mergeSurface(const Range surface)
merge surface entities
MoFEMErrorCode findLevelSetVolumes(Tag th, Range &vol_edges, const bool remove_adj_prims_edges, int verb=QUIET, const bool debug=false, const std::string edges_file_name=string())
Find level set volumes.
MoFEMErrorCode setVolume(const Range volume)
set volume entities
MoFEMErrorCode setCoords(Tag th, const BitRefLevel bit=BitRefLevel(), const BitRefLevel mask=BitRefLevel().set())
set coords from tag
CutMeshInterface(const MoFEM::Core &core)
MoFEMErrorCode makeFront(const bool debug=false)
Create front from the surface.
MoFEMErrorCode setTagData(Tag th, const BitRefLevel bit=BitRefLevel())
set coords to tag
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
const Range & getFront() const
MoFEMErrorCode copySurface(const Range surface, Tag th=NULL, double *shift=NULL, double *origin=NULL, double *transform=NULL, const std::string save_mesh="")
copy surface entities
const Range & getNewTrimVertices() const
double maxLength
Maximal edge length.
const Range & getNewCutVertices() const
map< EntityHandle, TreeData > verticesOnTrimEdges
MoFEMErrorCode trimSurface(Range *fixed_edge, Range *corner_nodes, const double tol=1e-4, Tag th=NULL, const bool debug=false)
Trim surface from faces beyond front.
MoFEMErrorCode getOptions()
Get options from command line.
const Range & getCutFrontVolumes() const
MoFEMErrorCode cutEdgesInMiddle(const BitRefLevel bit, Range &cut_vols, Range &cut_surf, Range &cut_verts, const bool debug=false)
cut edges
DeprecatedCoreInterface Interface
SaveData(moab::Interface &moab)
MoFEMErrorCode moveMidNodesOnCutEdges(Tag th=NULL)
projecting of mid edge nodes on new mesh on surface
LengthMapData(const double l, double q, const EntityHandle e)
MoFEMErrorCode clearMap()
MoFEMErrorCode mergeVolumes(const Range volume)
merge volume entities
#define CHKERR
Inline error check.
MoFEMErrorCode cutTrimAndMerge(int &first_bit, const int fraction_level, Tag th, const double tol_geometry, const double tol_cut_close, const double tol_trim_close, Range &fixed_edges, Range &corner_nodes, const bool update_meshsets=false, const bool debug=false)
Cut, trim and merge.
MoFEMErrorCode findEdgesToTrim(Range *fixed_edges, Range *corner_nodes, Tag th=NULL, const double tol=1e-4, const bool debug=false)
Find edges to trimEdges.
MoFEMErrorCode operator()(const std::string name, const Range &ents)
MoFEMErrorCode mergeBadEdges(const int fraction_level, const Range &tets, const Range &surface, const Range &fixed_edges, const Range &corner_nodes, Range &merged_nodes, Range &out_tets, Range &new_surf, Tag th, const bool update_meshsets=false, const BitRefLevel *bit_ptr=NULL, const bool debug=false)
Merge edges.
MoFEMErrorCode saveTrimEdges()
implementation of Data Operators for Forces and Sources
const Range & getTrimEdges() const
const Range & getMergedSurfaces() const
boost::shared_ptr< OrientedBoxTreeTool > treeSurfPtr
MoFEMErrorCode removePathologicalFrontTris(const BitRefLevel split_bit, Range &ents)
Remove pathological elements on surface internal front.
map< EntityHandle, TreeData > verticesOnCutEdges
MoFEMErrorCode setConstrainSurface(const Range surf)
Set the constrain surface object.
const Range & getMergedVolumes() const
MoFEMErrorCode snapSurfaceToEdges(const Range surface_edges, const Range fixed_edges, const double rel_tol, const double abs_tol, Tag th=nullptr, const bool debug=false)
MoFEMErrorCode setSurface(const Range surface)
set surface entities
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
MoFEMErrorCode findEdgesToCut(Range vol, int verb=QUIET, const bool debug=false)
find edges to cut
MoFEMErrorCode getSubInterfaceOptions()
MoFEMErrorCode createFrontLevelSets(Range vol, Tag th=nullptr, int verb=QUIET, const bool debug=false)
Calculate distance from mesh nodes to surface front.
MoFEMErrorCode buildTree()
build tree
const Range & getNewCutVolumes() const
base class for all interface classes
MoFEMErrorCode saveCutEdges(const std::string prefix="")
double projectEntitiesQualityTrashold
MoFEMErrorCode trimEdgesInTheMiddle(const BitRefLevel bit, const bool debug=false)
trim edges
MoFEMErrorCode snapSurfaceSkinToEdges(const Range fixed_edges, const double rel_tol, const double abs_tol, Tag th=nullptr, const bool debug=false)
const Range & getVolume() const
MoFEMErrorCode refineMesh(const int init_bit_level, const int surf_levels, const int front_levels, Range *fixed_edges=nullptr, int verb=QUIET, const bool debug=false)
Refine and set level sets.
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
const Range & getCutEdges() const
MoFEMErrorCode splitSides(const BitRefLevel split_bit, const BitRefLevel bit, const Range &ents, Tag th=NULL)
split sides
multi_index_container< LengthMapData, indexed_by< ordered_non_unique< member< LengthMapData, double, &LengthMapData::lEngth > >, hashed_unique< member< LengthMapData, EntityHandle, &LengthMapData::eDge > >, ordered_non_unique< member< LengthMapData, double, &LengthMapData::qUality > > > > LengthMapData_multi_index
Tensor2_Expr< transform_Tensor2< A, B, T, Dim0, Dim1, i, j >, T, Dim0, Dim1, i, j > transform(const Tensor2_Expr< A, T, Dim0, Dim1, i, j > &a, B function)
map< EntityHandle, TreeData > edgesToTrim
double aveLength
Average edge length.
const Range & getNewTrimSurfaces() const
MoFEMErrorCode createSurfaceLevelSets(int verb=QUIET, const bool debug=false)
Calculate distance from mesh nodes to cut surface.
map< EntityHandle, TreeData > edgesToCut
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
const Range & getNewTrimVolumes() const
MoFEMErrorCode cutAndTrim(int &first_bit, Tag th, const double tol_geometry, const double tol_cut_close, const double tol_trim_close, Range *fixed_edges=NULL, Range *corner_nodes=NULL, const bool update_meshsets=false, const bool debug=false)
Cut and trim.
boost::shared_ptr< OrientedBoxTreeTool > & getTreeSurfPtr()
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
const Range & getCutSurfaceVolumes() const
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
FTensor::Index< 'l', 3 > l
const Range & getSurface() const
MoFEMErrorCode moveMidNodesOnTrimmedEdges(Tag th=NULL)
move trimmed edges mid nodes
MoFEMErrorCode cutOnly(Range vol, const BitRefLevel cut_bit, Tag th, const double tol_geometry, const double tol_cut_close, Range *fixed_edges=NULL, Range *corner_nodes=NULL, const bool update_meshsets=false, const bool debug=false)
Cut mesh only.