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
78private:
79 int eleSense = 0;
80 int eleSide = 1;
81 EntityType eleType = MBENTITYSET;
84};
85
86template <typename OpBase> struct OpGetBrokenBaseSideData : public OpBase {
87
88 using OP = OpBase;
89
91 const std::string field_name,
92 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data);
93
94 MoFEMErrorCode doWork(int row_side, EntityType row_type,
96
97private:
98 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
99};
100
101template <typename OpBase>
103 const std::string field_name,
104 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data)
105 : OP(field_name, field_name, OP::OPROW),
106 brokenBaseSideData(broken_base_side_data) {}
107
108template <typename OpBase>
110OpGetBrokenBaseSideData<OpBase>::doWork(int row_side, EntityType row_type,
111 EntitiesFieldData::EntData &row_data) {
113 brokenBaseSideData->resize(OP::getLoopSize());
114
115 const auto n_in_the_loop = OP::getNinTheLoop();
116 const auto face_sense = OP::getSkeletonSense();
117
118#ifndef NDEBUG
119 if (face_sense != -1 && face_sense != 1)
120 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "face sense not set");
121#endif // NDEBUG
122
123 auto set_data = [&](auto &side_data) {
124 side_data.getSide() = row_side;
125 side_data.getType() = row_type;
126 side_data.getSense() = face_sense;
127 side_data.getData().sEnse = row_data.sEnse;
128 side_data.getData().sPace = row_data.sPace;
129 side_data.getData().bAse = row_data.bAse;
130 side_data.getData().iNdices = row_data.iNdices;
131 side_data.getData().localIndices = row_data.localIndices;
132 side_data.getData().dOfs = row_data.dOfs;
133 side_data.getData().fieldEntities = row_data.fieldEntities;
134 side_data.getData().fieldData = row_data.fieldData;
135 };
136
137 auto set_base = [&](auto &side_data) {
138 auto base = side_data.getData().getBase();
139 for (auto dd = 0; dd != BaseDerivatives::LastDerivative; ++dd) {
140 for (auto bb = 0; bb != LASTBASE; ++bb) {
141 side_data.getData().baseFunctionsAndBaseDerivatives[dd][bb].reset();
142 }
143 if (auto base_ptr = row_data.baseFunctionsAndBaseDerivatives[dd][base]) {
144 side_data.getData().baseFunctionsAndBaseDerivatives[dd][base] =
145 boost::make_shared<MatrixDouble>(*base_ptr);
146 }
147 }
148 };
149
150 set_data((*brokenBaseSideData)[n_in_the_loop]);
151 set_base((*brokenBaseSideData)[n_in_the_loop]);
152
154}
155
156template <typename OpBase> struct OpSetFlux : public OpBase {
157
158 using OP = OpBase;
159
160 OpSetFlux(
161 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
162 boost::shared_ptr<MatrixDouble> flux_ptr);
163
164 MoFEMErrorCode doWork(int row_side, EntityType row_type,
166
167private:
168 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
169 boost::shared_ptr<MatrixDouble> fluxPtr;
170};
171
172template <typename OpBase>
174 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
175 boost::shared_ptr<MatrixDouble> flux_ptr)
176 : OP(NOSPACE, OP::OPSPACE), brokenBaseSideData(broken_base_side_data),
177 fluxPtr(flux_ptr) {}
178
179template <typename OpBase>
183 auto swap_flux = [&](auto &side_data) { side_data.getFlux().swap(*fluxPtr); };
184 swap_flux((*brokenBaseSideData)[OP::getNinTheLoop()]);
186}
187
188template <typename OpBase> struct OpBrokenBaseImpl : public OpBase {
189
190 using OP = OpBase;
191
193 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
194 boost::shared_ptr<Range> ents_ptr = nullptr)
195 : OP(NOSPACE, OP::OPSPACE), brokenBaseSideData(broken_base_side_data) {
196 OP::entsPtr = ents_ptr;
197 }
198
200 const std::string row_field,
201 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
202 const bool assmb_transpose, const bool only_transpose,
203 boost::shared_ptr<Range> ents_ptr = nullptr)
204 : OP(row_field, row_field, OP::OPROW, ents_ptr),
205 brokenBaseSideData(broken_base_side_data) {
206 OP::entsPtr = ents_ptr;
207 OP::assembleTranspose = assmb_transpose;
208 OP::onlyTranspose = only_transpose;
209 OP::sYmm = false;
210 }
211
212 MoFEMErrorCode doWork(int row_side, EntityType row_type,
214
215protected:
216 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
217};
218
219template <typename OpBase>
221OpBrokenBaseImpl<OpBase>::doWork(int row_side, EntityType row_type,
222 EntitiesFieldData::EntData &row_data) {
224
225 if (OP::entsPtr) {
226 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
228 }
229
230#ifndef NDEBUG
231 if (!brokenBaseSideData) {
232 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not set");
233 }
234#endif // NDEBUG
235
236 auto do_work_rhs = [this](int row_side, EntityType row_type,
237 EntitiesFieldData::EntData &row_data, int sense) {
239 // get number of dofs on row
240 OP::nbRows = row_data.getIndices().size();
241 if (!OP::nbRows)
243 // get number of integration points
244 OP::nbIntegrationPts = OP::getGaussPts().size2();
245 // get row base functions
246 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
247 // resize and clear the right hand side vector
248 OP::locF.resize(OP::nbRows, false);
249 OP::locF.clear();
250 // integrate local vector
251 CHKERR this->iNtegrate(row_data);
252 // assemble local vector
253 OP::locF *= sense;
254 CHKERR this->aSsemble(row_data);
256 };
257
258 auto do_work_lhs = [this](int row_side, int col_side, EntityType row_type,
259 EntityType col_type,
261 EntitiesFieldData::EntData &col_data, int sense) {
263
264 auto check_if_assemble_transpose = [&] {
265 if (this->sYmm) {
266 if (OP::rowSide != OP::colSide || OP::rowType != OP::colType)
267 return true;
268 else
269 return false;
270 } else if (OP::assembleTranspose) {
271 return true;
272 }
273 return false;
274 };
275
276 OP::rowSide = row_side;
277 OP::rowType = row_type;
278 OP::colSide = col_side;
279 OP::colType = col_type;
280 OP::nbCols = col_data.getIndices().size();
281 OP::locMat.resize(OP::nbRows, OP::nbCols, false);
282 OP::locMat.clear();
283 CHKERR this->iNtegrate(row_data, col_data);
284 OP::locMat *= sense;
285 CHKERR this->aSsemble(row_data, col_data, check_if_assemble_transpose());
287 };
288
289 switch (OP::opType) {
290 case OP::OPROW:
291
292 OP::nbRows = row_data.getIndices().size();
293 if (!OP::nbRows)
295 OP::nbIntegrationPts = OP::getGaussPts().size2();
296 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
297
298 if (!OP::nbRows)
300
301 for (auto &bd : *brokenBaseSideData) {
302
303#ifndef NDEBUG
304 if (!bd.getData().getNSharedPtr(bd.getData().getBase())) {
305 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
306 "base functions not set");
307 }
308#endif
309
310 CHKERR do_work_lhs(
311
312 // side
313 row_side, bd.getSide(),
314
315 // type
316 row_type, bd.getType(),
317
318 // row_data
319 row_data, bd.getData(),
320
321 // sense
322 bd.getSense()
323
324 );
325 }
326
327 break;
328 case OP::OPSPACE:
329 for (auto &bd : *brokenBaseSideData) {
330 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData(),
331 bd.getSense());
332 }
333 break;
334 default:
336 (std::string("wrong op type ") +
338 .c_str());
339 }
340
342}
343
344template <int FIELD_DIM, IntegrationType I, typename OpBrokenBase>
346
347template <int FIELD_DIM, typename OpBrokenBase>
349 : public OpBrokenBase {
350
352
354 const std::string row_field,
355 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
356 boost::shared_ptr<double> beta_ptr, const bool assmb_transpose,
357 const bool only_transpose, boost::shared_ptr<Range> ents_ptr = nullptr)
358 : OP(row_field, broken_base_side_data, assmb_transpose, only_transpose,
359 ents_ptr),
360 scalarBetaPtr(beta_ptr) {}
361
363 const std::string row_field,
364 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
365 double beta, const bool assmb_transpose, const bool only_transpose,
366 boost::shared_ptr<Range> ents_ptr = nullptr)
367 : OpBrokenSpaceConstrainImpl(row_field, broken_base_side_data,
368 boost::make_shared<double>(beta),
369 assmb_transpose, only_transpose, ents_ptr) {}
370
371protected:
372 boost::shared_ptr<double> scalarBetaPtr;
375};
376
377template <int FIELD_DIM, typename OpBase>
380 EntitiesFieldData::EntData &col_data) {
382
383 auto nb_row_dofs = row_data.getIndices().size();
384 auto nb_col_dofs = col_data.getIndices().size();
385 if (!nb_row_dofs || !nb_col_dofs)
387
389 FTENSOR_INDEX(3, J);
390
391 auto t_w = this->getFTensor0IntegrationWeight();
392 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
393 size_t nb_base_functions = col_data.getN().size2() / 3;
394
395#ifndef NDEBUG
396 if (nb_row_dofs % FIELD_DIM != 0 || nb_col_dofs % FIELD_DIM != 0) {
397 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
398 "number of dofs not divisible by field dimension");
399 }
400 if (nb_row_dofs > row_data.getN().size2() * FIELD_DIM ||
401 nb_col_dofs > col_data.getN().size2() * FIELD_DIM) {
402 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
403 "number of dofs exceeds number of base functions");
404 }
405#endif // NDEBUG
406
407 auto t_col_base = col_data.getFTensor1N<3>();
408
409 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
410
411 int cc = 0;
412 for (; cc != nb_col_dofs / FIELD_DIM; cc++) {
413 auto t_row_base = row_data.getFTensor0N(gg, 0);
414 for (auto rr = 0; rr != nb_row_dofs / FIELD_DIM; ++rr) {
415 OP::locMat(FIELD_DIM * rr, FIELD_DIM * cc) +=
416 (0.5 * t_w * t_row_base) * (t_normal(J) * t_col_base(J));
417 ++t_row_base;
418 }
419
420 ++t_col_base;
421 }
422
423 for (; cc < nb_base_functions; ++cc)
424 ++t_col_base;
425
426 ++t_w;
427 ++t_normal;
428 }
429
430 for (auto rr = 0; rr != nb_row_dofs / FIELD_DIM; ++rr) {
431 for (auto cc = 0; cc != nb_col_dofs / FIELD_DIM; ++cc) {
432 for (auto dd = 1; dd < FIELD_DIM; ++dd) {
433 OP::locMat(FIELD_DIM * rr + dd, FIELD_DIM * cc + dd) =
434 OP::locMat(FIELD_DIM * rr, FIELD_DIM * cc);
435 }
436 }
437 }
438
439 if (scalarBetaPtr)
440 OP::locMat *= *scalarBetaPtr;
441
443}
444
445template <int FIELD_DIM, IntegrationType I, typename OpBase>
447
448template <int FIELD_DIM, IntegrationType I, typename OpBase>
450
451template <int FIELD_DIM, typename OpBrokenBase>
453 : public OpBrokenBase {
454
456
458 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
459 boost::shared_ptr<MatrixDouble> lagrange_ptr,
460 boost::shared_ptr<double> beta_ptr,
461 boost::shared_ptr<Range> ents_ptr = nullptr)
462 : OP(broken_base_side_data, ents_ptr), scalarBetaPtr(beta_ptr),
463 lagrangePtr(lagrange_ptr) {}
464
466 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
467 boost::shared_ptr<MatrixDouble> lagrange_ptr, double beta,
468 boost::shared_ptr<Range> ents_ptr = nullptr)
469 : OpBrokenSpaceConstrainDFluxImpl(broken_base_side_data, lagrange_ptr,
470 boost::make_shared<double>(beta),
471 ents_ptr) {}
472
473private:
475
476 boost::shared_ptr<double> scalarBetaPtr;
477 boost::shared_ptr<MatrixDouble> lagrangePtr;
478};
479
480template <int FIELD_DIM, typename OpBase>
483 EntitiesFieldData::EntData &row_data) {
485
487 FTENSOR_INDEX(3, J);
488
489 auto t_w = this->getFTensor0IntegrationWeight();
490 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
491 auto t_lagrange = getFTensor1FromMat<FIELD_DIM>(*lagrangePtr);
492
493 auto t_row_base = row_data.getFTensor1N<3>();
494 auto nb_base_functions = row_data.getN().size2() / 3;
495
496 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
497 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
498 size_t rr = 0;
499 for (; rr != row_data.getIndices().size() / FIELD_DIM; ++rr) {
500 t_vec(i) += (0.5 * t_w * (t_row_base(J) * t_normal(J))) * t_lagrange(i);
501 ++t_row_base;
502 ++t_vec;
503 }
504 for (; rr < nb_base_functions; ++rr)
505 ++t_row_base;
506 ++t_w;
507 ++t_normal;
508 ++t_lagrange;
509 }
510
511 if(scalarBetaPtr)
512 OP::locF *= *scalarBetaPtr;
513
515}
516
517template <int FIELD_DIM, typename OpBase>
519 : public OpBase {
520
521 using OP = OpBase;
522
524 const std::string row_field,
525 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
526 boost::shared_ptr<double> beta_ptr,
527 boost::shared_ptr<Range> ents_ptr = nullptr)
528 : OpBase(row_field, row_field, OpBase::OPROW, ents_ptr),
529 brokenSideDataPtr(broken_side_data_ptr), scalarBetaPtr(beta_ptr) {}
530
532 const std::string row_field,
533 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
534 double beta, boost::shared_ptr<Range> ents_ptr = nullptr)
535 : OpBrokenSpaceConstrainDHybridImpl(row_field, broken_side_data_ptr,
536 boost::make_shared<double>(beta),
537 ents_ptr) {}
538
539private:
540 boost::shared_ptr<double> scalarBetaPtr;
541 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenSideDataPtr;
543};
544
545template <int FIELD_DIM, typename OpBase>
548 EntitiesFieldData::EntData &row_data) {
550
552 FTENSOR_INDEX(3, J);
553
554 OP::locF.resize(row_data.getIndices().size(), false);
555 OP::locF.clear();
556
557 for (auto &bd : *brokenSideDataPtr) {
558 auto t_w = this->getFTensor0IntegrationWeight();
559 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
560 auto t_row_base = row_data.getFTensor0N();
561 auto t_flux = getFTensor2FromMat<FIELD_DIM, 3>(bd.getFlux());
562 auto nb_base_functions = row_data.getN().size2() / 3;
563 auto sense = bd.getSense();
564 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
565 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
566 size_t rr = 0;
567 for (; rr != row_data.getIndices().size() / FIELD_DIM; ++rr) {
568 t_vec(i) +=
569 (sense * 0.5 * t_w) * t_row_base * t_normal(J) * t_flux(i, J);
570 ++t_row_base;
571 ++t_vec;
572 }
573 for (; rr < nb_base_functions; ++rr)
574 ++t_row_base;
575 ++t_w;
576 ++t_normal;
577 ++t_flux;
578 }
579 }
580
581 if(scalarBetaPtr)
582 OP::locF *= *scalarBetaPtr;
583
585}
586
587template <int FIELD_DIM, IntegrationType I, typename OpBase>
589
590template <int FIELD_DIM, IntegrationType I, typename OpBase>
592
593template <typename OpBrokenBase> struct OpBrokenTopoBase : public OpBrokenBase {
595
597 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
598 boost::shared_ptr<MatrixDouble> adjoint_lambda_ptr,
599 boost::shared_ptr<MatrixDouble> tangent1_diff_ptr,
600 boost::shared_ptr<MatrixDouble> tangent2_diff_ptr,
601 boost::shared_ptr<double> beta_ptr, SmartPetscObj<Vec> assemble_vec,
602 Tag th, boost::shared_ptr<Range> ents_ptr = nullptr)
603 : OP(broken_base_side_data, ents_ptr),
604 adjointLambdaPtr(adjoint_lambda_ptr),
605 tangent1DiffPtr(tangent1_diff_ptr), tangent2DiffPtr(tangent2_diff_ptr),
606 assembleVec(assemble_vec), thGradTag(th) {}
607
609 const std::string row_field,
610 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
611 boost::shared_ptr<MatrixDouble> adjoint_lambda_ptr,
612 boost::shared_ptr<MatrixDouble> tangent1_ptr,
613 boost::shared_ptr<MatrixDouble> tangent2_ptr,
614 SmartPetscObj<Vec> assemble_vec, Tag th,
615 boost::shared_ptr<Range> ents_ptr = nullptr)
616 : OpBase(row_field, row_field, OpBase::OPROW, ents_ptr),
617 brokenSideDataPtr(broken_side_data_ptr),
618 adjointLambdaPtr(adjoint_lambda_ptr), tangent1Ptr(tangent1_ptr),
619 tangent2Ptr(tangent2_ptr), assembleVec(assemble_vec), thGradTag(th) {}
620
621protected:
623 if (!this->timeScalingFun.empty())
624 this->locF *= this->timeScalingFun(this->getFEMethod()->ts_t);
625 if (!this->feScalingFun.empty())
626 this->locF *= this->feScalingFun(this->getFEMethod());
627 if (assembleVec) {
628 auto *vec_ptr = this->locF.data().data();
629 const auto nb_dofs = row_data.getIndices().size();
630 auto *ind_ptr = row_data.getIndices().data().data();
631 return VecSetValues(assembleVec, nb_dofs, ind_ptr, vec_ptr, ADD_VALUES);
632 }
633 if (thGradTag) {
634 const auto field_ents = data.getFieldEntities();
635 std::vector<EntityHandle> ents(field_ents.size());
636 std::transform(field_ents.begin(), field_ents.end(), ents.begin(),
637 [](const auto *fe) { return fe->getEnt(); });
638 if (field_ents.empty())
640 if (type_from_handle(ents[0]) != MBVERTEX)
642 auto &moab = getMoab();
643 VectorDouble topo_values(this->locF.size());
644 CHKERR moab.tag_set_data(thGradTag, ents.data(), ents.size(),
645 topo_values.data().data());
646 topo_values += this->locF;
647 CHKERR moab.tag_set_data(thGradTag, ents.data(), ents.size(),
648 this->locF.data().data());
649 }
650 }
651
652 boost::shared_ptr<MatrixDouble> adjointLambdaPtr;
653 boost::shared_ptr<MatrixDouble> tangent1DiffPtr;
654 boost::shared_ptr<MatrixDouble> tangent2DiffPtr;
657}
658
659template <typename OpBrokenBase>
661 : public OpBrokenTopoBase<OpBrokenBase> {
662
664
666 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
667 boost::shared_ptr<MatrixDouble> lagrange_ptr,
668 boost::shared_ptr<MatrixDouble> adjoint_lambda_ptr,
669 boost::shared_ptr<MatrixDouble> tangent1_diff_ptr,
670 boost::shared_ptr<MatrixDouble> tangent2_diff_ptr,
671 boost::shared_ptr<double> beta_ptr, SmartPetscObj<Vec> assemble_vec,
672 Tag th, boost::shared_ptr<Range> ents_ptr = nullptr)
673 : OP(broken_base_side_data, adjoint_lambda_ptr, tangent1_diff_ptr,
674 tangent2_diff_ptr, beta_ptr, assemble_vec, th, ents_ptr),
675 lagrangePtr(lagrange_ptr) {}
676
677private:
679
680 boost::shared_ptr<double> scalarBetaPtr;
681 boost::shared_ptr<MatrixDouble> lagrangePtr;
682};
683
684template <typename OpBase>
686OpTopoDerivativeBrokenSpaceConstrainDFluxImpl<3, GAUSS, OpBase>::iNtegrate(
687 EntitiesFieldData::EntData &row_data) {
689
690 FTENSOR_INDEX(3, I);
691 FTENSOR_INDEX(3, J);
692 FTENSOR_INDEX(3, i);
693 FTENSOR_INDEX(3, j);
694 FTENSOR_INDEX(3, k);
695 FTENSOR_INDEX(3, m);
698
699 auto t_w = this->getFTensor0IntegrationWeight();
700 auto t_lagrange = getFTensor1FromMat<3>(*lagrangePtr);
701 auto t_adjoint_lambda =
702 getFTensor2FromMat<3, 3>(*adjointLambdaPtr);
703 auto t_t1 = OpBase::getFTensor1Tangent1AtGaussPts();
704 auto t_t2 = OpBase::getFTensor1Tangent2AtGaussPts();
705
706 auto t_diff_base = row_data.getFTensor1DiffN<2>();
707 auto nb_base_functions = row_data.getN().size2();
708
709 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
710
711 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
712
713 FTensor::Tensor2<double, 3, 3> t_normal_diff;
714 t_normal_diff(j, m) =
715 FTensor::levi_civita(i, j, k) * t_t1(k) *
716 (t_kd(i, m) * t_diff_base(N1));
717 t_normal_diff(j, m) +=
718 FTensor::levi_civita(i, j, k) * (t_kd(k, m) * t_diff_base(N0)) *
719 t_t2(i);
720 auto t_vec = getFTensor1FromPtr<3>(&*OP::locF.data().begin());
721 size_t rr = 0;
722 for (; rr != row_data.getIndices().size() / 3; ++rr) {
723 t_vec(m) += (0.5 * t_w) * (t_adjoint_lambda(i, J) * t_lagrange(i)) *
724 t_normal_diff(J, m);
725 ++t_diff_base;
726 ++t_vec;
727 }
728 for (; rr < nb_base_functions; ++rr)
729 ++t_diff_base;
730 ++t_w;
731 ++t_lagrange;
732 ++t_adjoint_lambda;
733 ++t_t1;
734 ++t_t2;
735 }
736
737 if (scalarBetaPtr)
738 OP::locF *= *scalarBetaPtr;
739
741}
742
743template <typename OpBase>
745 : public OpBrokenTopoBase<OpBase> {
746
747 using OP = typename OpBrokenTopoBase<OpBase>;
748
750 const std::string row_field,
751 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
752 boost::shared_ptr<MatrixDouble> adjoint_lambda_ptr,
753 boost::shared_ptr<MatrixDouble> tangent1_ptr,
754 boost::shared_ptr<MatrixDouble> tangent2_ptr,
755 boost::shared_ptr<double> beta_ptr, SmartPetscObj<Vec> assemble_vec,
756 Tag th, boost::shared_ptr<Range> ents_ptr = nullptr)
757 : OpBrokenTopoBase<OpBase>(row_field, row_field, OpBase::OPROW, ents_ptr),
758 brokenSideDataPtr(broken_side_data_ptr), scalarBetaPtr(beta_ptr) {}
759
760private:
761 boost::shared_ptr<double> scalarBetaPtr;
762 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenSideDataPtr;
764};
765
766template <typename OpBase>
770
771 FTENSOR_INDEX(3, i);
772 FTENSOR_INDEX(3, j);
773 FTENSOR_INDEX(3, k);
774 FTENSOR_INDEX(3, m);
775 FTENSOR_INDEX(3, J);
776 FTENSOR_INDEX(3, K);
777 FTENSOR_INDEX(3, L);
778 FTENSOR_INDEX(3, M);
779
782
783 OP::locF.resize(row_data.getIndices().size(), false);
784 OP::locF.clear();
785
786 for (auto &bd : *brokenSideDataPtr) {
787 auto t_w = this->getFTensor0IntegrationWeight();
788 auto t_t1 = getFTensor1FromMat<3>(*tangent1Ptr);
789 auto t_t2 = getFTensor1FromMat<3>(*tangent2Ptr);
790 auto t_flux = getFTensor2FromMat<3, 3>(bd.getFlux());
791 auto t_adjoint_lambda = getFTensor1FromMat<3>(*adjointLambdaPtr);
792
793 auto nb_base_functions = row_data.getN().size2();
794 auto sense = bd.getSense();
795 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
796
797 auto t_diff_base = row_data.getFTensor1DiffN<2>();
798 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
799 auto t_vec = getFTensor1FromPtr<3>(&*OP::locF.data().begin());
800 size_t rr = 0;
801 for (; rr != row_data.getIndices().size() / 3; ++rr) {
802 FTensor::Tensor2<double, 3, 3> t_normal_diff;
803 t_normal_diff(j, m) =
804 FTensor::levi_civita(i, j, k) * t_t1(k) *
805 (t_kd(i, m) * t_diff_base(N1));
806 t_normal_diff(j, m) +=
807 FTensor::levi_civita(i, j, k) * (t_kd(k, m) * t_diff_base(N0)) *
808 t_t2(i);
809
810 t_vec(m) += (sense * 0.5 * t_w) * (t_adjoint_lambda(i) * t_flux(i, J)) *
811 t_normal_diff(J, m);
812
813 ++t_diff_base;
814 ++t_vec;
815 }
816 for (; rr < nb_base_functions; ++rr)
817 ++t_diff_base;
818 ++t_w;
819 ++t_t1;
820 ++t_t2;
821 ++t_flux;
822 ++t_adjoint_lambda;
823 }
824 }
825
826 if (scalarBetaPtr)
827 OP::locF *= *scalarBetaPtr;
828
830}
831
832} // namespace MoFEM
833
834#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
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
boost::shared_ptr< double > scalarBetaPtr
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data)
boost::shared_ptr< MatrixDouble > lagrangePtr
MoFEM::OpBrokenTopoBase GAUSS
Gaussian quadrature integration.
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.
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 > adjointLambdaPtr
OpBrokenTopoBase(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< MatrixDouble > adjoint_lambda_ptr, boost::shared_ptr< MatrixDouble > tangent1_diff_ptr, boost::shared_ptr< MatrixDouble > tangent2_diff_ptr, boost::shared_ptr< double > beta_ptr, SmartPetscObj< Vec > assemble_vec, Tag th, boost::shared_ptr< Range > ents_ptr=nullptr)
boost::shared_ptr< MatrixDouble > tangent2DiffPtr
boost::shared_ptr< MatrixDouble > tangent1DiffPtr
OpBrokenTopoBase(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_side_data_ptr, boost::shared_ptr< MatrixDouble > adjoint_lambda_ptr, boost::shared_ptr< MatrixDouble > tangent1_ptr, boost::shared_ptr< MatrixDouble > tangent2_ptr, SmartPetscObj< Vec > assemble_vec, Tag th, boost::shared_ptr< Range > ents_ptr=nullptr)
MoFEMErrorCode aSsemble(EntitiesFieldData::EntData &row_data)
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
OpTopoDerivativeBrokenSpaceConstrainDHybridImpl(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_side_data_ptr, boost::shared_ptr< MatrixDouble > adjoint_lambda_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