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
15template <typename E> struct OpBrokenLoopSide : public OpLoopSide<E> {
16
18 using OpLoopSide<E>::OpLoopSide;
19
20 MoFEMErrorCode doWork(int side, EntityType type,
23
24 auto prev_side_fe_ptr = OP::getSidePtrFE();
25 if (OP::sideFEName == prev_side_fe_ptr->getFEName()) {
26 auto prev_fe_uid =
27 prev_side_fe_ptr->numeredEntFiniteElementPtr->getFEUId();
29 CHKERR OP::sideFEPtr->setSideFEPtr(OP::getPtrFE());
30 CHKERR OP::sideFEPtr->copyBasicMethod(*OP::getFEMethod());
31 CHKERR OP::sideFEPtr->copyPetscData(*OP::getFEMethod());
35 OP::sideFEPtr->cacheWeakPtr = prev_side_fe_ptr->cacheWeakPtr;
36 OP::sideFEPtr->loopSize = 1;
37 CHKERR OP::sideFEPtr->preProcess();
38 OP::sideFEPtr->nInTheLoop = 0;
39 OP::sideFEPtr->numeredEntFiniteElementPtr =
40 prev_side_fe_ptr->numeredEntFiniteElementPtr;
42 CHKERR OP::sideFEPtr->postProcess();
43 } else {
44 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
45 "sideFEName is different");
46 }
47
49 };
50};
51
54 inline auto &getSense() { return eleSense; }
55 inline auto &getSide() { return eleSide; }
56 inline auto &getType() { return eleType; }
57 inline auto &getData() { return entData; }
58 inline auto &getFlux() { return fluxMat; }
59
60private:
61 int eleSense = 0;
62 int eleSide = 1;
63 EntityType eleType = MBENTITYSET;
66};
67
68template <typename OpBase> struct OpGetBrokenBaseSideData : public OpBase {
69
70 using OP = OpBase;
71
73 const std::string field_name,
74 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data);
75
76 MoFEMErrorCode doWork(int row_side, EntityType row_type,
78
79private:
80 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
81};
82
83template <typename OpBase>
85 const std::string field_name,
86 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data)
87 : OP(field_name, field_name, OP::OPROW),
88 brokenBaseSideData(broken_base_side_data) {}
89
90template <typename OpBase>
92OpGetBrokenBaseSideData<OpBase>::doWork(int row_side, EntityType row_type,
95
96 if (row_data.getIndices().size() == 0)
98
99 brokenBaseSideData->resize(OP::getLoopSize());
100
101 const auto n_in_the_loop = OP::getNinTheLoop();
102 const auto face_sense = OP::getSkeletonSense();
103
104#ifndef NDEBUG
105 if (face_sense != -1 && face_sense != 1)
106 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "face sense not set");
107#endif // NDEBUG
108
109 auto set_data = [&](auto &side_data) {
110 side_data.getSide() = row_side;
111 side_data.getType() = row_type;
112 side_data.getSense() = face_sense;
113 side_data.getData().sEnse = row_data.sEnse;
114 side_data.getData().sPace = row_data.sPace;
115 side_data.getData().bAse = row_data.bAse;
116 side_data.getData().iNdices = row_data.iNdices;
117 side_data.getData().localIndices = row_data.localIndices;
118 side_data.getData().dOfs = row_data.dOfs;
119 side_data.getData().fieldEntities = row_data.fieldEntities;
120 side_data.getData().fieldData = row_data.fieldData;
121 };
122
123 auto set_base = [&](auto &side_data) {
124 auto base = side_data.getData().getBase();
125 for (auto dd = 0; dd != BaseDerivatives::LastDerivative; ++dd) {
126 if (auto base_ptr = row_data.baseFunctionsAndBaseDerivatives[dd][base]) {
127 side_data.getData().baseFunctionsAndBaseDerivatives[dd][base] =
128 boost::make_shared<MatrixDouble>(*base_ptr);
129 }
130 }
131 };
132
133 set_data((*brokenBaseSideData)[n_in_the_loop]);
134 set_base((*brokenBaseSideData)[n_in_the_loop]);
135
137}
138
139template <typename OpBase> struct OpSetFlux : public OpBase {
140
141 using OP = OpBase;
142
143 OpSetFlux(
144 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
145 boost::shared_ptr<MatrixDouble> flux_ptr);
146
147 MoFEMErrorCode doWork(int row_side, EntityType row_type,
149
150private:
151 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
152 boost::shared_ptr<MatrixDouble> fluxPtr;
153};
154
155template <typename OpBase>
157 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
158 boost::shared_ptr<MatrixDouble> flux_ptr)
159 : OP(NOSPACE, OP::OPSPACE), brokenBaseSideData(broken_base_side_data),
160 fluxPtr(flux_ptr) {}
161
162template <typename OpBase>
166 auto swap_flux = [&](auto &side_data) { side_data.getFlux().swap(*fluxPtr); };
167 swap_flux((*brokenBaseSideData)[OP::getNinTheLoop()]);
169}
170
171template <typename OpBase> struct OpBrokenBaseImpl : public OpBase {
172
173 using OP = OpBase;
174
176 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
177 boost::shared_ptr<Range> ents_ptr = nullptr)
178 : OP(NOSPACE, OP::OPSPACE), brokenBaseSideData(broken_base_side_data) {
179 OP::entsPtr = ents_ptr;
180 }
181
183 const std::string row_field,
184 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
185 const bool assmb_transpose, const bool only_transpose,
186 boost::shared_ptr<Range> ents_ptr = nullptr)
187 : OP(row_field, row_field, OP::OPROW, ents_ptr),
188 brokenBaseSideData(broken_base_side_data) {
189 OP::entsPtr = ents_ptr;
190 OP::assembleTranspose = assmb_transpose;
191 OP::onlyTranspose = only_transpose;
192 OP::sYmm = false;
193 }
194
195 MoFEMErrorCode doWork(int row_side, EntityType row_type,
197
198protected:
199 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideData;
200};
201
202template <typename OpBase>
204OpBrokenBaseImpl<OpBase>::doWork(int row_side, EntityType row_type,
205 EntitiesFieldData::EntData &row_data) {
207
208 if (OP::entsPtr) {
209 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
211 }
212
213#ifndef NDEBUG
214 if (!brokenBaseSideData) {
215 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not set");
216 }
217#endif // NDEBUG
218
219 auto do_work_rhs = [this](int row_side, EntityType row_type,
220 EntitiesFieldData::EntData &row_data, int sense) {
222 // get number of dofs on row
223 OP::nbRows = row_data.getIndices().size();
224 if (!OP::nbRows)
226 // get number of integration points
227 OP::nbIntegrationPts = OP::getGaussPts().size2();
228 // get row base functions
229 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
230 // resize and clear the right hand side vector
231 OP::locF.resize(OP::nbRows, false);
232 OP::locF.clear();
233 // integrate local vector
234 CHKERR this->iNtegrate(row_data);
235 // assemble local vector
236 OP::locF *= sense;
237 CHKERR this->aSsemble(row_data);
239 };
240
241 auto do_work_lhs = [this](int row_side, int col_side, EntityType row_type,
242 EntityType col_type,
244 EntitiesFieldData::EntData &col_data, int sense) {
246
247 auto check_if_assemble_transpose = [&] {
248 if (this->sYmm) {
249 if (OP::rowSide != OP::colSide || OP::rowType != OP::colType)
250 return true;
251 else
252 return false;
253 } else if (OP::assembleTranspose) {
254 return true;
255 }
256 return false;
257 };
258
259 OP::rowSide = row_side;
260 OP::rowType = row_type;
261 OP::colSide = col_side;
262 OP::colType = col_type;
263 OP::nbCols = col_data.getIndices().size();
264 OP::locMat.resize(OP::nbRows, OP::nbCols, false);
265 OP::locMat.clear();
266 CHKERR this->iNtegrate(row_data, col_data);
267 OP::locMat *= sense;
268 CHKERR this->aSsemble(row_data, col_data, check_if_assemble_transpose());
270 };
271
272 switch (OP::opType) {
273 case OP::OPROW:
274
275 OP::nbRows = row_data.getIndices().size();
276 if (!OP::nbRows)
278 OP::nbIntegrationPts = OP::getGaussPts().size2();
279 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
280
281 if (!OP::nbRows)
283
284 for (auto &bd : *brokenBaseSideData) {
285
286#ifndef NDEBUG
287 if (bd.getData().getSpace() != HDIV && bd.getData().getSpace() != HCURL) {
288 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "%s",
289 (std::string("Expect Hdiv or Hcurl space but received ") +
290 FieldSpaceNames[bd.getData().getSpace()])
291 .c_str());
292 }
293 if (!bd.getData().getNSharedPtr(bd.getData().getBase())) {
294 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
295 "base functions not set");
296 }
297#endif
298
299 CHKERR do_work_lhs(
300
301 // side
302 row_side, bd.getSide(),
303
304 // type
305 row_type, bd.getType(),
306
307 // row_data
308 row_data, bd.getData(),
309
310 // sense
311 bd.getSense()
312
313 );
314 }
315
316 break;
317 case OP::OPSPACE:
318 for (auto &bd : *brokenBaseSideData) {
319 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData(),
320 bd.getSense());
321 }
322 break;
323 default:
325 (std::string("wrong op type ") +
327 .c_str());
328 }
329
331}
332
333template <int FIELD_DIM, IntegrationType I, typename OpBrokenBase>
335
336template <int FIELD_DIM, typename OpBrokenBase>
338 : public OpBrokenBase {
339
341
343 const std::string row_field,
344 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
345 boost::shared_ptr<double> beta_ptr, const bool assmb_transpose,
346 const bool only_transpose, boost::shared_ptr<Range> ents_ptr = nullptr)
347 : OP(row_field, broken_base_side_data, assmb_transpose, only_transpose,
348 ents_ptr),
349 scalarBetaPtr(beta_ptr) {}
350
352 const std::string row_field,
353 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
354 double beta, const bool assmb_transpose, const bool only_transpose,
355 boost::shared_ptr<Range> ents_ptr = nullptr)
356 : OpBrokenSpaceConstrainImpl(row_field, broken_base_side_data,
357 boost::make_shared<double>(beta),
358 assmb_transpose, only_transpose, ents_ptr) {}
359
360protected:
361 boost::shared_ptr<double> scalarBetaPtr;
364};
365
366template <int FIELD_DIM, typename OpBase>
369 EntitiesFieldData::EntData &col_data) {
371
372 auto nb_row_dofs = row_data.getIndices().size();
373 auto nb_col_dofs = col_data.getIndices().size();
374 if (!nb_row_dofs || !nb_col_dofs)
376
378 FTENSOR_INDEX(3, J);
379
380 auto t_w = this->getFTensor0IntegrationWeight();
381 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
382 size_t nb_base_functions = col_data.getN().size2() / 3;
383
384#ifndef NDEBUG
385 if (nb_row_dofs % FIELD_DIM != 0 || nb_col_dofs % FIELD_DIM != 0) {
386 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
387 "number of dofs not divisible by field dimension");
388 }
389 if (nb_row_dofs > row_data.getN().size2() * FIELD_DIM ||
390 nb_col_dofs > col_data.getN().size2() * FIELD_DIM) {
391 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
392 "number of dofs exceeds number of base functions");
393 }
394#endif // NDEBUG
395
396 auto t_col_base = col_data.getFTensor1N<3>();
397
398 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
399
400 int cc = 0;
401 for (; cc != nb_col_dofs / FIELD_DIM; cc++) {
402 auto t_row_base = row_data.getFTensor0N(gg, 0);
403 for (auto rr = 0; rr != nb_row_dofs / FIELD_DIM; ++rr) {
404 OP::locMat(FIELD_DIM * rr, FIELD_DIM * cc) +=
405 (t_w * t_row_base) * (t_normal(J) * t_col_base(J));
406 ++t_row_base;
407 }
408
409 ++t_col_base;
410 }
411
412 for (; cc < nb_base_functions; ++cc)
413 ++t_col_base;
414
415 ++t_w;
416 ++t_normal;
417 }
418
419 for (auto rr = 0; rr != nb_row_dofs / FIELD_DIM; ++rr) {
420 for (auto cc = 0; cc != nb_col_dofs / FIELD_DIM; ++cc) {
421 for (auto dd = 1; dd < FIELD_DIM; ++dd) {
422 OP::locMat(FIELD_DIM * rr + dd, FIELD_DIM * cc + dd) =
423 OP::locMat(FIELD_DIM * rr, FIELD_DIM * cc);
424 }
425 }
426 }
427
428 if (scalarBetaPtr)
429 OP::locMat *= *scalarBetaPtr;
430
432}
433
434template <int FIELD_DIM, IntegrationType I, typename OpBase>
436
437template <int FIELD_DIM, IntegrationType I, typename OpBase>
439
440template <int FIELD_DIM, typename OpBrokenBase>
442 : public OpBrokenBase {
443
445
447 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
448 boost::shared_ptr<MatrixDouble> lagrange_ptr,
449 boost::shared_ptr<double> beta_ptr,
450 boost::shared_ptr<Range> ents_ptr = nullptr)
451 : OP(broken_base_side_data, ents_ptr), scalarBetaPtr(beta_ptr),
452 lagrangePtr(lagrange_ptr) {}
453
455 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
456 boost::shared_ptr<MatrixDouble> lagrange_ptr, double beta,
457 boost::shared_ptr<Range> ents_ptr = nullptr)
458 : OpBrokenSpaceConstrainDFluxImpl(broken_base_side_data, lagrange_ptr,
459 boost::make_shared<double>(beta),
460 ents_ptr) {}
461
462private:
464
465 boost::shared_ptr<double> scalarBetaPtr;
466 boost::shared_ptr<MatrixDouble> lagrangePtr;
467};
468
469template <int FIELD_DIM, typename OpBase>
472 EntitiesFieldData::EntData &row_data) {
474
476 FTENSOR_INDEX(3, J);
477
478 auto t_w = this->getFTensor0IntegrationWeight();
479 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
480 auto t_lagrange = getFTensor1FromMat<FIELD_DIM>(*lagrangePtr);
481
482 auto t_row_base = row_data.getFTensor1N<3>();
483 auto nb_base_functions = row_data.getN().size2() / 3;
484
485 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
486 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
487 size_t rr = 0;
488 for (; rr != row_data.getIndices().size() / FIELD_DIM; ++rr) {
489 t_vec(i) += (t_w * (t_row_base(J) * t_normal(J))) * t_lagrange(i);
490 ++t_row_base;
491 ++t_vec;
492 }
493 for (; rr < nb_base_functions; ++rr)
494 ++t_row_base;
495 ++t_w;
496 ++t_normal;
497 ++t_lagrange;
498 }
499
500 if(scalarBetaPtr)
501 OP::locF *= *scalarBetaPtr;
502
504}
505
506template <int FIELD_DIM, typename OpBase>
508 : public OpBase {
509
510 using OP = OpBase;
511
513 const std::string row_field,
514 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
515 boost::shared_ptr<double> beta_ptr,
516 boost::shared_ptr<Range> ents_ptr = nullptr)
517 : OpBase(row_field, row_field, OpBase::OPROW, ents_ptr),
518 brokenSideDataPtr(broken_side_data_ptr), scalarBetaPtr(beta_ptr) {}
519
521 const std::string row_field,
522 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
523 double beta, boost::shared_ptr<Range> ents_ptr = nullptr)
524 : OpBrokenSpaceConstrainDHybridImpl(row_field, broken_side_data_ptr,
525 boost::make_shared<double>(beta),
526 ents_ptr) {}
527
528private:
529 boost::shared_ptr<double> scalarBetaPtr;
530 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenSideDataPtr;
531 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data);
532};
533
534template <int FIELD_DIM, typename OpBase>
537 EntitiesFieldData::EntData &row_data) {
539
541 FTENSOR_INDEX(3, J);
542
543 OP::locF.resize(row_data.getIndices().size(), false);
544 OP::locF.clear();
545
546 for (auto &bd : *brokenSideDataPtr) {
547 auto t_w = this->getFTensor0IntegrationWeight();
548 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
549 auto t_row_base = row_data.getFTensor0N();
550 auto t_flux = getFTensor2FromMat<FIELD_DIM, 3>(bd.getFlux());
551 auto nb_base_functions = row_data.getN().size2() / 3;
552 auto sense = bd.getSense();
553 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
554 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
555 size_t rr = 0;
556 for (; rr != row_data.getIndices().size() / FIELD_DIM; ++rr) {
557 t_vec(i) += (sense * t_w) * t_row_base * t_normal(J) * t_flux(i, J);
558 ++t_row_base;
559 ++t_vec;
560 }
561 for (; rr < nb_base_functions; ++rr)
562 ++t_row_base;
563 ++t_w;
564 ++t_normal;
565 ++t_flux;
566 }
567 }
568
569 if(scalarBetaPtr)
570 OP::locF *= *scalarBetaPtr;
571
573}
574
575} // namespace MoFEM
576
577#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
@ HCURL
field with continuous tangents
Definition definitions.h:86
@ HDIV
field with continuous normal traction
Definition definitions.h:87
static const char *const FieldSpaceNames[]
Definition definitions.h:92
#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.
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, need to calculate base functions with 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 dofs 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
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)
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
OpBrokenBaseImpl(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< Range > ents_ptr=nullptr)
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