|
| v0.14.0
|
Go to the documentation of this file.
11 namespace bio = boost::iostreams;
13 using bio::tee_device;
15 using namespace MoFEM;
17 static char help[] =
"...\n\n";
19 int main(
int argc,
char *argv[]) {
28 PetscBool flg = PETSC_TRUE;
30 #if PETSC_VERSION_GE(3, 6, 4)
37 if (flg != PETSC_TRUE) {
38 SETERRQ(PETSC_COMM_SELF, 1,
"*** ERROR -my_file (MESH FILE NEEDED)");
54 std::vector<BitRefLevel> bit_levels;
61 CHKERR PetscPrintf(PETSC_COMM_WORLD,
"Insert Interface %d\n",
67 CHKERR moab.create_meshset(MESHSET_SET, ref_level_meshset);
70 ->getEntitiesByTypeAndRefLevel(bit_levels.back(),
74 ->getEntitiesByTypeAndRefLevel(bit_levels.back(),
78 CHKERR moab.get_entities_by_handle(ref_level_meshset, ref_level_tets,
82 CHKERR interface->
getSides(cubit_meshset, bit_levels.back(),
true, 0);
87 cubit_meshset,
true,
true, 0);
89 CHKERR moab.delete_entities(&ref_level_meshset, 1);
96 ->updateMeshsetByEntitiesChildren(cubit_meshset, bit_levels.back(),
97 cubit_meshset, MBMAXTYPE,
true);
109 const double coords[] = {0, 0, 0};
113 Range range_no_field_vertex;
114 range_no_field_vertex.insert(no_field_vertex);
118 CHKERR m_field.
get_moab().add_entities(meshset, range_no_field_vertex);
131 "MESH_NODE_POSITIONS");
157 "MESH_NODE_POSITIONS");
187 m_field,
"MESH_NODE_POSITIONS");
188 CHKERR m_field.
loop_dofs(
"MESH_NODE_POSITIONS", ent_method_mesh_positions);
207 typedef tee_device<std::ostream, std::ofstream>
TeeDevice;
210 std::ofstream ofs(
"forces_and_sources_testing_flat_prism_element.txt");
220 "FIELD1",
"FIELD1",
type),
230 const double eps = 1e-4;
231 for (DoubleAllocator::iterator it = getNormal().data().begin();
232 it != getNormal().data().end(); it++) {
233 *it = fabs(*it) <
eps ? 0.0 : *it;
235 for (DoubleAllocator::iterator it =
236 getNormalsAtGaussPtsF3().data().begin();
237 it != getNormalsAtGaussPtsF3().data().end(); it++) {
238 *it = fabs(*it) <
eps ? 0.0 : *it;
240 for (DoubleAllocator::iterator it =
241 getTangent1AtGaussPtF3().data().begin();
242 it != getTangent1AtGaussPtF3().data().end(); it++) {
243 *it = fabs(*it) <
eps ? 0.0 : *it;
245 for (DoubleAllocator::iterator it =
246 getTangent2AtGaussPtF3().data().begin();
247 it != getTangent2AtGaussPtF3().data().end(); it++) {
248 *it = fabs(*it) <
eps ? 0.0 : *it;
250 for (DoubleAllocator::iterator it =
251 getNormalsAtGaussPtsF4().data().begin();
252 it != getNormalsAtGaussPtsF4().data().end(); it++) {
253 *it = fabs(*it) <
eps ? 0.0 : *it;
255 for (DoubleAllocator::iterator it =
256 getTangent1AtGaussPtF4().data().begin();
257 it != getTangent1AtGaussPtF4().data().end(); it++) {
258 *it = fabs(*it) <
eps ? 0.0 : *it;
260 for (DoubleAllocator::iterator it =
261 getTangent2AtGaussPtF4().data().begin();
262 it != getTangent2AtGaussPtF4().data().end(); it++) {
263 *it = fabs(*it) <
eps ? 0.0 : *it;
266 mySplit <<
"NH1" << std::endl;
267 mySplit <<
"side: " << side <<
" type: " <<
type << std::endl;
268 mySplit << data << std::endl;
269 mySplit << std::setprecision(3) << getCoords() << std::endl;
270 mySplit << std::setprecision(3) << getCoordsAtGaussPts() << std::endl;
271 mySplit << std::setprecision(3) << getArea(0) << std::endl;
272 mySplit << std::setprecision(3) << getArea(1) << std::endl;
273 mySplit << std::setprecision(3) <<
"normal F3 " << getNormalF3()
275 mySplit << std::setprecision(3) <<
"normal F4 " << getNormalF4()
277 mySplit << std::setprecision(3) <<
"normal at Gauss pt F3 "
278 << getNormalsAtGaussPtsF3() << std::endl;
279 mySplit << std::setprecision(3) << getTangent1AtGaussPtF3()
281 mySplit << std::setprecision(3) << getTangent2AtGaussPtF3()
283 mySplit << std::setprecision(3) <<
"normal at Gauss pt F4 "
284 << getNormalsAtGaussPtsF4() << std::endl;
285 mySplit << std::setprecision(3) << getTangent1AtGaussPtF4()
287 mySplit << std::setprecision(3) << getTangent2AtGaussPtF4()
292 MoFEMErrorCode doWork(
int row_side,
int col_side, EntityType row_type,
301 mySplit <<
"NH1NH1" << std::endl;
302 mySplit <<
"row side: " << row_side <<
" row_type: " << row_type
304 mySplit << row_data << std::endl;
305 mySplit <<
"NH1NH1" << std::endl;
306 mySplit <<
"col side: " << col_side <<
" col_type: " << col_type
308 mySplit << row_data << std::endl;
320 "FIELD1",
"FIELD2",
type),
327 if (
type != MBENTITYSET)
330 mySplit <<
"NOFIELD" << std::endl;
331 mySplit <<
"side: " << side <<
" type: " <<
type << std::endl;
332 mySplit << data << std::endl;
336 MoFEMErrorCode doWork(
int row_side,
int col_side, EntityType row_type,
344 if (col_type != MBENTITYSET)
347 mySplit <<
"NOFILEDH1" << std::endl;
348 mySplit <<
"row side: " << row_side <<
" row_type: " << row_type
350 mySplit << row_data << std::endl;
351 mySplit <<
"col side: " << col_side <<
" col_type: " << col_type
353 mySplit << col_data << std::endl;
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
friend class UserDataOperator
Data on single entity (This is passed as argument to DataOperator::doWork)
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Problem manager is used to build and partition problems.
Create interface from given surface and insert flat prisms in-between.
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
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
const VectorDouble & getFieldData() const
get dofs values
Projection of edge entities with one mid-node on hierarchical basis.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
virtual EntityHandle get_field_meshset(const std::string name) const =0
get field meshset
@ OPROWCOL
operator doWork is executed on FE rows &columns
int main(int argc, char *argv[])
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.
Deprecated interface functions.
DeprecatedCoreInterface Interface
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
#define CHKERR
Inline error check.
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual moab::Interface & get_moab()=0
implementation of Data Operators for Forces and Sources
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 build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
@ OPCOL
operator doWork function is executed on FE columns
MoFEMErrorCode splitSides(const EntityHandle meshset, const BitRefLevel &bit, const int msId, const CubitBCType cubit_bc_type, const bool add_interface_entities, const bool recursive=false, int verb=QUIET)
Split nodes and other entities of tetrahedra on both sides of the interface and insert flat prisms in...
default operator for Flat Prism element
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
FlatPrism finite element.
#define CATCH_ERRORS
Catch errors.
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode getSides(const int msId, const CubitBCType cubit_bc_type, const BitRefLevel mesh_bit_level, const bool recursive, int verb=QUIET)
Store tetrahedra from each side of the interface separately in two child meshsets of the parent meshs...
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Operator used to check consistency between local coordinates and global cooridnates for integrated po...
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
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.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
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.
@ NOFIELD
scalar or vector of scalars describe (no true field)
@ OPROW
operator doWork function is executed on FE rows