|
| v0.14.0
|
Go to the documentation of this file.
5 namespace bio = boost::iostreams;
9 typedef tee_device<std::ostream, std::ofstream>
TeeDevice;
12 using namespace MoFEM;
14 static char help[] =
"...\n\n";
21 "FIELD1",
"FIELD1",
type, face_type),
22 faceType(face_type) {}
31 MOFEM_LOG(
"ATOM_TEST", Sev::inform) <<
"NH1";
33 <<
"side: " << side <<
" type: " <<
type;
34 MOFEM_LOG(
"ATOM_TEST", Sev::inform) << data;
36 if (faceType == FACEMASTER) {
38 <<
"coords Master " << std::fixed << std::setprecision(2)
41 <<
"area Master " << std::fixed << std::setprecision(2)
44 <<
"normal Master " << std::fixed << std::setprecision(2)
47 <<
"coords at Gauss Pts Master " << std::fixed
48 << std::setprecision(2) << getCoordsAtGaussPtsMaster();
51 <<
"coords Slave " << std::fixed << std::setprecision(2)
54 <<
"area Slave " << std::fixed << std::setprecision(2)
57 <<
"normal Slave " << std::fixed << std::setprecision(2)
60 <<
"coords at Gauss Pts Slave " << std::fixed
61 << std::setprecision(2) << getCoordsAtGaussPtsSlave();
75 MOFEM_LOG(
"ATOM_TEST", Sev::inform) <<
"NH1NH1";
77 <<
"row side: " << row_side <<
" row_type: " << row_type;
78 MOFEM_LOG(
"ATOM_TEST", Sev::inform) << row_data;
79 MOFEM_LOG(
"ATOM_TEST", Sev::inform) <<
"NH1NH1";
81 <<
"col side: " << col_side <<
" col_type: " << col_type;
82 MOFEM_LOG(
"ATOM_TEST", Sev::inform) << col_data;
100 MOFEM_LOG(
"ATOM_TEST", Sev::inform) <<
"Calling Operator NH1";
102 <<
"side: " << side <<
" type: " <<
type;
103 MOFEM_LOG(
"ATOM_TEST", Sev::inform) << data;
117 MOFEM_LOG(
"ATOM_TEST", Sev::inform) <<
"Calling Operator NH1NH1";
119 <<
"row side: " << row_side <<
" row_type: " << row_type;
120 MOFEM_LOG(
"ATOM_TEST", Sev::inform) << row_data;
121 MOFEM_LOG(
"ATOM_TEST", Sev::inform) <<
"NH1NH1";
123 <<
"col side: " << col_side <<
" col_type: " << col_type;
124 MOFEM_LOG(
"ATOM_TEST", Sev::inform) << col_data;
135 "FIELD1",
"FIELD2",
type, face_type) {}
141 if (
type != MBENTITYSET)
144 MOFEM_LOG(
"ATOM_TEST", Sev::inform) <<
"NPFIELD";
146 <<
"side: " << side <<
" type: " <<
type;
147 MOFEM_LOG(
"ATOM_TEST", Sev::inform) << data;
159 if (col_type != MBENTITYSET)
162 MOFEM_LOG(
"ATOM_TEST", Sev::inform) <<
"NOFILEDH1";
164 <<
"row side: " << row_side <<
" row_type: " << row_type;
165 MOFEM_LOG(
"ATOM_TEST", Sev::inform) << row_data;
167 <<
"col side: " << col_side <<
" col_type: " << col_type;
168 MOFEM_LOG(
"ATOM_TEST", Sev::inform) << col_data;
174 int main(
int argc,
char *argv[]) {
183 PetscBool flg = PETSC_TRUE;
187 if (flg != PETSC_TRUE)
190 enum spaces { MYH1, MYHDIV, MYLASTSPACEOP };
191 const char *list_spaces[] = {
"h1",
"hdiv"};
192 PetscInt choice_space_value =
H1;
194 MYLASTSPACEOP, &choice_space_value, &flg);
195 bool is_hdiv = (choice_space_value == MYHDIV) ?
true :
false;
202 auto add_atom_logging = [is_hdiv]() {
204 auto get_log_file_name = [is_hdiv]() {
206 return "forces_and_sources_testing_contact_prism_element_HDIV.txt";
208 return "forces_and_sources_testing_contact_prism_element.txt";
211 auto core_log = logging::core::get();
218 logging::add_file_log(keywords::file_name = get_log_file_name(),
220 MoFEM::LogKeywords::channel ==
"ATOM_TEST");
224 CHKERR add_atom_logging();
235 std::vector<BitRefLevel> bit_levels;
241 CHKERR PetscPrintf(PETSC_COMM_WORLD,
"Insert Interface %d\n",
242 cit->getMeshsetId());
247 CHKERR moab.create_meshset(MESHSET_SET, ref_level_meshset);
249 ->getEntitiesByTypeAndRefLevel(bit_levels.back(),
253 ->getEntitiesByTypeAndRefLevel(bit_levels.back(),
256 Range ref_level_tets;
257 CHKERR moab.get_entities_by_handle(ref_level_meshset, ref_level_tets,
260 CHKERR interface->
getSides(cubit_meshset, bit_levels.back(),
true, 0);
265 cubit_meshset,
true,
true, 0);
267 CHKERR moab.delete_entities(&ref_level_meshset, 1);
273 ->updateMeshsetByEntitiesChildren(cubit_meshset, bit_levels.back(),
274 cubit_meshset, MBMAXTYPE,
true);
288 auto set_no_field_vertex = [&]() {
291 const double coords[] = {0, 0, 0};
294 Range range_no_field_vertex;
295 range_no_field_vertex.insert(no_field_vertex);
299 CHKERR m_field.
get_moab().add_entities(meshset, range_no_field_vertex);
303 CHKERR set_no_field_vertex();
314 "MESH_NODE_POSITIONS");
338 "MESH_NODE_POSITIONS");
346 constexpr
int order = 3;
373 m_field,
"MESH_NODE_POSITIONS");
374 CHKERR m_field.
loop_dofs(
"MESH_NODE_POSITIONS", ent_method_mesh_positions);
393 using ContactDataOp =
398 new MyOp(UMDataOp::OPROW, ContactDataOp::FACEMASTER));
400 new MyOp(UMDataOp::OPROW, ContactDataOp::FACESLAVE));
402 new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACEMASTERMASTER));
404 new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACEMASTERSLAVE));
406 new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACESLAVEMASTER));
408 new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACESLAVESLAVE));
416 new MyOp2(UMDataOp::OPCOL, ContactDataOp::FACEMASTER));
418 new MyOp2(UMDataOp::OPCOL, ContactDataOp::FACESLAVE));
420 new MyOp2(UMDataOp::OPROWCOL, ContactDataOp::FACEMASTERMASTER));
422 new MyOp2(UMDataOp::OPROWCOL, ContactDataOp::FACEMASTERSLAVE));
424 new MyOp2(UMDataOp::OPROWCOL, ContactDataOp::FACESLAVEMASTER));
426 new MyOp2(UMDataOp::OPROWCOL, ContactDataOp::FACESLAVESLAVE));
#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.
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
MyOp2(const char type, const char face_type)
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
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
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
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
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.
MyOp(const char type, const char face_type)
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...
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
friend class UserDataOperator
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
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
structure to get information form mofem into EntitiesFieldData
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
#define MOFEM_LOG(channel, severity)
Log.
#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)
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
#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 ...
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
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.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
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.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
CallingOp(const char type)
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HDIV
field with continuous normal traction
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.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
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.
@ NOFIELD
scalar or vector of scalars describe (no true field)