v0.15.0
Loading...
Searching...
No Matches
forces_and_sources_testing_contact_prism_element.cpp
Go to the documentation of this file.
1
2
3#include <MoFEM.hpp>
4
5namespace bio = boost::iostreams;
6using bio::stream;
7using bio::tee_device;
8
9typedef tee_device<std::ostream, std::ofstream> TeeDevice;
10typedef stream<TeeDevice> TeeStream;
11
12using namespace MoFEM;
13
14static char help[] = "...\n\n";
15
17
18 const char faceType;
19 MyOp(const char type, const char face_type)
21 "FIELD1", "FIELD1", type, face_type),
22 faceType(face_type) {}
23
24 MoFEMErrorCode doWork(int side, EntityType type,
27
28 if (data.getFieldData().empty())
30
31 MOFEM_LOG("ATOM_TEST", Sev::inform) << "NH1";
32 MOFEM_LOG("ATOM_TEST", Sev::inform)
33 << "side: " << side << " type: " << type;
34 MOFEM_LOG("ATOM_TEST", Sev::inform) << data;
35
36 if (faceType == FACEMASTER) {
37 MOFEM_LOG("ATOM_TEST", Sev::inform)
38 << "coords Master " << std::fixed << std::setprecision(2)
39 << getCoordsMaster();
40 MOFEM_LOG("ATOM_TEST", Sev::inform)
41 << "area Master " << std::fixed << std::setprecision(2)
42 << getAreaMaster();
43 MOFEM_LOG("ATOM_TEST", Sev::inform)
44 << "normal Master " << std::fixed << std::setprecision(2)
45 << getNormalMaster();
46 MOFEM_LOG("ATOM_TEST", Sev::inform)
47 << "coords at Gauss Pts Master " << std::fixed
48 << std::setprecision(2) << getCoordsAtGaussPtsMaster();
49 } else {
50 MOFEM_LOG("ATOM_TEST", Sev::inform)
51 << "coords Slave " << std::fixed << std::setprecision(2)
52 << getCoordsSlave();
53 MOFEM_LOG("ATOM_TEST", Sev::inform)
54 << "area Slave " << std::fixed << std::setprecision(2)
55 << getAreaSlave();
56 MOFEM_LOG("ATOM_TEST", Sev::inform)
57 << "normal Slave " << std::fixed << std::setprecision(2)
58 << getNormalSlave();
59 MOFEM_LOG("ATOM_TEST", Sev::inform)
60 << "coords at Gauss Pts Slave " << std::fixed
61 << std::setprecision(2) << getCoordsAtGaussPtsSlave();
62 }
64 }
65
66 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
67 EntityType col_type,
71
72 if (row_data.getFieldData().empty())
74
75 MOFEM_LOG("ATOM_TEST", Sev::inform) << "NH1NH1";
76 MOFEM_LOG("ATOM_TEST", Sev::inform)
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";
80 MOFEM_LOG("ATOM_TEST", Sev::inform)
81 << "col side: " << col_side << " col_type: " << col_type;
82 MOFEM_LOG("ATOM_TEST", Sev::inform) << col_data;
83
85 }
86};
87
88struct CallingOp : public ForcesAndSourcesCore::UserDataOperator {
89
90 CallingOp(const char type)
91 : ForcesAndSourcesCore::UserDataOperator("FIELD1", "FIELD1", type) {}
92
93 MoFEMErrorCode doWork(int side, EntityType type,
96
97 if (data.getFieldData().empty())
99
100 MOFEM_LOG("ATOM_TEST", Sev::inform) << "Calling Operator NH1";
101 MOFEM_LOG("ATOM_TEST", Sev::inform)
102 << "side: " << side << " type: " << type;
103 MOFEM_LOG("ATOM_TEST", Sev::inform) << data;
104
106 }
107
108 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
109 EntityType col_type,
111 EntitiesFieldData::EntData &col_data) {
113
114 if (row_data.getFieldData().empty())
116
117 MOFEM_LOG("ATOM_TEST", Sev::inform) << "Calling Operator NH1NH1";
118 MOFEM_LOG("ATOM_TEST", Sev::inform)
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";
122 MOFEM_LOG("ATOM_TEST", Sev::inform)
123 << "col side: " << col_side << " col_type: " << col_type;
124 MOFEM_LOG("ATOM_TEST", Sev::inform) << col_data;
125
127 }
128};
129
130int main(int argc, char *argv[]) {
131
132 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
133
134 try {
135
136 moab::Core mb_instance;
137 moab::Interface &moab = mb_instance;
138
139 PetscBool flg = PETSC_TRUE;
140 char mesh_file_name[255];
141 CHKERR PetscOptionsGetString(PETSC_NULLPTR, "", "-my_file", mesh_file_name,
142 255, &flg);
143 if (flg != PETSC_TRUE)
144 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "error -my_file (mesh file not given");
145
146 enum spaces { MYH1, MYHDIV, MYLASTSPACEOP };
147 const char *list_spaces[] = {"h1", "hdiv"};
148 PetscInt choice_space_value = H1;
149 CHKERR PetscOptionsGetEList(PETSC_NULLPTR, NULL, "-space", list_spaces,
150 MYLASTSPACEOP, &choice_space_value, &flg);
151 bool is_hdiv = (choice_space_value == MYHDIV) ? true : false;
152
153 const char *option;
154 option = "";
155 CHKERR moab.load_file(mesh_file_name, 0, option);
156
157 // Create channel "ATOM_TEST" and add sink to file for that channel
158 auto add_atom_logging = [is_hdiv]() {
160 auto get_log_file_name = [is_hdiv]() {
161 if (is_hdiv)
162 return "forces_and_sources_testing_contact_prism_element_HDIV.txt";
163 else
164 return "forces_and_sources_testing_contact_prism_element.txt";
165 };
166
167 auto core_log = logging::core::get();
168 core_log->add_sink(
170 LogManager::setLog("ATOM_TEST");
171 MOFEM_LOG_TAG("ATOM_TEST", "atom test");
172
173 // Add log to file
174 logging::add_file_log(keywords::file_name = get_log_file_name(),
175 keywords::filter =
176 MoFEM::LogKeywords::channel == "ATOM_TEST");
178 };
179
180 CHKERR add_atom_logging();
181
182 // Create MoFEM (Joseph) database
183 MoFEM::Core core(moab);
184 MoFEM::Interface &m_field = core;
185 PrismInterface *interface;
186 CHKERR m_field.getInterface(interface);
187
188 // set entities bit level
189 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
190 0, 3, BitRefLevel().set(0));
191 std::vector<BitRefLevel> bit_levels;
192 bit_levels.push_back(BitRefLevel().set(0));
193
194 int ll = 1;
196
197 CHKERR PetscPrintf(PETSC_COMM_WORLD, "Insert Interface %d\n",
198 cit->getMeshsetId());
199 EntityHandle cubit_meshset = cit->getMeshset();
200 {
201 // get tet enties form back bit_level
202 EntityHandle ref_level_meshset = 0;
203 CHKERR moab.create_meshset(MESHSET_SET, ref_level_meshset);
205 ->getEntitiesByTypeAndRefLevel(bit_levels.back(),
206 BitRefLevel().set(), MBTET,
207 ref_level_meshset);
209 ->getEntitiesByTypeAndRefLevel(bit_levels.back(),
210 BitRefLevel().set(), MBPRISM,
211 ref_level_meshset);
212 Range ref_level_tets;
213 CHKERR moab.get_entities_by_handle(ref_level_meshset, ref_level_tets,
214 true);
215 // get faces and test to split
216 CHKERR interface->getSides(cubit_meshset, bit_levels.back(), true, 0);
217 // set new bit level
218 bit_levels.push_back(BitRefLevel().set(ll++));
219 // split faces and
220 CHKERR interface->splitSides(ref_level_meshset, bit_levels.back(),
221 cubit_meshset, true, true, 0);
222 // clean meshsets
223 CHKERR moab.delete_entities(&ref_level_meshset, 1);
224 }
225 // update cubit meshsets
226 for (_IT_CUBITMESHSETS_FOR_LOOP_(m_field, ciit)) {
227 EntityHandle cubit_meshset = ciit->meshset;
229 ->updateMeshsetByEntitiesChildren(cubit_meshset, bit_levels.back(),
230 cubit_meshset, MBMAXTYPE, true);
231 }
232 }
233
234 // Fields
235 if (is_hdiv)
236 CHKERR m_field.add_field("FIELD1", HDIV, DEMKOWICZ_JACOBI_BASE, 3);
237 else
238 CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 3);
239
240 CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
241 3);
242
243 // FE
244 CHKERR m_field.add_finite_element("TEST_FE1");
245
246 // Define rows/cols and element data
247 CHKERR m_field.modify_finite_element_add_field_row("TEST_FE1", "FIELD1");
248 CHKERR m_field.modify_finite_element_add_field_col("TEST_FE1", "FIELD1");
249 CHKERR m_field.modify_finite_element_add_field_data("TEST_FE1", "FIELD1");
251 "MESH_NODE_POSITIONS");
252
253 // Problem
254 CHKERR m_field.add_problem("TEST_PROBLEM");
255
256 // set finite elements for problem
257 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
258 "TEST_FE1");
259 // set refinement level for problem
260 CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM",
261 bit_levels.back());
262
263 // meshset consisting all entities in mesh
264 EntityHandle root_set = moab.get_root_set();
265 // add entities to field
266 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD1");
267 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET,
268 "MESH_NODE_POSITIONS");
269 // add entities to finite element
270 CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBPRISM,
271 "TEST_FE1");
272
273 // set app. order
274 constexpr int order = 3;
275
276 if (is_hdiv) {
277 CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", 0);
278 CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
279 } else {
280 CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", order);
281 CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
282 CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", order);
283 CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
284 }
285
286 CHKERR m_field.set_field_order(root_set, MBTET, "MESH_NODE_POSITIONS", 2);
287 CHKERR m_field.set_field_order(root_set, MBTRI, "MESH_NODE_POSITIONS", 2);
288 CHKERR m_field.set_field_order(root_set, MBEDGE, "MESH_NODE_POSITIONS", 2);
289 CHKERR m_field.set_field_order(root_set, MBVERTEX, "MESH_NODE_POSITIONS",
290 1);
291
292 // build field
293 CHKERR m_field.build_fields();
294
295 if (!is_hdiv) {
296 // set FIELD1 from positions of 10 node tets
297 Projection10NodeCoordsOnField ent_method_field1(m_field, "FIELD1");
298 CHKERR m_field.loop_dofs("FIELD1", ent_method_field1);
299 }
300 Projection10NodeCoordsOnField ent_method_mesh_positions(
301 m_field, "MESH_NODE_POSITIONS");
302 CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method_mesh_positions);
303
304 // build finite elemnts
306 // build adjacencies
307 CHKERR m_field.build_adjacencies(bit_levels.back());
308 // build problem
309 ProblemsManager *prb_mng_ptr;
310 CHKERR m_field.getInterface(prb_mng_ptr);
311 CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", false);
312
313 // mesh partitioning
314 // partition
315 CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
316 CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
317 // what are ghost nodes, see Petsc Manual
318 CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
319
320 using UMDataOp = ForcesAndSourcesCore::UserDataOperator;
321 using ContactDataOp =
323
325 fe1.getOpPtrVector().push_back(
326 new MyOp(UMDataOp::OPROW, ContactDataOp::FACEMASTER));
327 fe1.getOpPtrVector().push_back(
328 new MyOp(UMDataOp::OPROW, ContactDataOp::FACESLAVE));
329 fe1.getOpPtrVector().push_back(
330 new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACEMASTERMASTER));
331 fe1.getOpPtrVector().push_back(
332 new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACEMASTERSLAVE));
333 fe1.getOpPtrVector().push_back(
334 new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACESLAVEMASTER));
335 fe1.getOpPtrVector().push_back(
336 new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACESLAVESLAVE));
337 fe1.getOpPtrVector().push_back(new CallingOp(UMDataOp::OPCOL));
338 fe1.getOpPtrVector().push_back(new CallingOp(UMDataOp::OPROW));
339 fe1.getOpPtrVector().push_back(new CallingOp(UMDataOp::OPROWCOL));
340 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE1", fe1);
341
342 }
344
346
347 return 0;
348}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
int main()
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ H1
continuous field
Definition definitions.h:85
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ SIDESET
@ MOFEM_INVALID_DATA
Definition definitions.h:36
#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 ...
constexpr int order
tee_device< std::ostream, std::ofstream > TeeDevice
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 build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
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.
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.
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.
#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 partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
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 add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
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
char mesh_file_name[255]
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Managing BitRefLevels.
MatrixDouble & getCoordsAtGaussPtsMaster()
get coordinates at Gauss pts on full prism.
MatrixDouble & getCoordsAtGaussPtsSlave()
get coordinates at Gauss pts on full prism.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
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.
Core (interface) class.
Definition Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:118
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorDouble & getFieldData() const
get dofs values
structure to get information form mofem into EntitiesFieldData
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Create interface from given surface and insert flat prisms in-between.
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...
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...
Problem manager is used to build and partition problems.
Projection of edge entities with one mid-node on hierarchical basis.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Operator used to check consistency between local coordinates and global cooridnates for integrated po...
MyOp(const char type, const char face_type)
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 doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.