v0.15.0
Loading...
Searching...
No Matches
SolutionMapping.hpp
Go to the documentation of this file.
1
2
3/** \file SolutionMapping.hpp
4 * \example SolutionMapping.hpp
5 */
6
7#include <DGProjection.hpp>
9
11
12 MoFEMErrorCode operator()(ForcesAndSourcesCore *fe_raw_ptr, int order_row,
13 int order_col, int order_data);
14
15private:
16 struct Fe : public ForcesAndSourcesCore {
17 using ForcesAndSourcesCore::dataOnElement;
18
19 private:
20 using ForcesAndSourcesCore::ForcesAndSourcesCore;
21 };
22
23 static inline std::map<long int, MatrixDouble> mapRefCoords;
24};
25
27
29
30 MoFEMErrorCode mapFields(SmartPetscObj<DM> &intermediateDM,
31 BitRefLevel parentBit, BitRefLevel childBit,
32 PetscInt numVecs = 1,
33 Vec vecsToMapFrom[] = PETSC_NULLPTR,
34 Vec vecsToMapTo[] = PETSC_NULLPTR);
35 MoFEMErrorCode reSetUp(std::vector<std::string> L2Fields,
36 std::vector<int> L2Oders,
37 std::vector<std::string> H1Fields,
38 std::vector<int> H1Orders,
39 std::vector<std::string> HdivFields,
40 std::vector<int> HdivOrders, BitRefLevel bIt);
41
42private:
44 BitRefLevel parentBit, childBit;
45 SmartPetscObj<DM> subDM, intermediateDM;
46
47 boost::shared_ptr<DomainEle> feRhs;
48 boost::shared_ptr<DomainEle> feLhs;
49
50 MoFEMErrorCode setUp(SmartPetscObj<DM> &subDM,
51 SmartPetscObj<DM> &intermediateDM,
52 BitRefLevel parent_bit, BitRefLevel child_bit);
53 MoFEMErrorCode projectResults(BitRefLevel parent_bit, BitRefLevel child_bit);
54 MoFEMErrorCode solveMap(SmartPetscObj<DM> &subDM,
55 SmartPetscObj<DM> &intermediateDM,
56 BitRefLevel parent_bit, BitRefLevel child_bit,
57 PetscInt numVecs, Vec vecsToMapFrom[],
58 Vec vecsToMapTo[]);
59 MoFEMErrorCode postProcess();
60};
61
63 feRhs = boost::make_shared<DomainEle>(m_field);
64 feLhs = boost::make_shared<DomainEle>(m_field);
65}
66
67MoFEMErrorCode SolutionMapping::mapFields(SmartPetscObj<DM> &intermediate_dm,
68 BitRefLevel parent_bit,
69 BitRefLevel child_bit,
70 PetscInt num_vecs,
71 Vec vecs_to_map_from[],
72 Vec vecs_to_map_to[]) {
74
75 CHKERR setUp(subDM, intermediate_dm, parent_bit, child_bit);
76 CHKERR solveMap(subDM, intermediate_dm, parent_bit, child_bit, num_vecs,
77 vecs_to_map_from, vecs_to_map_to);
78
80};
81
82MoFEMErrorCode SolutionMapping::setUp(SmartPetscObj<DM> &sub_dm,
83 SmartPetscObj<DM> &intermediate_dm,
84 BitRefLevel parent_bit,
85 BitRefLevel child_bit) {
87
88 MOFEM_LOG("REMESHING", Sev::inform) << "Set up sub DM for mapping";
89
90 Simple *simple = m_field.getInterface<Simple>();
91
92 Range sub_dm_elems;
93
94 auto create_sub_dm = [&]() {
96 sub_dm = createDM(m_field.get_comm(), "DMMOFEM");
97 auto create_impl = [&]() {
99 CHKERR DMSetType(sub_dm, "DMMOFEM");
100 CHKERR DMMoFEMCreateSubDM(sub_dm, intermediate_dm, "SUB_MAPPING");
101 CHKERR DMMoFEMAddElement(sub_dm, simple->getDomainFEName());
102 CHKERR DMMoFEMSetSquareProblem(sub_dm, PETSC_TRUE);
103 CHKERR DMMoFEMSetDestroyProblem(sub_dm, PETSC_TRUE);
104
105 MOFEM_LOG("REMESHING", Sev::inform) << "Created sub DM";
106
107 // get only refinement bit DOFs
108 auto ref_entities_ptr = boost::make_shared<Range>();
109 auto child_mask = parent_bit;
110 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByDimAndRefLevel(
111 child_bit, child_mask.flip(), SPACE_DIM, *ref_entities_ptr);
112
113 Range edges;
114 CHKERR m_field.get_moab().get_adjacencies(*ref_entities_ptr,
115 SPACE_DIM - 1, true, edges,
116 moab::Interface::UNION);
117 ref_entities_ptr->merge(edges);
118
119 Range verts;
120 CHKERR m_field.get_moab().get_adjacencies(*ref_entities_ptr, 0, true,
121 verts, moab::Interface::UNION);
122 ref_entities_ptr->merge(verts);
123
124 for (auto f : {"U", "FLUX"}) {
125 CHKERR DMMoFEMAddSubFieldRow(sub_dm, f, ref_entities_ptr);
126 CHKERR DMMoFEMAddSubFieldCol(sub_dm, f, ref_entities_ptr);
127 }
128
130 parent_bit | child_bit);
132 BitRefLevel().set());
133
134 CHKERR DMMoFEMSetIsPartitioned(sub_dm, is_distributed_mesh);
135 CHKERR DMSubDMSetUp_MoFEM(sub_dm);
136
138 };
139
140 CHKERR create_impl();
142 };
143
144 CHK_THROW_MESSAGE(create_sub_dm(), "failed to create sub dm");
145
146 MOFEM_LOG("REMESHING", Sev::inform) << "Set up sub DM for mapping <= done";
148};
149
150MoFEMErrorCode SolutionMapping::reSetUp(std::vector<std::string> L2_fields,
151 std::vector<int> L2_orders,
152 std::vector<std::string> H1_fields,
153 std::vector<int> H1_orders,
154 std::vector<std::string> Hdiv_fields,
155 std::vector<int> Hdiv_orders,
156 BitRefLevel bit) {
158
159 Simple *simple = m_field.getInterface<Simple>();
160
161 auto add_ents_order_to_field =
163 const std::initializer_list<std::string> &f,
164 const std::initializer_list<EntityType> &t, int _order) {
166 for (auto _f : f) {
167 for (auto _t : t) {
169 CHKERR m_field.set_field_order(0, _t, _f, _order);
170 }
171 }
173 };
174
175 for (size_t i = 0; i < L2_fields.size(); ++i) {
176 CHKERR add_ents_order_to_field(m_field, {L2_fields[i]}, {MBTRI},
177 L2_orders[i]);
178 }
179 for (size_t i = 0; i < H1_fields.size(); ++i) {
180 CHKERR add_ents_order_to_field(m_field, {H1_fields[i]}, {MBTRI, MBEDGE},
181 H1_orders[i]);
182 CHKERR m_field.set_field_order(0, MBVERTEX, {H1_fields[i]},
183 1); // H1 only supports 1st order on vertices
184 }
185 for (size_t i = 0; i < Hdiv_fields.size(); ++i) {
186 CHKERR add_ents_order_to_field(m_field, {Hdiv_fields[i]}, {MBTRI, MBEDGE},
187 Hdiv_orders[i]);
188 }
189
190 Skinner skin(&m_field.get_moab());
191 Range faces;
192 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
193 BitRefLevel().set(FLIPPED_BIT + num_refinement_levels),
194 BitRefLevel().set(), MBTRI, faces);
195 Range skin_edges;
196 CHKERR skin.find_skin(0, faces, false, skin_edges);
197
198#ifdef ADD_CONTACT
199
200 auto filter_blocks = [&](auto skin) {
201 bool is_contact_block = false;
202 Range contact_range;
203 for (auto m :
204 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
205
206 (boost::format("%s(.*)") % "CONTACT").str()
207
208 ))
209
210 ) {
211 is_contact_block =
212 true; ///< blocs interation is collective, so that is set irrespective
213 ///< if there are entities in given rank or not in the block
214 MOFEM_LOG("CONTACT", Sev::inform)
215 << "Find contact block set: " << m->getName();
216 auto meshset = m->getMeshset();
217 Range contact_meshset_range;
218 CHKERR m_field.get_moab().get_entities_by_dimension(
219 meshset, SPACE_DIM - 1, contact_meshset_range, true);
220
221 CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
222 contact_meshset_range);
223 contact_range.merge(contact_meshset_range);
224 }
225 if (is_contact_block) {
226 MOFEM_LOG("SYNC", Sev::inform)
227 << "Nb entities in contact surface: " << contact_range.size();
229 skin = intersect(skin, contact_range);
230 }
231 return skin;
232 };
233
234 auto sigma_trace_ents = filter_blocks(skin_edges);
235 CHKERR m_field.add_ents_to_field_by_type(0, MBTRI, "SIGMA");
236 CHKERR m_field.set_field_order(0, MBTRI, "SIGMA", 0);
237 CHKERR m_field.add_ents_to_field_by_type(sigma_trace_ents, MBEDGE, "SIGMA");
238 CHKERR m_field.set_field_order(sigma_trace_ents, "SIGMA", order - 1);
239#endif // ADD_CONTACT
240
242 simple->getDomainFEName());
244 skin_edges, MBEDGE, simple->getBoundaryFEName());
245
249
251};
252
253MoFEMErrorCode SolutionMapping::projectResults(BitRefLevel parent_bit,
254 BitRefLevel child_bit) {
256
257 feRhs->getOpPtrVector().clear();
258 feLhs->getOpPtrVector().clear();
259
260 MOFEM_LOG("REMESHING", Sev::inform) << "Cleared FEs";
261
262 Simple *simple = m_field.getInterface<Simple>();
263
264 auto vol_rule = [](int, int, int ao) { return 2 * ao + geom_order; };
265
266 feRhs->getRuleHook = vol_rule;
267 feLhs->getRuleHook = vol_rule;
268
269 // OpLoopThis, is child operator, and is use to execute
270 // fe_child_ptr, only on bit ref level and mask
271 // for child elements
272 auto get_child_op = [&](auto &pip) {
273 auto op_this_child =
274 new OpLoopThis<DomainEle>(m_field, simple->getDomainFEName(), child_bit,
275 parent_bit.flip(), Sev::noisy);
276 auto fe_child_ptr = op_this_child->getThisFEPtr();
277 fe_child_ptr->getRuleHook = [](int, int, int p) { return -1; };
278 Range child_edges;
279 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
280 child_bit, parent_bit.flip(), MBEDGE, child_edges);
281 // set integration rule, such that integration points are not on flipped edge
282 CHKERR setDGSetIntegrationPoints<2>( // Currently only implemented for 2D
283 fe_child_ptr->setRuleHook, [](int, int, int p) { return 2 * p; },
284 boost::make_shared<Range>(child_edges));
285 pip.push_back(op_this_child);
286 return fe_child_ptr;
287 };
288
289 // Use field evaluator to calculate field values on parent bitref level,
290 // i.e. elements which were flipped.
291 auto get_field_eval_op = [&](auto fe_child_ptr) {
292 auto field_eval_ptr = m_field.getInterface<FieldEvaluatorInterface>();
293
294 // Get pointer of FieldEvaluator data. Note finite element and method
295 // set integrating points is destroyed when this pointer is releases
296 auto field_eval_data = field_eval_ptr->getData<DomainEle>();
297 // Build tree for particular element
298 CHKERR field_eval_ptr->buildTree<SPACE_DIM>(
299 field_eval_data, simple->getDomainFEName(), parent_bit,
300 BitRefLevel().set());
301
302 // You can add more fields here
303 auto new_T_ptr = boost::make_shared<MatrixDouble>();
304 auto old_T_ptr = boost::make_shared<MatrixDouble>();
305 auto new_TAU_ptr = boost::make_shared<MatrixDouble>();
306 auto old_TAU_ptr = boost::make_shared<MatrixDouble>();
307 auto new_EP_ptr = boost::make_shared<MatrixDouble>();
308 auto old_EP_ptr = boost::make_shared<MatrixDouble>();
309
310 auto new_U_ptr = boost::make_shared<MatrixDouble>();
311 auto old_U_ptr = boost::make_shared<MatrixDouble>();
312 auto new_flux_ptr = boost::make_shared<MatrixDouble>();
313 auto old_flux_ptr = boost::make_shared<MatrixDouble>();
314
315 if (auto fe_eval_ptr = field_eval_data->feMethodPtr.lock()) {
316 fe_eval_ptr->getRuleHook = [](int, int, int p) { return -1; };
317 fe_eval_ptr->getOpPtrVector().push_back(
318 new OpCalculateVectorFieldValues<1>("T", old_T_ptr));
319 fe_eval_ptr->getOpPtrVector().push_back(
320 new OpCalculateVectorFieldValues<1>("TAU", old_TAU_ptr));
321 constexpr int size_symm = (SPACE_DIM * (SPACE_DIM + 1)) / 2;
322 fe_eval_ptr->getOpPtrVector().push_back(
323 new OpCalculateVectorFieldValues<size_symm>("EP", old_EP_ptr));
324 fe_eval_ptr->getOpPtrVector().push_back(
325 new OpCalculateVectorFieldValues<SPACE_DIM>("U", old_U_ptr));
326 fe_eval_ptr->getOpPtrVector().push_back(
327 new OpCalculateHVecVectorField<3>("FLUX", old_flux_ptr));
328
329 } else {
331 "Field evaluator method pointer is expired");
332 }
333
334 auto op_ptr = field_eval_ptr->getDataOperator<SPACE_DIM>(
335 {
336
337 {old_T_ptr, new_T_ptr},
338 {old_TAU_ptr, new_TAU_ptr},
339 {old_EP_ptr, new_EP_ptr},
340 {old_U_ptr, new_U_ptr},
341 {old_flux_ptr, new_flux_ptr}
342
343 },
344 simple->getDomainFEName(), field_eval_data, 0, m_field.get_comm_size(),
345 parent_bit, child_bit.flip(), MF_EXIST, QUIET);
346
347 fe_child_ptr->getOpPtrVector().push_back(op_ptr);
348
349 using FieldTupleVec = std::vector<
350 std::tuple<std::string, boost::shared_ptr<MatrixDouble>, int>>;
351
352 return std::make_pair(
353
354 std::vector<std::pair<std::string, boost::shared_ptr<MatrixDouble>>>{
355 {"T", new_T_ptr}, {"TAU", new_TAU_ptr}, {"EP", new_EP_ptr}},
356
357 std::vector<std::pair<std::string, boost::shared_ptr<MatrixDouble>>>{
358 {"U", new_U_ptr}, {"FLUX", new_flux_ptr}}
359
360 );
361 };
362
363 // calculate coefficients on child (flipped) elements
364 auto dg_projection_base = [&](auto fe_child_ptr, auto vec_data_ptr) {
366 // constexpr int projection_order = order;
367 auto entity_data_l2 = boost::make_shared<EntitiesFieldData>(MBENTITYSET);
368 auto mass_ptr = boost::make_shared<MatrixDouble>();
369 auto coeffs_ptr = boost::make_shared<MatrixDouble>();
370
371 constexpr int projection_order = 2; // use order 2 for DG projection
372
373 for (auto &[field_name, data_ptr] : vec_data_ptr.first) {
374
375 fe_child_ptr->getOpPtrVector().push_back(new OpDGProjectionMassMatrix(
376 projection_order, mass_ptr, entity_data_l2, AINSWORTH_LEGENDRE_BASE,
377 L2));
378 fe_child_ptr->getOpPtrVector().push_back(new OpDGProjectionCoefficients(
379 data_ptr, coeffs_ptr, mass_ptr, entity_data_l2,
381
382 // set coefficients to flipped element
383 auto op_set_coeffs = new DomainEleOp(field_name, DomainEleOp::OPROW);
384 op_set_coeffs->doWorkRhsHook =
385 [coeffs_ptr](DataOperator *base_op_ptr, int side, EntityType type,
386 EntitiesFieldData::EntData &data) -> MoFEMErrorCode {
388 auto field_ents = data.getFieldEntities();
389 auto nb_dofs = data.getIndices().size();
390 if (!field_ents.size())
392 if (auto e_ptr = field_ents[0]) {
393 auto field_ent_data = e_ptr->getEntFieldData();
394 std::copy(coeffs_ptr->data().data(),
395 coeffs_ptr->data().data() + nb_dofs,
396 field_ent_data.begin());
397 }
399 };
400 fe_child_ptr->getOpPtrVector().push_back(op_set_coeffs);
401 }
402
403 for (auto &p : vec_data_ptr.second) {
404 auto field_name = p.first;
405 auto data_ptr = p.second;
406
407 if (field_name == "U") {
408
409 using OpDomainMass = FormsIntegrators<DomainEleOp>::Assembly<
411
412 // assemble to global matrix, since this is H1 (you will do the shame for Hcurl of Hdiv)
413 auto beta = [](const double, const double, const double) { return 1; };
414 fe_child_ptr->getOpPtrVector().push_back(
416 using OpVec = FormsIntegrators<DomainEleOp>::Assembly<
418 fe_child_ptr->getOpPtrVector().push_back(
419 new OpVec(field_name, data_ptr, beta));
420 }
421
422 if (field_name == "FLUX") {
423
424 using OpDomainMass = FormsIntegrators<DomainEleOp>::Assembly<
426
427 // assemble to global matrix, since this is H1 (you will do the shame for Hcurl of Hdiv)
428 auto beta = [](const double, const double, const double) { return 1; };
429 fe_child_ptr->getOpPtrVector().push_back(
431 using OpVec = FormsIntegrators<DomainEleOp>::Assembly<
433 fe_child_ptr->getOpPtrVector().push_back(
434 new OpVec(field_name, data_ptr, beta));
435 }
436 }
437
439 };
440
441 // create child operator, and fe_child_ptr element in it
442 auto fe_child_ptr = get_child_op(feRhs->getOpPtrVector());
443 // run dg projection, note that get_field_eval_op,
444 // pass data_ptr values used to project and calculate coefficients
445 CHKERR dg_projection_base(fe_child_ptr, get_field_eval_op(fe_child_ptr));
446
447 MOFEM_LOG("REMESHING", Sev::inform) << "Before looping FEs";
448 CHKERR DMoFEMLoopFiniteElementsUpAndLowRank(
449 subDM, simple->getDomainFEName(), feRhs, 0, m_field.get_comm_size());
450 // CHKERR DMoFEMLoopFiniteElementsUpAndLowRank(subDM, simple->getDomainFEName(),
451 // feRhs, m_field.get_comm_rank(),
452 // m_field.get_comm_rank());
453 MOFEM_LOG("REMESHING", Sev::inform) << "After looping RHS";
454 // CHKERR DMoFEMLoopFiniteElementsUpAndLowRank(
455 // subDM, simple->getDomainFEName(), feLhs, 0, m_field.get_comm_size());
456 // MOFEM_LOG("REMESHING", Sev::inform) << "After looping LHS";
457
459};
460
462 SmartPetscObj<DM> &sub_dm, SmartPetscObj<DM> &intermediate_dm,
463 BitRefLevel parent_bit, BitRefLevel child_bit, PetscInt num_vecs,
464 Vec vecs_to_map_from[], Vec vecs_to_map_to[]) {
466
467 auto simple = m_field.getInterface<Simple>();
468
469 auto post_proc = [&](auto dm, std::string file_name) {
471
472 auto T_ptr = boost::make_shared<VectorDouble>();
473 auto TAU_ptr = boost::make_shared<VectorDouble>();
474 auto EP_ptr = boost::make_shared<MatrixDouble>();
475
476 auto U_ptr = boost::make_shared<MatrixDouble>();
477 auto FLUX_ptr = boost::make_shared<MatrixDouble>();
478
479 auto post_proc_fe = boost::make_shared<PostProcEle>(m_field);
480
481 using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
482
483 post_proc_fe->getOpPtrVector().push_back(
484 new OpCalculateScalarFieldValues("T", T_ptr));
485 post_proc_fe->getOpPtrVector().push_back(
486 new OpCalculateScalarFieldValues("TAU", TAU_ptr));
487 post_proc_fe->getOpPtrVector().push_back(
488 new OpCalculateTensor2SymmetricFieldValues<SPACE_DIM>("EP", EP_ptr));
489 post_proc_fe->getOpPtrVector().push_back(
490 new OpCalculateVectorFieldValues<SPACE_DIM>("U", U_ptr));
491 post_proc_fe->getOpPtrVector().push_back(
492 new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", FLUX_ptr));
493
494 post_proc_fe->getOpPtrVector().push_back(new OpPPMap(
495 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
496
497 {{"T", T_ptr}, {"TAU", TAU_ptr}},
498
499 {{"U", U_ptr}, {"FLUX", FLUX_ptr}},
500
501 {},
502
503 {{"EP", EP_ptr}}
504
505 )
506
507 );
508
509 MOFEM_LOG("WORLD", Sev::inform) << "Just before looping post_proc_fe";
510 {
511 CHKERR DMoFEMLoopFiniteElementsUpAndLowRank(dm, simple->getDomainFEName(),
512 post_proc_fe, 0,
514 auto &post_proc_moab = post_proc_fe->getPostProcMesh();
515 auto output_name = "out_all_rank_" +
516 std::to_string(m_field.get_comm_rank()) + "_" +
517 file_name + ".h5m";
518 CHKERR post_proc_moab.write_file(output_name.c_str());
519 }
520 {
521 CHKERR DMoFEMLoopFiniteElementsUpAndLowRank(
522 dm, simple->getDomainFEName(), post_proc_fe, m_field.get_comm_rank(),
524 auto &post_proc_moab = post_proc_fe->getPostProcMesh();
525 auto output_name = "out_only_on_rank" +
526 std::to_string(m_field.get_comm_rank()) + "_" +
527 file_name + ".h5m";
528 CHKERR post_proc_moab.write_file(output_name.c_str());
529 }
530
532 };
533
534 auto null_fe = boost::shared_ptr<FEMethod>();
535
536 // CHKERR DMMoFEMKSPSetComputeOperators(sub_dm, simple->getDomainFEName(), feLhs,
537 // null_fe, null_fe);
538 // CHKERR DMMoFEMKSPSetComputeRHS(sub_dm, simple->getDomainFEName(), feRhs,
539 // null_fe, null_fe);
540
541
542
543 auto vec = createDMVector(sub_dm);
544 int size;
545 CHKERR VecGetSize(vec, &size);
546 SmartPetscObj<Mat> mat;
547 if (size)
548 mat = createDMMatrix(sub_dm);
549 auto sol = createDMVector(sub_dm);
550
551 auto ksp = createKSP(m_field.get_comm());
552 CHKERR KSPSetOperators(ksp, mat, mat);
553 CHKERR KSPSetFromOptions(ksp);
554
555 feRhs->ksp_A = mat;
556 feRhs->ksp_B = mat;
557 feRhs->ksp_f = vec;
558 feRhs->data_ctx =
559 PetscData::CtxSetF | PetscData::CtxSetA | PetscData::CtxSetB;
560
561 // We create a temp_D to store current state of the mesh. We doing to preserve
562 // whatever is currently stored on entities, we will restore it to keep it as
563 // it was before.
564
565 auto tmp_D = createDMVector(intermediate_dm);
566 CHKERR DMoFEMMeshToLocalVector(intermediate_dm, tmp_D, INSERT_VALUES,
567 SCATTER_FORWARD);
568 CHKERR VecGhostUpdateBegin(tmp_D, INSERT_VALUES, SCATTER_FORWARD);
569 CHKERR VecGhostUpdateEnd(tmp_D, INSERT_VALUES, SCATTER_FORWARD);
570
571 for (int i = 0; i < num_vecs; ++i) {
572 // load data to intermediate dm from vecs_to_map_from[i]
573 CHKERR VecGhostUpdateBegin(vecs_to_map_from[i], INSERT_VALUES,
574 SCATTER_FORWARD);
575 CHKERR VecGhostUpdateEnd(vecs_to_map_from[i], INSERT_VALUES,
576 SCATTER_FORWARD);
577 CHKERR DMoFEMMeshToGlobalVector(intermediate_dm, vecs_to_map_from[i],
578 INSERT_VALUES, SCATTER_REVERSE);
579
580 if (size)
581 CHKERR MatZeroEntries(mat);
582 CHKERR VecZeroEntries(vec);
583 CHKERR VecGhostUpdateBegin(vec, INSERT_VALUES, SCATTER_FORWARD);
584 CHKERR VecGhostUpdateEnd(vec, INSERT_VALUES, SCATTER_FORWARD);
585 // map
586 CHKERR projectResults(parent_bit, child_bit);
587
588 if (size) {
589 CHKERR VecAssemblyBegin(vec);
590 CHKERR VecAssemblyEnd(vec);
591 CHKERR MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);
592 CHKERR MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);
593
594 CHKERR KSPSolve(ksp, vec, sol);
595 CHKERR VecGhostUpdateBegin(sol, INSERT_VALUES, SCATTER_FORWARD);
596 CHKERR VecGhostUpdateEnd(sol, INSERT_VALUES, SCATTER_FORWARD);
597 CHKERR DMoFEMMeshToGlobalVector(sub_dm, sol, INSERT_VALUES,
598 SCATTER_REVERSE);
599 }
600
601 // update data
602 CHKERR DMoFEMMeshToLocalVector(intermediate_dm, vecs_to_map_to[i],
603 INSERT_VALUES, SCATTER_FORWARD);
604 CHKERR VecGhostUpdateBegin(vecs_to_map_to[i], INSERT_VALUES,
605 SCATTER_FORWARD);
606 CHKERR VecGhostUpdateEnd(vecs_to_map_to[i], INSERT_VALUES, SCATTER_FORWARD);
607
608 CHKERR post_proc(intermediate_dm, "map_" + std::to_string(i));
609 }
610
611 CHKERR DMoFEMMeshToLocalVector(intermediate_dm, tmp_D, INSERT_VALUES,
612 SCATTER_REVERSE);
613
614 MOFEM_LOG("REMESHING", Sev::inform) << "Mapped thermoplastic state variables";
615
617};
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
constexpr int SPACE_DIM
[Define dimension]
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
@ QUIET
@ MF_EXIST
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ L2
field with C-1 continuity
Definition definitions.h:88
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr int order
Order displacement.
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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 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.
#define MOFEM_LOG(channel, severity)
Log.
virtual MoFEMErrorCode modify_problem_mask_ref_level_set_bit(const std::string &name_problem, const BitRefLevel &bit)=0
set dof mask ref level for problem
virtual MoFEMErrorCode modify_problem_ref_level_set_bit(const std::string &name_problem, const BitRefLevel &bit)=0
set ref level for problem
auto bit
set bit
FTensor::Index< 'i', SPACE_DIM > i
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
constexpr double t
plate stiffness
Definition plate.cpp:58
constexpr auto field_name
FTensor::Index< 'm', 3 > m
virtual int get_comm_size() const =0
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Deprecated interface functions.
Post post-proc data at points from hash maps.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
MoFEMErrorCode operator()(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)
static std::map< long int, MatrixDouble > mapRefCoords
MoFEMErrorCode solveMap(SmartPetscObj< DM > &subDM, SmartPetscObj< DM > &intermediateDM, BitRefLevel parent_bit, BitRefLevel child_bit, PetscInt numVecs, Vec vecsToMapFrom[], Vec vecsToMapTo[])
MoFEMErrorCode postProcess()
SmartPetscObj< DM > intermediateDM
MoFEMErrorCode mapFields(SmartPetscObj< DM > &intermediateDM, BitRefLevel parentBit, BitRefLevel childBit, PetscInt numVecs=1, Vec vecsToMapFrom[]=PETSC_NULLPTR, Vec vecsToMapTo[]=PETSC_NULLPTR)
SolutionMapping(MoFEM::Interface &m_field)
boost::shared_ptr< DomainEle > feRhs
MoFEMErrorCode reSetUp(std::vector< std::string > L2Fields, std::vector< int > L2Oders, std::vector< std::string > H1Fields, std::vector< int > H1Orders, std::vector< std::string > HdivFields, std::vector< int > HdivOrders, BitRefLevel bIt)
MoFEMErrorCode setUp(SmartPetscObj< DM > &subDM, SmartPetscObj< DM > &intermediateDM, BitRefLevel parent_bit, BitRefLevel child_bit)
MoFEMErrorCode projectResults(BitRefLevel parent_bit, BitRefLevel child_bit)
MoFEM::Interface & m_field
boost::shared_ptr< DomainEle > feLhs
SmartPetscObj< DM > subDM
BitRefLevel parentBit
@ FLIPPED_BIT
PetscBool is_distributed_mesh
int num_refinement_levels
constexpr auto size_symm
Definition plastic.cpp:42
int geom_order
Order if fixed.
Definition plastic.cpp:141