v0.15.5
Loading...
Searching...
No Matches
FormsBrokenSpaceConstraintImpl.hpp
Go to the documentation of this file.
1/**
2 * @file FormsBrokenSpaceConstraintImpl.hpp
3 * @brief Integrator for broken space constraints
4 * @date 2024-07-01
5 *
6 * @copyright Copyright (c) 2024
7 *
8 */
9
10#ifndef __FORMSBROKENSPACECONSTRAINTIMPL_HPP__
11#define __FORMSBROKENSPACECONSTRAINTIMPL_HPP__
12
13namespace MoFEM {
14
15/**
16 * @brief Operator for broken loop side
17 *
18 * That is used follow pattern,
19 *
20 * First: Iterate over skeleton FEs adjacent to Domain FEs
21 * Note: BoundaryEle, i.e. uses skeleton interation rule
22 *
23 * ans then
24 *
25 * Iterate over domain FEs adjacent to skelton, particularly one
26 * domain element.
27 *
28 * You use this to integrate over skelton, but field which is on adjacent domain
29 * element.
30 *
31 * @tparam E
32 */
33template <typename E> struct OpBrokenLoopSide : public OpLoopSide<E> {
34
36 using OpLoopSide<E>::OpLoopSide;
37
38 MoFEMErrorCode doWork(int side, EntityType type,
41
42 auto prev_side_fe_ptr = OP::getSidePtrFE();
43 if (OP::sideFEName == prev_side_fe_ptr->getFEName()) {
44 auto prev_fe_uid =
45 prev_side_fe_ptr->numeredEntFiniteElementPtr->getFEUId();
47 CHKERR OP::sideFEPtr->setSideFEPtr(OP::getPtrFE());
48 CHKERR OP::sideFEPtr->copyBasicMethod(*OP::getFEMethod());
49 CHKERR OP::sideFEPtr->copyPetscData(*OP::getFEMethod());
53 OP::sideFEPtr->cacheWeakPtr = prev_side_fe_ptr->cacheWeakPtr;
54 OP::sideFEPtr->loopSize = 1;
55 CHKERR OP::sideFEPtr->preProcess();
56 OP::sideFEPtr->nInTheLoop = 0;
57 OP::sideFEPtr->numeredEntFiniteElementPtr =
58 prev_side_fe_ptr->numeredEntFiniteElementPtr;
60 CHKERR OP::sideFEPtr->postProcess();
61 } else {
62 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
63 "sideFEName is different");
64 }
65
67 };
68};
69
72 inline auto &getSense() { return eleSense; }
73 inline auto &getSide() { return eleSide; }
74 inline auto &getType() { return eleType; }
75 inline auto &getData() { return entData; }
76 inline auto &getFlux() { return fluxMat; }
77 inline auto &getVarFlux() { return fluxVarMat; }
78
79private:
80 int eleSense = 0;
81 int eleSide = 1;
82 EntityType eleType = MBENTITYSET;
86};
87
88template <typename OpBase> struct OpGetBrokenBaseSideData : public OpBase {
89
90 using OP = OpBase;
91
93 const std::string field_name,
94 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data);
95
96 MoFEMErrorCode doWork(int row_side, EntityType row_type,
98
99private:
100 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
101};
102
103template <typename OpBase>
105 const std::string field_name,
106 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data)
107 : OP(field_name, field_name, OP::OPROW),
108 brokenBaseSideData(broken_base_side_data) {}
109
110template <typename OpBase>
112OpGetBrokenBaseSideData<OpBase>::doWork(int row_side, EntityType row_type,
113 EntitiesFieldData::EntData &row_data) {
115 brokenBaseSideData->resize(OP::getLoopSize());
116
117 const auto n_in_the_loop = OP::getNinTheLoop();
118 const auto face_sense = OP::getSkeletonSense();
119
120#ifndef NDEBUG
121 if (face_sense != -1 && face_sense != 1)
122 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "face sense not set");
123#endif // NDEBUG
124
125 auto set_data = [&](auto &side_data) {
126 side_data.getSide() = row_side;
127 side_data.getType() = row_type;
128 side_data.getSense() = face_sense;
129 side_data.getData().sEnse = row_data.sEnse;
130 side_data.getData().sPace = row_data.sPace;
131 side_data.getData().bAse = row_data.bAse;
132 side_data.getData().iNdices = row_data.iNdices;
133 side_data.getData().localIndices = row_data.localIndices;
134 side_data.getData().dOfs = row_data.dOfs;
135 side_data.getData().fieldEntities = row_data.fieldEntities;
136 side_data.getData().fieldData = row_data.fieldData;
137 };
138
139 auto set_base = [&](auto &side_data) {
140 auto base = side_data.getData().getBase();
141 for (auto dd = 0; dd != BaseDerivatives::LastDerivative; ++dd) {
142 for (auto bb = 0; bb != LASTBASE; ++bb) {
143 side_data.getData().baseFunctionsAndBaseDerivatives[dd][bb].reset();
144 }
145 if (auto base_ptr = row_data.baseFunctionsAndBaseDerivatives[dd][base]) {
146 side_data.getData().baseFunctionsAndBaseDerivatives[dd][base] =
147 boost::make_shared<MatrixDouble>(*base_ptr);
148 }
149 }
150 };
151
152 set_data((*brokenBaseSideData)[n_in_the_loop]);
153 set_base((*brokenBaseSideData)[n_in_the_loop]);
154
156}
157
158template <typename OpBase> struct OpSetFlux : public OpBase {
159
160 using OP = OpBase;
161
162 OpSetFlux(
163 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
164 boost::shared_ptr<MatrixDouble> flux_ptr);
165
166 MoFEMErrorCode doWork(int row_side, EntityType row_type,
168
169private:
170 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
171 boost::shared_ptr<MatrixDouble> fluxPtr;
172};
173
174template <typename OpBase>
176 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
177 boost::shared_ptr<MatrixDouble> flux_ptr)
178 : OP(NOSPACE, OP::OPSPACE), brokenBaseSideData(broken_base_side_data),
179 fluxPtr(flux_ptr) {}
180
181template <typename OpBase>
185 auto swap_flux = [&](auto &side_data) { side_data.getFlux().swap(*fluxPtr); };
186 swap_flux((*brokenBaseSideData)[OP::getNinTheLoop()]);
188}
189template <typename OpBase> struct OpSetVarFlux : public OpBase {
190
191 using OP = OpBase;
192
194 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
195 boost::shared_ptr<MatrixDouble> flux_var_ptr);
196
197 MoFEMErrorCode doWork(int row_side, EntityType row_type,
199
200private:
201 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
202 boost::shared_ptr<MatrixDouble> fluxVarPtr;
203};
204
205template <typename OpBase>
207 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
208 boost::shared_ptr<MatrixDouble> flux_var_ptr)
209 : OP(NOSPACE, OP::OPSPACE), brokenBaseSideData(broken_base_side_data),
210 fluxVarPtr(flux_var_ptr) {}
211
212template <typename OpBase>
216 auto swap_flux = [&](auto &side_data) {
217 side_data.getVarFlux().swap(*fluxVarPtr);
218 };
219 swap_flux((*brokenBaseSideData)[OP::getNinTheLoop()]);
221}
222
223
224template <typename OpBase> struct OpBrokenBaseImpl : public OpBase {
225
226 using OP = OpBase;
227
229 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
230 boost::shared_ptr<Range> ents_ptr = nullptr)
231 : OP(NOSPACE, OP::OPSPACE), brokenBaseSideData(broken_base_side_data) {
232 OP::entsPtr = ents_ptr;
233 OP::assembleTranspose = false;
234 OP::onlyTranspose = false;
235 OP::sYmm = false;
236 }
237
239 const std::string row_field,
240 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
241 const bool assmb_transpose, const bool only_transpose,
242 boost::shared_ptr<Range> ents_ptr = nullptr)
243 : OP(row_field, row_field, OP::OPROW, ents_ptr),
244 brokenBaseSideData(broken_base_side_data) {
245 OP::entsPtr = ents_ptr;
246 OP::assembleTranspose = assmb_transpose;
247 OP::onlyTranspose = only_transpose;
248 OP::sYmm = false;
249 }
250
251 MoFEMErrorCode doWork(int row_side, EntityType row_type,
253
254protected:
255 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
256};
257
258template <typename OpBase>
260OpBrokenBaseImpl<OpBase>::doWork(int row_side, EntityType row_type,
261 EntitiesFieldData::EntData &row_data) {
263
264 if (OP::entsPtr) {
265 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
267 }
268
269#ifndef NDEBUG
270 if (!brokenBaseSideData) {
271 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not set");
272 }
273#endif // NDEBUG
274
275 auto do_work_rhs = [this](int row_side, EntityType row_type,
276 EntitiesFieldData::EntData &row_data, int sense) {
278 // get number of dofs on row
279 OP::nbRows = row_data.getIndices().size();
280 if (!OP::nbRows)
282 // get number of integration points
283 OP::nbIntegrationPts = OP::getGaussPts().size2();
284 // get row base functions
285 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
286 // resize and clear the right hand side vector
287 OP::locF.resize(OP::nbRows, false);
288 OP::locF.clear();
289 // integrate local vector
290 CHKERR this->iNtegrate(row_data);
291 // assemble local vector
292 OP::locF *= sense;
293 CHKERR this->aSsemble(row_data);
295 };
296
297 auto do_work_lhs = [this](int row_side, int col_side, EntityType row_type,
298 EntityType col_type,
300 EntitiesFieldData::EntData &col_data, int sense) {
302
303 auto check_if_assemble_transpose = [&] {
304 if (this->sYmm) {
305 if (OP::rowSide != OP::colSide || OP::rowType != OP::colType)
306 return true;
307 else
308 return false;
309 } else if (OP::assembleTranspose) {
310 return true;
311 }
312 return false;
313 };
314
315 OP::rowSide = row_side;
316 OP::rowType = row_type;
317 OP::colSide = col_side;
318 OP::colType = col_type;
319 OP::nbCols = col_data.getIndices().size();
320 OP::locMat.resize(OP::nbRows, OP::nbCols, false);
321 OP::locMat.clear();
322 CHKERR this->iNtegrate(row_data, col_data);
323 OP::locMat *= sense;
324 CHKERR this->aSsemble(row_data, col_data, check_if_assemble_transpose());
326 };
327
328 switch (OP::opType) {
329 case OP::OPROW:
330
331 OP::nbRows = row_data.getIndices().size();
332 if (!OP::nbRows)
334 OP::nbIntegrationPts = OP::getGaussPts().size2();
335 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
336
337 if (!OP::nbRows)
339
340 for (auto &bd : *brokenBaseSideData) {
341
342#ifndef NDEBUG
343 if (!bd.getData().getNSharedPtr(bd.getData().getBase())) {
344 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
345 "base functions not set");
346 }
347#endif
348
349 CHKERR do_work_lhs(
350
351 // side
352 row_side, bd.getSide(),
353
354 // type
355 row_type, bd.getType(),
356
357 // row_data
358 row_data, bd.getData(),
359
360 // sense
361 bd.getSense()
362
363 );
364 }
365
366 break;
367 case OP::OPSPACE:
368 for (auto &bd : *brokenBaseSideData) {
369 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData(),
370 bd.getSense());
371 }
372 break;
373 default:
375 (std::string("wrong op type ") +
377 .c_str());
378 }
379
381}
382
383template <int FIELD_DIM, IntegrationType I, typename OpBrokenBase>
385
386template <int FIELD_DIM, typename OpBrokenBase>
388 : public OpBrokenBase {
389
391
393 const std::string row_field,
394 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
395 boost::shared_ptr<double> beta_ptr, const bool assmb_transpose,
396 const bool only_transpose, boost::shared_ptr<Range> ents_ptr = nullptr)
397 : OP(row_field, broken_base_side_data, assmb_transpose, only_transpose,
398 ents_ptr),
399 scalarBetaPtr(beta_ptr) {}
400
402 const std::string row_field,
403 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
404 double beta, const bool assmb_transpose, const bool only_transpose,
405 boost::shared_ptr<Range> ents_ptr = nullptr)
406 : OpBrokenSpaceConstrainImpl(row_field, broken_base_side_data,
407 boost::make_shared<double>(beta),
408 assmb_transpose, only_transpose, ents_ptr) {}
409
410protected:
411 boost::shared_ptr<double> scalarBetaPtr;
414};
415
416template <int FIELD_DIM, typename OpBase>
419 EntitiesFieldData::EntData &col_data) {
421
422 auto nb_row_dofs = row_data.getIndices().size();
423 auto nb_col_dofs = col_data.getIndices().size();
424 if (!nb_row_dofs || !nb_col_dofs)
426
428 FTENSOR_INDEX(3, J);
429
430 auto t_w = this->getFTensor0IntegrationWeight();
431 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
432 size_t nb_base_functions = col_data.getN().size2() / 3;
433
434#ifndef NDEBUG
435 if (nb_row_dofs % FIELD_DIM != 0 || nb_col_dofs % FIELD_DIM != 0) {
436 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
437 "number of dofs not divisible by field dimension");
438 }
439 if (nb_row_dofs > row_data.getN().size2() * FIELD_DIM ||
440 nb_col_dofs > col_data.getN().size2() * FIELD_DIM) {
441 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
442 "number of dofs exceeds number of base functions");
443 }
444#endif // NDEBUG
445
446 auto t_col_base = col_data.getFTensor1N<3>();
447
448 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
449
450 int cc = 0;
451 for (; cc != nb_col_dofs / FIELD_DIM; cc++) {
452 auto t_row_base = row_data.getFTensor0N(gg, 0);
453 for (auto rr = 0; rr != nb_row_dofs / FIELD_DIM; ++rr) {
454 OP::locMat(FIELD_DIM * rr, FIELD_DIM * cc) +=
455 (0.5 * t_w * t_row_base) * (t_normal(J) * t_col_base(J));
456 ++t_row_base;
457 }
458
459 ++t_col_base;
460 }
461
462 for (; cc < nb_base_functions; ++cc)
463 ++t_col_base;
464
465 ++t_w;
466 ++t_normal;
467 }
468
469 for (auto rr = 0; rr != nb_row_dofs / FIELD_DIM; ++rr) {
470 for (auto cc = 0; cc != nb_col_dofs / FIELD_DIM; ++cc) {
471 for (auto dd = 1; dd < FIELD_DIM; ++dd) {
472 OP::locMat(FIELD_DIM * rr + dd, FIELD_DIM * cc + dd) =
473 OP::locMat(FIELD_DIM * rr, FIELD_DIM * cc);
474 }
475 }
476 }
477
478 if (scalarBetaPtr)
479 OP::locMat *= *scalarBetaPtr;
480
482}
483
484template <int FIELD_DIM, IntegrationType I, typename OpBase>
486
487template <int FIELD_DIM, IntegrationType I, typename OpBase>
489
490template <int FIELD_DIM, typename OpBrokenBase>
492 : public OpBrokenBase {
493
495
497 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
498 boost::shared_ptr<MatrixDouble> lagrange_ptr,
499 boost::shared_ptr<double> beta_ptr,
500 boost::shared_ptr<Range> ents_ptr = nullptr)
501 : OP(broken_base_side_data, ents_ptr), scalarBetaPtr(beta_ptr),
502 lagrangePtr(lagrange_ptr) {}
503
505 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
506 boost::shared_ptr<MatrixDouble> lagrange_ptr, double beta,
507 boost::shared_ptr<Range> ents_ptr = nullptr)
508 : OpBrokenSpaceConstrainDFluxImpl(broken_base_side_data, lagrange_ptr,
509 boost::make_shared<double>(beta),
510 ents_ptr) {}
511
512private:
514
515 boost::shared_ptr<double> scalarBetaPtr;
516 boost::shared_ptr<MatrixDouble> lagrangePtr;
517};
518
519template <int FIELD_DIM, typename OpBase>
522 EntitiesFieldData::EntData &row_data) {
524
526 FTENSOR_INDEX(3, J);
527
528 auto t_w = this->getFTensor0IntegrationWeight();
529 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
530 auto t_lagrange = getFTensor1FromMat<FIELD_DIM>(*lagrangePtr);
531
532 auto t_row_base = row_data.getFTensor1N<3>();
533 auto nb_base_functions = row_data.getN().size2() / 3;
534
535 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
536 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
537 size_t rr = 0;
538 for (; rr != row_data.getIndices().size() / FIELD_DIM; ++rr) {
539 t_vec(i) += (0.5 * t_w * (t_row_base(J) * t_normal(J))) * t_lagrange(i);
540 ++t_row_base;
541 ++t_vec;
542 }
543 for (; rr < nb_base_functions; ++rr)
544 ++t_row_base;
545 ++t_w;
546 ++t_normal;
547 ++t_lagrange;
548 }
549
550 if(scalarBetaPtr)
551 OP::locF *= *scalarBetaPtr;
552
554}
555
556template <int FIELD_DIM, typename OpBase>
558 : public OpBase {
559
560 using OP = OpBase;
561
563 const std::string row_field,
564 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
565 boost::shared_ptr<double> beta_ptr,
566 boost::shared_ptr<Range> ents_ptr = nullptr)
567 : OpBase(row_field, row_field, OpBase::OPROW, ents_ptr),
568 brokenSideDataPtr(broken_side_data_ptr), scalarBetaPtr(beta_ptr) {}
569
571 const std::string row_field,
572 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
573 double beta, boost::shared_ptr<Range> ents_ptr = nullptr)
574 : OpBrokenSpaceConstrainDHybridImpl(row_field, broken_side_data_ptr,
575 boost::make_shared<double>(beta),
576 ents_ptr) {}
577
578private:
579 boost::shared_ptr<double> scalarBetaPtr;
580 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenSideDataPtr;
581 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data);
582};
583
584template <int FIELD_DIM, typename OpBase>
587 EntitiesFieldData::EntData &row_data) {
589
591 FTENSOR_INDEX(3, J);
592
593 OP::locF.resize(row_data.getIndices().size(), false);
594 OP::locF.clear();
595
596 for (auto &bd : *brokenSideDataPtr) {
597 auto t_w = this->getFTensor0IntegrationWeight();
598 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
599 auto t_row_base = row_data.getFTensor0N();
600 auto t_flux = getFTensor2FromMat<FIELD_DIM, 3>(bd.getFlux());
601 auto nb_base_functions = row_data.getN().size2() / 3;
602 auto sense = bd.getSense();
603 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
604 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
605 size_t rr = 0;
606 for (; rr != row_data.getIndices().size() / FIELD_DIM; ++rr) {
607 t_vec(i) +=
608 (sense * 0.5 * t_w) * t_row_base * t_normal(J) * t_flux(i, J);
609 ++t_row_base;
610 ++t_vec;
611 }
612 for (; rr < nb_base_functions; ++rr)
613 ++t_row_base;
614 ++t_w;
615 ++t_normal;
616 ++t_flux;
617 }
618 }
619
620 if(scalarBetaPtr)
621 OP::locF *= *scalarBetaPtr;
622
624}
625
626template <int FIELD_DIM, IntegrationType I, typename OpBase>
628
629template <int FIELD_DIM, IntegrationType I, typename OpBase>
631
632template <typename OpBrokenBase> struct OpBrokenTopoBase : public OpBrokenBase {
634
635protected:
636
638 const std::string row_field,
639 boost::shared_ptr<MatrixDouble> tangent1_diff_ptr,
640 boost::shared_ptr<MatrixDouble> tangent2_diff_ptr,
641 SmartPetscObj<Vec> assemble_vec, Tag th,
642 boost::shared_ptr<Range> ents_ptr = nullptr)
643 : OP(row_field, row_field, OP::OPROW, ents_ptr),
644 tangent1DiffPtr(tangent1_diff_ptr),
645 tangent2DiffPtr(tangent2_diff_ptr), assembleVec(assemble_vec),
646 thGradTag(th) {}
647
650 if (!this->timeScalingFun.empty())
651 this->locF *= this->timeScalingFun(this->getFEMethod()->ts_t);
652 if (!this->feScalingFun.empty())
653 this->locF *= this->feScalingFun(this->getFEMethod());
654 if (assembleVec) {
655 auto *vec_ptr = this->locF.data().data();
656 const auto nb_dofs = row_data.getIndices().size();
657 auto *ind_ptr = row_data.getIndices().data().data();
658 return VecSetValues(assembleVec, nb_dofs, ind_ptr, vec_ptr, ADD_VALUES);
659 }
660 if (thGradTag) {
661 const auto field_ents = row_data.getFieldEntities();
662 std::vector<EntityHandle> ents(field_ents.size());
663 std::transform(field_ents.begin(), field_ents.end(), ents.begin(),
664 [](const auto *fe) { return fe->getEnt(); });
665 if (field_ents.empty())
667 if (type_from_handle(ents[0]) != MBVERTEX)
669 auto &moab = this->getMoab();
670 VectorDouble topo_values(this->locF.size());
671 CHKERR moab.tag_get_data(thGradTag, ents.data(), ents.size(),
672 topo_values.data().data());
673 topo_values += this->locF;
674 CHKERR moab.tag_set_data(thGradTag, ents.data(), ents.size(),
675 topo_values.data().data());
676 }
678 }
679
680 boost::shared_ptr<MatrixDouble> tangent1DiffPtr;
681 boost::shared_ptr<MatrixDouble> tangent2DiffPtr;
684};
685
686template <typename OpBase>
688 : public OpBrokenTopoBase<OpBase> {
689
691
693 const std::string row_field,
694 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
695 boost::shared_ptr<MatrixDouble> lamgrange_ptr,
696 boost::shared_ptr<MatrixDouble> tangent1_ptr,
697 boost::shared_ptr<MatrixDouble> tangent2_ptr,
698 boost::shared_ptr<double> beta_ptr, SmartPetscObj<Vec> assemble_vec,
699 Tag th, boost::shared_ptr<Range> ents_ptr = nullptr)
700 : OP(row_field, tangent1_ptr, tangent2_ptr, assemble_vec, th, ents_ptr),
701 brokenSideDataPtr(broken_side_data_ptr), lagrangePtr(lamgrange_ptr),
702 scalarBetaPtr(beta_ptr) {}
703
704private:
705 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data);
706 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenSideDataPtr;
707 boost::shared_ptr<MatrixDouble> lagrangePtr;
708 boost::shared_ptr<double> scalarBetaPtr;
709};
710
711template <typename OpBase>
714 EntitiesFieldData::EntData &row_data) {
716
717 FTENSOR_INDEX(3, I);
718 FTENSOR_INDEX(3, J);
719 FTENSOR_INDEX(3, i);
720 FTENSOR_INDEX(3, j);
721 FTENSOR_INDEX(3, k);
722 FTENSOR_INDEX(3, m);
725
726 auto get_ftensor1 = [](MatrixDouble &m) {
728 &m(0, 0), &m(0, 1), &m(0, 2));
729 };
730
731 for (auto &bd : *brokenSideDataPtr) {
732
733 auto t_w = this->getFTensor0IntegrationWeight();
734 auto t_lagrange = getFTensor1FromMat<3>(*lagrangePtr);
735 auto t_adjoint_lambda = getFTensor2FromMat<3, 3>(bd.getVarFlux());
736 auto t_t1 = get_ftensor1(*this->tangent1DiffPtr);
737 auto t_t2 = get_ftensor1(*this->tangent2DiffPtr);
738
739 auto t_diff_base = row_data.getFTensor1DiffN<2>();
740 auto nb_base_functions = row_data.getN().size2();
741
742 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
743
744 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
745
746 FTensor::Tensor2<double, 3, 3> t_normal_diff;
747 t_normal_diff(j, m) = FTensor::levi_civita(i, j, k) * t_t1(k) *
748 (t_kd(i, m) * t_diff_base(N1));
749 t_normal_diff(j, m) += FTensor::levi_civita(i, j, k) *
750 (t_kd(k, m) * t_diff_base(N0)) * t_t2(i);
751 auto t_vec = getFTensor1FromPtr<3>(&*OP::locF.data().begin());
752 size_t rr = 0;
753 for (; rr != row_data.getIndices().size() / 3; ++rr) {
754 t_vec(m) += (0.5 * t_w) * (t_adjoint_lambda(i, J) * t_lagrange(i)) *
755 t_normal_diff(J, m);
756 ++t_diff_base;
757 ++t_vec;
758 }
759 for (; rr < nb_base_functions; ++rr)
760 ++t_diff_base;
761 ++t_w;
762 ++t_lagrange;
763 ++t_adjoint_lambda;
764 ++t_t1;
765 ++t_t2;
766 }
767 }
768
769 if (scalarBetaPtr)
770 OP::locF *= *scalarBetaPtr;
771
773}
774
775template <typename OpBase>
777 : public OpBrokenTopoBase<OpBase> {
778
780
782 const std::string row_field,
783 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
784 boost::shared_ptr<MatrixDouble> adjoint_hybrid_ptr,
785 boost::shared_ptr<MatrixDouble> tangent1_ptr,
786 boost::shared_ptr<MatrixDouble> tangent2_ptr,
787 boost::shared_ptr<double> beta_ptr, SmartPetscObj<Vec> assemble_vec,
788 Tag th, boost::shared_ptr<Range> ents_ptr = nullptr)
789 : OP(row_field, tangent1_ptr, tangent2_ptr, assemble_vec, th, ents_ptr),
790 brokenSideDataPtr(broken_side_data_ptr),
791 adjointHybridPtr(adjoint_hybrid_ptr), scalarBetaPtr(beta_ptr) {}
792
793private:
794 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenSideDataPtr;
795 boost::shared_ptr<MatrixDouble> adjointHybridPtr;
796 boost::shared_ptr<double> scalarBetaPtr;
797 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data);
798};
799
800template <typename OpBase>
802 3, GAUSS, OpBase>::iNtegrate(EntitiesFieldData::EntData &row_data) {
804
805 FTENSOR_INDEX(3, i);
806 FTENSOR_INDEX(3, j);
807 FTENSOR_INDEX(3, k);
808 FTENSOR_INDEX(3, m);
809 FTENSOR_INDEX(3, J);
810 FTENSOR_INDEX(3, K);
811 FTENSOR_INDEX(3, L);
812 FTENSOR_INDEX(3, M);
813
816
817 OP::locF.resize(row_data.getIndices().size(), false);
818 OP::locF.clear();
819
820 auto get_ftensor1 = [](MatrixDouble &m) {
822 &m(0, 0), &m(0, 1), &m(0, 2));
823 };
824
825 for (auto &bd : *brokenSideDataPtr) {
826 auto t_w = this->getFTensor0IntegrationWeight();
827 auto t_t1 = get_ftensor1(*this->tangent1DiffPtr);
828 auto t_t2 = get_ftensor1(*this->tangent2DiffPtr);
829 auto t_flux = getFTensor2FromMat<3, 3>(bd.getFlux());
830 auto t_adjoint_lambda = getFTensor1FromMat<3>(*adjointHybridPtr);
831
832 auto nb_base_functions = row_data.getN().size2();
833 auto sense = bd.getSense();
834 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
835
836 auto t_diff_base = row_data.getFTensor1DiffN<2>();
837 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
838 auto t_vec = getFTensor1FromPtr<3>(&*OP::locF.data().begin());
839 size_t rr = 0;
840 for (; rr != row_data.getIndices().size() / 3; ++rr) {
841 FTensor::Tensor2<double, 3, 3> t_normal_diff;
842 t_normal_diff(j, m) =
843 FTensor::levi_civita(i, j, k) * t_t1(k) *
844 (t_kd(i, m) * t_diff_base(N1));
845 t_normal_diff(j, m) +=
846 FTensor::levi_civita(i, j, k) * (t_kd(k, m) * t_diff_base(N0)) *
847 t_t2(i);
848
849 t_vec(m) += (sense * 0.5 * t_w) * (t_adjoint_lambda(i) * t_flux(i, J)) *
850 t_normal_diff(J, m);
851
852 ++t_diff_base;
853 ++t_vec;
854 }
855 for (; rr < nb_base_functions; ++rr)
856 ++t_diff_base;
857 ++t_w;
858 ++t_t1;
859 ++t_t2;
860 ++t_flux;
861 ++t_adjoint_lambda;
862 }
863 }
864
865 if (scalarBetaPtr)
866 OP::locF *= *scalarBetaPtr;
867
869}
870
871} // namespace MoFEM
872
873#endif // __FORMSBROKENSPACECONSTRAINTIMPL_HPP__
#define FTENSOR_INDEX(DIM, I)
constexpr int FIELD_DIM
Kronecker Delta class.
@ LASTBASE
Definition definitions.h:69
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ NOSPACE
Definition definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
@ MOFEM_IMPOSSIBLE_CASE
Definition definitions.h:35
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr auto t_kd
@ GAUSS
Gaussian quadrature integration.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition ddTensor0.hpp:33
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto type_from_handle(const EntityHandle h)
get type from entity handle
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
constexpr IntegrationType I
constexpr auto field_name
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Definition radiation.cpp:29
FTensor::Index< 'm', 3 > m
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
int sEnse
Entity sense (orientation)
VectorDouble fieldData
Field data on entity.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
VectorInt localIndices
Local indices on entity.
const VectorFieldEntities & getFieldEntities() const
Get field entities (const version)
VectorInt iNdices
Global indices on entity.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
virtual int getSense() const
Get entity sense for conforming approximation fields.
VectorFieldEntities fieldEntities
Field entities.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
std::array< std::array< boost::shared_ptr< MatrixDouble >, LASTBASE >, LastDerivative > baseFunctionsAndBaseDerivatives
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
FieldApproximationBase bAse
Field approximation base.
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
boost::shared_ptr< Range > entsPtr
Entities on which element is run.
int nbIntegrationPts
number of integration points
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideData
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
OpBrokenBaseImpl(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, const bool assmb_transpose, const bool only_transpose, boost::shared_ptr< Range > ents_ptr=nullptr)
OpBrokenBaseImpl(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< Range > ents_ptr=nullptr)
Operator for broken loop side.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data)
OpBrokenSpaceConstrainDFluxImpl(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< MatrixDouble > lagrange_ptr, boost::shared_ptr< double > beta_ptr, boost::shared_ptr< Range > ents_ptr=nullptr)
OpBrokenSpaceConstrainDFluxImpl(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< MatrixDouble > lagrange_ptr, double beta, boost::shared_ptr< Range > ents_ptr=nullptr)
OpBrokenSpaceConstrainDHybridImpl(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_side_data_ptr, boost::shared_ptr< double > beta_ptr, boost::shared_ptr< Range > ents_ptr=nullptr)
OpBrokenSpaceConstrainDHybridImpl(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_side_data_ptr, double beta, boost::shared_ptr< Range > ents_ptr=nullptr)
OpBrokenSpaceConstrainImpl(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< double > beta_ptr, const bool assmb_transpose, const bool only_transpose, boost::shared_ptr< Range > ents_ptr=nullptr)
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
OpBrokenSpaceConstrainImpl(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, double beta, const bool assmb_transpose, const bool only_transpose, boost::shared_ptr< Range > ents_ptr=nullptr)
boost::shared_ptr< MatrixDouble > tangent2DiffPtr
boost::shared_ptr< MatrixDouble > tangent1DiffPtr
MoFEMErrorCode aSsemble(EntitiesFieldData::EntData &row_data)
OpBrokenTopoBase(const std::string row_field, boost::shared_ptr< MatrixDouble > tangent1_diff_ptr, boost::shared_ptr< MatrixDouble > tangent2_diff_ptr, SmartPetscObj< Vec > assemble_vec, Tag th, boost::shared_ptr< Range > ents_ptr=nullptr)
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideData
OpGetBrokenBaseSideData(const std::string field_name, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data)
Element used to execute operators on side of the element.
const std::string sideFEName
boost::shared_ptr< E > sideFEPtr
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
OpSetFlux(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< MatrixDouble > flux_ptr)
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideData
boost::shared_ptr< MatrixDouble > fluxPtr
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
boost::shared_ptr< MatrixDouble > fluxVarPtr
OpSetVarFlux(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< MatrixDouble > flux_var_ptr)
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideData
OpTopoDerivativeBrokenSpaceConstrainDFluxImpl(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_side_data_ptr, boost::shared_ptr< MatrixDouble > lamgrange_ptr, boost::shared_ptr< MatrixDouble > tangent1_ptr, boost::shared_ptr< MatrixDouble > tangent2_ptr, boost::shared_ptr< double > beta_ptr, SmartPetscObj< Vec > assemble_vec, Tag th, boost::shared_ptr< Range > ents_ptr=nullptr)
OpTopoDerivativeBrokenSpaceConstrainDHybridImpl(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_side_data_ptr, boost::shared_ptr< MatrixDouble > adjoint_hybrid_ptr, boost::shared_ptr< MatrixDouble > tangent1_ptr, boost::shared_ptr< MatrixDouble > tangent2_ptr, boost::shared_ptr< double > beta_ptr, SmartPetscObj< Vec > assemble_vec, Tag th, boost::shared_ptr< Range > ents_ptr=nullptr)
intrusive_ptr for managing petsc objects