10#ifndef __FORMSBROKENSPACECONSTRAINTIMPL_HPP__
11#define __FORMSBROKENSPACECONSTRAINTIMPL_HPP__
45 prev_side_fe_ptr->numeredEntFiniteElementPtr->getFEUId();
58 prev_side_fe_ptr->numeredEntFiniteElementPtr;
63 "sideFEName is different");
92 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data);
101template <
typename OpBase>
104 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data)
106 brokenBaseSideData(broken_base_side_data) {}
108template <
typename OpBase>
113 brokenBaseSideData->resize(OP::getLoopSize());
115 const auto n_in_the_loop = OP::getNinTheLoop();
116 const auto face_sense = OP::getSkeletonSense();
119 if (face_sense != -1 && face_sense != 1)
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;
134 side_data.getData().fieldData = row_data.
fieldData;
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();
144 side_data.getData().baseFunctionsAndBaseDerivatives[dd][base] =
145 boost::make_shared<MatrixDouble>(*base_ptr);
150 set_data((*brokenBaseSideData)[n_in_the_loop]);
151 set_base((*brokenBaseSideData)[n_in_the_loop]);
161 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
162 boost::shared_ptr<MatrixDouble> flux_ptr);
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),
179template <
typename OpBase>
183 auto swap_flux = [&](
auto &side_data) { side_data.getFlux().swap(*fluxPtr); };
184 swap_flux((*brokenBaseSideData)[OP::getNinTheLoop()]);
193 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
194 boost::shared_ptr<Range> ents_ptr =
nullptr)
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),
219template <
typename OpBase>
226 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
231 if (!brokenBaseSideData) {
236 auto do_work_rhs = [
this](
int row_side, EntityType row_type,
244 OP::nbIntegrationPts = OP::getGaussPts().size2();
246 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
248 OP::locF.resize(OP::nbRows,
false);
254 CHKERR this->aSsemble(row_data);
258 auto do_work_lhs = [
this](
int row_side,
int col_side, EntityType row_type,
264 auto check_if_assemble_transpose = [&] {
266 if (OP::rowSide != OP::colSide || OP::rowType != OP::colType)
270 }
else if (OP::assembleTranspose) {
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);
285 CHKERR this->aSsemble(row_data, col_data, check_if_assemble_transpose());
289 switch (OP::opType) {
295 OP::nbIntegrationPts = OP::getGaussPts().size2();
296 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
301 for (
auto &bd : *brokenBaseSideData) {
304 if (!bd.getData().getNSharedPtr(bd.getData().getBase())) {
306 "base functions not set");
313 row_side, bd.getSide(),
316 row_type, bd.getType(),
319 row_data, bd.getData(),
329 for (
auto &bd : *brokenBaseSideData) {
330 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData(),
336 (std::string(
"wrong op type ") +
344template <
int FIELD_DIM, IntegrationType I,
typename OpBrokenBase>
347template <
int FIELD_DIM,
typename OpBrokenBase>
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,
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)
369 assmb_transpose, only_transpose, ents_ptr) {}
377template <
int FIELD_DIM,
typename OpBase>
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)
391 auto t_w = this->getFTensor0IntegrationWeight();
392 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
393 size_t nb_base_functions = col_data.
getN().size2() / 3;
398 "number of dofs not divisible by field dimension");
403 "number of dofs exceeds number of base functions");
412 for (; cc != nb_col_dofs /
FIELD_DIM; cc++) {
414 for (
auto rr = 0; rr != nb_row_dofs /
FIELD_DIM; ++rr) {
416 (0.5 * t_w * t_row_base) * (t_normal(
J) * t_col_base(
J));
423 for (; cc < nb_base_functions; ++cc)
430 for (
auto rr = 0; rr != nb_row_dofs /
FIELD_DIM; ++rr) {
431 for (
auto cc = 0; cc != nb_col_dofs /
FIELD_DIM; ++cc) {
445template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
448template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
451template <
int FIELD_DIM,
typename OpBrokenBase>
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)
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)
480template <
int FIELD_DIM,
typename OpBase>
489 auto t_w = this->getFTensor0IntegrationWeight();
490 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
491 auto t_lagrange = getFTensor1FromMat<FIELD_DIM>(*
lagrangePtr);
494 auto nb_base_functions = row_data.
getN().size2() / 3;
497 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
500 t_vec(
i) += (0.5 * t_w * (t_row_base(
J) * t_normal(
J))) * t_lagrange(
i);
504 for (; rr < nb_base_functions; ++rr)
517template <
int FIELD_DIM,
typename OpBase>
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) {}
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)
545template <
int FIELD_DIM,
typename OpBase>
554 OP::locF.resize(row_data.
getIndices().size(),
false);
557 for (
auto &bd : *brokenSideDataPtr) {
558 auto t_w = this->getFTensor0IntegrationWeight();
559 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
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();
565 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
569 (sense * 0.5 * t_w) * t_row_base * t_normal(
J) * t_flux(
i,
J);
573 for (; rr < nb_base_functions; ++rr)
587template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
590template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
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,
602 Tag th, boost::shared_ptr<Range> ents_ptr =
nullptr)
603 :
OP(broken_base_side_data, ents_ptr),
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,
615 boost::shared_ptr<Range> ents_ptr =
nullptr)
616 :
OpBase(row_field, row_field,
OpBase::OPROW, ents_ptr),
617 brokenSideDataPtr(broken_side_data_ptr),
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());
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();
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())
642 auto &moab = getMoab();
645 topo_values.data().data());
646 topo_values += this->locF;
648 this->locF.data().data());
659template <
typename OpBrokenBase>
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,
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),
684template <
typename OpBase>
686OpTopoDerivativeBrokenSpaceConstrainDFluxImpl<3, GAUSS, OpBase>::iNtegrate(
687 EntitiesFieldData::EntData &row_data) {
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();
706 auto t_diff_base = row_data.getFTensor1DiffN<2>();
707 auto nb_base_functions = row_data.getN().size2();
714 t_normal_diff(
j,
m) =
716 (
t_kd(
i,
m) * t_diff_base(N1));
717 t_normal_diff(
j,
m) +=
720 auto t_vec = getFTensor1FromPtr<3>(&*OP::locF.data().begin());
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)) *
728 for (; rr < nb_base_functions; ++rr)
743template <
typename OpBase>
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,
756 Tag th, boost::shared_ptr<Range> ents_ptr =
nullptr)
758 brokenSideDataPtr(broken_side_data_ptr),
scalarBetaPtr(beta_ptr) {}
766template <
typename OpBase>
783 OP::locF.resize(row_data.
getIndices().size(),
false);
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);
793 auto nb_base_functions = row_data.
getN().size2();
794 auto sense = bd.getSense();
799 auto t_vec = getFTensor1FromPtr<3>(&*OP::locF.data().begin());
801 for (; rr != row_data.
getIndices().size() / 3; ++rr) {
803 t_normal_diff(
j,
m) =
805 (
t_kd(
i,
m) * t_diff_base(N1));
806 t_normal_diff(
j,
m) +=
810 t_vec(
m) += (sense * 0.5 * t_w) * (t_adjoint_lambda(
i) * t_flux(
i,
J)) *
816 for (; rr < nb_base_functions; ++rr)
#define FTENSOR_INDEX(DIM, I)
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#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_DATA_INCONSISTENCY
#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
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)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
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
FTensor::Index< 'm', 3 > m
EntitiesFieldData::EntData entData
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.
FieldSpace sPace
Entity space.
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
VectorDofs dOfs
DoFs on entity.
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
FieldApproximationBase bAse
Field approximation base.
static const char *const OpTypeNames[]
ForcesAndSourcesCore * getPtrFE() const
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
ForcesAndSourcesCore * getSidePtrFE() const
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.
boost::shared_ptr< double > scalarBetaPtr
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)
boost::shared_ptr< MatrixDouble > lagrangePtr
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)
boost::shared_ptr< double > scalarBetaPtr
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenSideDataPtr
boost::shared_ptr< double > scalarBetaPtr
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)
SmartPetscObj< Vec > assembleVec
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
boost::shared_ptr< double > scalarBetaPtr
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenSideDataPtr
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