v0.15.0
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
114 if (row_data.getIndices().size() == 0)
116
117 brokenBaseSideData->resize(OP::getLoopSize());
118
119 const auto n_in_the_loop = OP::getNinTheLoop();
120 const auto face_sense = OP::getSkeletonSense();
121
122#ifndef NDEBUG
123 if (face_sense != -1 && face_sense != 1)
124 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "face sense not set");
125#endif // NDEBUG
126
127 auto set_data = [&](auto &side_data) {
128 side_data.getSide() = row_side;
129 side_data.getType() = row_type;
130 side_data.getSense() = face_sense;
131 side_data.getData().sEnse = row_data.sEnse;
132 side_data.getData().sPace = row_data.sPace;
133 side_data.getData().bAse = row_data.bAse;
134 side_data.getData().iNdices = row_data.iNdices;
135 side_data.getData().localIndices = row_data.localIndices;
136 side_data.getData().dOfs = row_data.dOfs;
137 side_data.getData().fieldEntities = row_data.fieldEntities;
138 side_data.getData().fieldData = row_data.fieldData;
139 };
140
141 auto set_base = [&](auto &side_data) {
142 auto base = side_data.getData().getBase();
143 for (auto dd = 0; dd != BaseDerivatives::LastDerivative; ++dd) {
144 side_data.getData().baseFunctionsAndBaseDerivatives[dd][base].reset();
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}
189
190template <typename OpBase> struct OpBrokenBaseImpl : public OpBase {
191
192 using OP = OpBase;
193
195 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
196 boost::shared_ptr<Range> ents_ptr = nullptr)
197 : OP(NOSPACE, OP::OPSPACE), brokenBaseSideData(broken_base_side_data) {
198 OP::entsPtr = ents_ptr;
199 }
200
202 const std::string row_field,
203 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
204 const bool assmb_transpose, const bool only_transpose,
205 boost::shared_ptr<Range> ents_ptr = nullptr)
206 : OP(row_field, row_field, OP::OPROW, ents_ptr),
207 brokenBaseSideData(broken_base_side_data) {
208 OP::entsPtr = ents_ptr;
209 OP::assembleTranspose = assmb_transpose;
210 OP::onlyTranspose = only_transpose;
211 OP::sYmm = false;
212 }
213
214 MoFEMErrorCode doWork(int row_side, EntityType row_type,
216
217protected:
218 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
219};
220
221template <typename OpBase>
223OpBrokenBaseImpl<OpBase>::doWork(int row_side, EntityType row_type,
224 EntitiesFieldData::EntData &row_data) {
226
227 if (OP::entsPtr) {
228 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
230 }
231
232#ifndef NDEBUG
233 if (!brokenBaseSideData) {
234 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not set");
235 }
236#endif // NDEBUG
237
238 auto do_work_rhs = [this](int row_side, EntityType row_type,
239 EntitiesFieldData::EntData &row_data, int sense) {
241 // get number of dofs on row
242 OP::nbRows = row_data.getIndices().size();
243 if (!OP::nbRows)
245 // get number of integration points
246 OP::nbIntegrationPts = OP::getGaussPts().size2();
247 // get row base functions
248 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
249 // resize and clear the right hand side vector
250 OP::locF.resize(OP::nbRows, false);
251 OP::locF.clear();
252 // integrate local vector
253 CHKERR this->iNtegrate(row_data);
254 // assemble local vector
255 OP::locF *= sense;
256 CHKERR this->aSsemble(row_data);
258 };
259
260 auto do_work_lhs = [this](int row_side, int col_side, EntityType row_type,
261 EntityType col_type,
263 EntitiesFieldData::EntData &col_data, int sense) {
265
266 auto check_if_assemble_transpose = [&] {
267 if (this->sYmm) {
268 if (OP::rowSide != OP::colSide || OP::rowType != OP::colType)
269 return true;
270 else
271 return false;
272 } else if (OP::assembleTranspose) {
273 return true;
274 }
275 return false;
276 };
277
278 OP::rowSide = row_side;
279 OP::rowType = row_type;
280 OP::colSide = col_side;
281 OP::colType = col_type;
282 OP::nbCols = col_data.getIndices().size();
283 OP::locMat.resize(OP::nbRows, OP::nbCols, false);
284 OP::locMat.clear();
285 CHKERR this->iNtegrate(row_data, col_data);
286 OP::locMat *= sense;
287 CHKERR this->aSsemble(row_data, col_data, check_if_assemble_transpose());
289 };
290
291 switch (OP::opType) {
292 case OP::OPROW:
293
294 OP::nbRows = row_data.getIndices().size();
295 if (!OP::nbRows)
297 OP::nbIntegrationPts = OP::getGaussPts().size2();
298 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
299
300 if (!OP::nbRows)
302
303 for (auto &bd : *brokenBaseSideData) {
304
305#ifndef NDEBUG
306 if (!bd.getData().getNSharedPtr(bd.getData().getBase())) {
307 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
308 "base functions not set");
309 }
310#endif
311
312 CHKERR do_work_lhs(
313
314 // side
315 row_side, bd.getSide(),
316
317 // type
318 row_type, bd.getType(),
319
320 // row_data
321 row_data, bd.getData(),
322
323 // sense
324 bd.getSense()
325
326 );
327 }
328
329 break;
330 case OP::OPSPACE:
331 for (auto &bd : *brokenBaseSideData) {
332 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData(),
333 bd.getSense());
334 }
335 break;
336 default:
338 (std::string("wrong op type ") +
340 .c_str());
341 }
342
344}
345
346template <int FIELD_DIM, IntegrationType I, typename OpBrokenBase>
348
349template <int FIELD_DIM, typename OpBrokenBase>
351 : public OpBrokenBase {
352
354
356 const std::string row_field,
357 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
358 boost::shared_ptr<double> beta_ptr, const bool assmb_transpose,
359 const bool only_transpose, boost::shared_ptr<Range> ents_ptr = nullptr)
360 : OP(row_field, broken_base_side_data, assmb_transpose, only_transpose,
361 ents_ptr),
362 scalarBetaPtr(beta_ptr) {}
363
365 const std::string row_field,
366 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
367 double beta, const bool assmb_transpose, const bool only_transpose,
368 boost::shared_ptr<Range> ents_ptr = nullptr)
369 : OpBrokenSpaceConstrainImpl(row_field, broken_base_side_data,
370 boost::make_shared<double>(beta),
371 assmb_transpose, only_transpose, ents_ptr) {}
372
373protected:
374 boost::shared_ptr<double> scalarBetaPtr;
377};
378
379template <int FIELD_DIM, typename OpBase>
382 EntitiesFieldData::EntData &col_data) {
384
385 auto nb_row_dofs = row_data.getIndices().size();
386 auto nb_col_dofs = col_data.getIndices().size();
387 if (!nb_row_dofs || !nb_col_dofs)
389
391 FTENSOR_INDEX(3, J);
392
393 auto t_w = this->getFTensor0IntegrationWeight();
394 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
395 size_t nb_base_functions = col_data.getN().size2() / 3;
396
397#ifndef NDEBUG
398 if (nb_row_dofs % FIELD_DIM != 0 || nb_col_dofs % FIELD_DIM != 0) {
399 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
400 "number of dofs not divisible by field dimension");
401 }
402 if (nb_row_dofs > row_data.getN().size2() * FIELD_DIM ||
403 nb_col_dofs > col_data.getN().size2() * FIELD_DIM) {
404 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
405 "number of dofs exceeds number of base functions");
406 }
407#endif // NDEBUG
408
409 auto t_col_base = col_data.getFTensor1N<3>();
410
411 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
412
413 int cc = 0;
414 for (; cc != nb_col_dofs / FIELD_DIM; cc++) {
415 auto t_row_base = row_data.getFTensor0N(gg, 0);
416 for (auto rr = 0; rr != nb_row_dofs / FIELD_DIM; ++rr) {
417 OP::locMat(FIELD_DIM * rr, FIELD_DIM * cc) +=
418 (0.5 * t_w * t_row_base) * (t_normal(J) * t_col_base(J));
419 ++t_row_base;
420 }
421
422 ++t_col_base;
423 }
424
425 for (; cc < nb_base_functions; ++cc)
426 ++t_col_base;
427
428 ++t_w;
429 ++t_normal;
430 }
431
432 for (auto rr = 0; rr != nb_row_dofs / FIELD_DIM; ++rr) {
433 for (auto cc = 0; cc != nb_col_dofs / FIELD_DIM; ++cc) {
434 for (auto dd = 1; dd < FIELD_DIM; ++dd) {
435 OP::locMat(FIELD_DIM * rr + dd, FIELD_DIM * cc + dd) =
436 OP::locMat(FIELD_DIM * rr, FIELD_DIM * cc);
437 }
438 }
439 }
440
441 if (scalarBetaPtr)
442 OP::locMat *= *scalarBetaPtr;
443
445}
446
447template <int FIELD_DIM, IntegrationType I, typename OpBase>
449
450template <int FIELD_DIM, IntegrationType I, typename OpBase>
452
453template <int FIELD_DIM, typename OpBrokenBase>
455 : public OpBrokenBase {
456
458
460 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
461 boost::shared_ptr<MatrixDouble> lagrange_ptr,
462 boost::shared_ptr<double> beta_ptr,
463 boost::shared_ptr<Range> ents_ptr = nullptr)
464 : OP(broken_base_side_data, ents_ptr), scalarBetaPtr(beta_ptr),
465 lagrangePtr(lagrange_ptr) {}
466
468 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
469 boost::shared_ptr<MatrixDouble> lagrange_ptr, double beta,
470 boost::shared_ptr<Range> ents_ptr = nullptr)
471 : OpBrokenSpaceConstrainDFluxImpl(broken_base_side_data, lagrange_ptr,
472 boost::make_shared<double>(beta),
473 ents_ptr) {}
474
475private:
477
478 boost::shared_ptr<double> scalarBetaPtr;
479 boost::shared_ptr<MatrixDouble> lagrangePtr;
480};
481
482template <int FIELD_DIM, typename OpBase>
485 EntitiesFieldData::EntData &row_data) {
487
489 FTENSOR_INDEX(3, J);
490
491 auto t_w = this->getFTensor0IntegrationWeight();
492 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
493 auto t_lagrange = getFTensor1FromMat<FIELD_DIM>(*lagrangePtr);
494
495 auto t_row_base = row_data.getFTensor1N<3>();
496 auto nb_base_functions = row_data.getN().size2() / 3;
497
498 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
499 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
500 size_t rr = 0;
501 for (; rr != row_data.getIndices().size() / FIELD_DIM; ++rr) {
502 t_vec(i) += (0.5 * t_w * (t_row_base(J) * t_normal(J))) * t_lagrange(i);
503 ++t_row_base;
504 ++t_vec;
505 }
506 for (; rr < nb_base_functions; ++rr)
507 ++t_row_base;
508 ++t_w;
509 ++t_normal;
510 ++t_lagrange;
511 }
512
513 if(scalarBetaPtr)
514 OP::locF *= *scalarBetaPtr;
515
517}
518
519template <int FIELD_DIM, typename OpBase>
521 : public OpBase {
522
523 using OP = OpBase;
524
526 const std::string row_field,
527 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
528 boost::shared_ptr<double> beta_ptr,
529 boost::shared_ptr<Range> ents_ptr = nullptr)
530 : OpBase(row_field, row_field, OpBase::OPROW, ents_ptr),
531 brokenSideDataPtr(broken_side_data_ptr), scalarBetaPtr(beta_ptr) {}
532
534 const std::string row_field,
535 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
536 double beta, boost::shared_ptr<Range> ents_ptr = nullptr)
537 : OpBrokenSpaceConstrainDHybridImpl(row_field, broken_side_data_ptr,
538 boost::make_shared<double>(beta),
539 ents_ptr) {}
540
541private:
542 boost::shared_ptr<double> scalarBetaPtr;
543 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenSideDataPtr;
544 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data);
545};
546
547template <int FIELD_DIM, typename OpBase>
550 EntitiesFieldData::EntData &row_data) {
552
554 FTENSOR_INDEX(3, J);
555
556 OP::locF.resize(row_data.getIndices().size(), false);
557 OP::locF.clear();
558
559 for (auto &bd : *brokenSideDataPtr) {
560 auto t_w = this->getFTensor0IntegrationWeight();
561 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
562 auto t_row_base = row_data.getFTensor0N();
563 auto t_flux = getFTensor2FromMat<FIELD_DIM, 3>(bd.getFlux());
564 auto nb_base_functions = row_data.getN().size2() / 3;
565 auto sense = bd.getSense();
566 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
567 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
568 size_t rr = 0;
569 for (; rr != row_data.getIndices().size() / FIELD_DIM; ++rr) {
570 t_vec(i) +=
571 (sense * 0.5 * t_w) * t_row_base * t_normal(J) * t_flux(i, J);
572 ++t_row_base;
573 ++t_vec;
574 }
575 for (; rr < nb_base_functions; ++rr)
576 ++t_row_base;
577 ++t_w;
578 ++t_normal;
579 ++t_flux;
580 }
581 }
582
583 if(scalarBetaPtr)
584 OP::locF *= *scalarBetaPtr;
585
587}
588
589} // namespace MoFEM
590
591#endif // __FORMSBROKENSPACECONSTRAINTIMPL_HPP__
#define FTENSOR_INDEX(DIM, I)
constexpr int FIELD_DIM
#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.
@ GAUSS
Gaussian quadrature integration.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
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
constexpr auto field_name
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Definition radiation.cpp:29
Data on single entity (This is passed as argument to DataOperator::doWork)
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)
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