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");
94 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data);
103template <
typename OpBase>
106 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data)
108 brokenBaseSideData(broken_base_side_data) {}
110template <
typename OpBase>
115 brokenBaseSideData->resize(OP::getLoopSize());
117 const auto n_in_the_loop = OP::getNinTheLoop();
118 const auto face_sense = OP::getSkeletonSense();
121 if (face_sense != -1 && face_sense != 1)
125 auto set_data = [&](
auto &side_data) {
126 side_data.getSide() = row_side;
127 side_data.getType() = row_type;
128 side_data.getSense() = face_sense;
129 side_data.getData().sEnse = row_data.
sEnse;
130 side_data.getData().sPace = row_data.
sPace;
131 side_data.getData().bAse = row_data.
bAse;
132 side_data.getData().iNdices = row_data.
iNdices;
133 side_data.getData().localIndices = row_data.
localIndices;
134 side_data.getData().dOfs = row_data.
dOfs;
136 side_data.getData().fieldData = row_data.
fieldData;
139 auto set_base = [&](
auto &side_data) {
140 auto base = side_data.getData().getBase();
141 for (
auto dd = 0; dd != BaseDerivatives::LastDerivative; ++dd) {
142 for (
auto bb = 0; bb !=
LASTBASE; ++bb) {
143 side_data.getData().baseFunctionsAndBaseDerivatives[dd][bb].reset();
146 side_data.getData().baseFunctionsAndBaseDerivatives[dd][base] =
147 boost::make_shared<MatrixDouble>(*base_ptr);
152 set_data((*brokenBaseSideData)[n_in_the_loop]);
153 set_base((*brokenBaseSideData)[n_in_the_loop]);
163 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
164 boost::shared_ptr<MatrixDouble> flux_ptr);
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),
181template <
typename OpBase>
185 auto swap_flux = [&](
auto &side_data) { side_data.getFlux().swap(*fluxPtr); };
186 swap_flux((*brokenBaseSideData)[OP::getNinTheLoop()]);
194 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
195 boost::shared_ptr<MatrixDouble> flux_var_ptr);
205template <
typename OpBase>
207 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
208 boost::shared_ptr<MatrixDouble> flux_var_ptr)
209 :
OP(
NOSPACE,
OP::OPSPACE), brokenBaseSideData(broken_base_side_data),
210 fluxVarPtr(flux_var_ptr) {}
212template <
typename OpBase>
216 auto swap_flux = [&](
auto &side_data) {
217 side_data.getVarFlux().swap(*fluxVarPtr);
219 swap_flux((*brokenBaseSideData)[OP::getNinTheLoop()]);
229 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
230 boost::shared_ptr<Range> ents_ptr =
nullptr)
239 const std::string row_field,
240 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
241 const bool assmb_transpose,
const bool only_transpose,
242 boost::shared_ptr<Range> ents_ptr =
nullptr)
243 :
OP(row_field, row_field,
OP::OPROW, ents_ptr),
258template <
typename OpBase>
265 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
270 if (!brokenBaseSideData) {
275 auto do_work_rhs = [
this](
int row_side, EntityType row_type,
283 OP::nbIntegrationPts = OP::getGaussPts().size2();
285 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
287 OP::locF.resize(OP::nbRows,
false);
290 CHKERR this->iNtegrate(row_data);
293 CHKERR this->aSsemble(row_data);
297 auto do_work_lhs = [
this](
int row_side,
int col_side, EntityType row_type,
303 auto check_if_assemble_transpose = [&] {
305 if (OP::rowSide != OP::colSide || OP::rowType != OP::colType)
309 }
else if (OP::assembleTranspose) {
315 OP::rowSide = row_side;
316 OP::rowType = row_type;
317 OP::colSide = col_side;
318 OP::colType = col_type;
319 OP::nbCols = col_data.getIndices().size();
320 OP::locMat.resize(OP::nbRows, OP::nbCols,
false);
322 CHKERR this->iNtegrate(row_data, col_data);
324 CHKERR this->aSsemble(row_data, col_data, check_if_assemble_transpose());
328 switch (OP::opType) {
334 OP::nbIntegrationPts = OP::getGaussPts().size2();
335 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
340 for (
auto &bd : *brokenBaseSideData) {
343 if (!bd.getData().getNSharedPtr(bd.getData().getBase())) {
345 "base functions not set");
352 row_side, bd.getSide(),
355 row_type, bd.getType(),
358 row_data, bd.getData(),
368 for (
auto &bd : *brokenBaseSideData) {
369 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData(),
375 (std::string(
"wrong op type ") +
383template <
int FIELD_DIM, IntegrationType I,
typename OpBrokenBase>
386template <
int FIELD_DIM,
typename OpBrokenBase>
393 const std::string row_field,
394 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
395 boost::shared_ptr<double> beta_ptr,
const bool assmb_transpose,
396 const bool only_transpose, boost::shared_ptr<Range> ents_ptr =
nullptr)
397 :
OP(row_field, broken_base_side_data, assmb_transpose, only_transpose,
399 scalarBetaPtr(beta_ptr) {}
402 const std::string row_field,
403 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
404 double beta,
const bool assmb_transpose,
const bool only_transpose,
405 boost::shared_ptr<Range> ents_ptr =
nullptr)
408 assmb_transpose, only_transpose, ents_ptr) {}
416template <
int FIELD_DIM,
typename OpBase>
422 auto nb_row_dofs = row_data.
getIndices().size();
423 auto nb_col_dofs = col_data.
getIndices().size();
424 if (!nb_row_dofs || !nb_col_dofs)
430 auto t_w = this->getFTensor0IntegrationWeight();
431 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
432 size_t nb_base_functions = col_data.
getN().size2() / 3;
437 "number of dofs not divisible by field dimension");
442 "number of dofs exceeds number of base functions");
451 for (; cc != nb_col_dofs /
FIELD_DIM; cc++) {
453 for (
auto rr = 0; rr != nb_row_dofs /
FIELD_DIM; ++rr) {
455 (0.5 * t_w * t_row_base) * (t_normal(
J) * t_col_base(
J));
462 for (; cc < nb_base_functions; ++cc)
469 for (
auto rr = 0; rr != nb_row_dofs /
FIELD_DIM; ++rr) {
470 for (
auto cc = 0; cc != nb_col_dofs /
FIELD_DIM; ++cc) {
479 OP::locMat *= *scalarBetaPtr;
484template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
487template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
490template <
int FIELD_DIM,
typename OpBrokenBase>
497 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
498 boost::shared_ptr<MatrixDouble> lagrange_ptr,
499 boost::shared_ptr<double> beta_ptr,
500 boost::shared_ptr<Range> ents_ptr =
nullptr)
501 :
OP(broken_base_side_data, ents_ptr), scalarBetaPtr(beta_ptr),
502 lagrangePtr(lagrange_ptr) {}
505 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
506 boost::shared_ptr<MatrixDouble> lagrange_ptr,
double beta,
507 boost::shared_ptr<Range> ents_ptr =
nullptr)
519template <
int FIELD_DIM,
typename OpBase>
528 auto t_w = this->getFTensor0IntegrationWeight();
529 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
530 auto t_lagrange = getFTensor1FromMat<FIELD_DIM>(*lagrangePtr);
533 auto nb_base_functions = row_data.
getN().size2() / 3;
536 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
539 t_vec(
i) += (0.5 * t_w * (t_row_base(
J) * t_normal(
J))) * t_lagrange(
i);
543 for (; rr < nb_base_functions; ++rr)
551 OP::locF *= *scalarBetaPtr;
556template <
int FIELD_DIM,
typename OpBase>
563 const std::string row_field,
564 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
565 boost::shared_ptr<double> beta_ptr,
566 boost::shared_ptr<Range> ents_ptr =
nullptr)
567 :
OpBase(row_field, row_field,
OpBase::OPROW, ents_ptr),
568 brokenSideDataPtr(broken_side_data_ptr), scalarBetaPtr(beta_ptr) {}
571 const std::string row_field,
572 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
573 double beta, boost::shared_ptr<Range> ents_ptr =
nullptr)
584template <
int FIELD_DIM,
typename OpBase>
593 OP::locF.resize(row_data.
getIndices().size(),
false);
596 for (
auto &bd : *brokenSideDataPtr) {
597 auto t_w = this->getFTensor0IntegrationWeight();
598 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
600 auto t_flux = getFTensor2FromMat<FIELD_DIM, 3>(bd.getFlux());
601 auto nb_base_functions = row_data.
getN().size2() / 3;
602 auto sense = bd.getSense();
604 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
608 (sense * 0.5 * t_w) * t_row_base * t_normal(
J) * t_flux(
i,
J);
612 for (; rr < nb_base_functions; ++rr)
621 OP::locF *= *scalarBetaPtr;
626template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
629template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
638 const std::string row_field,
639 boost::shared_ptr<MatrixDouble> tangent1_diff_ptr,
640 boost::shared_ptr<MatrixDouble> tangent2_diff_ptr,
642 boost::shared_ptr<Range> ents_ptr =
nullptr)
643 :
OP(row_field, row_field,
OP::OPROW, ents_ptr),
650 if (!this->timeScalingFun.empty())
651 this->locF *= this->timeScalingFun(this->getFEMethod()->ts_t);
652 if (!this->feScalingFun.empty())
653 this->locF *= this->feScalingFun(this->getFEMethod());
655 auto *vec_ptr = this->locF.data().data();
656 const auto nb_dofs = row_data.
getIndices().size();
657 auto *ind_ptr = row_data.
getIndices().data().data();
662 std::vector<EntityHandle> ents(field_ents.size());
663 std::transform(field_ents.begin(), field_ents.end(), ents.begin(),
664 [](
const auto *fe) { return fe->getEnt(); });
665 if (field_ents.empty())
669 auto &moab = this->getMoab();
672 topo_values.data().data());
673 topo_values += this->locF;
675 topo_values.data().data());
686template <
typename OpBase>
693 const std::string row_field,
694 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
695 boost::shared_ptr<MatrixDouble> lamgrange_ptr,
696 boost::shared_ptr<MatrixDouble> tangent1_ptr,
697 boost::shared_ptr<MatrixDouble> tangent2_ptr,
699 Tag th, boost::shared_ptr<Range> ents_ptr =
nullptr)
700 :
OP(row_field, tangent1_ptr, tangent2_ptr, assemble_vec,
th, ents_ptr),
701 brokenSideDataPtr(broken_side_data_ptr), lagrangePtr(lamgrange_ptr),
702 scalarBetaPtr(beta_ptr) {}
711template <
typename OpBase>
728 &
m(0, 0), &
m(0, 1), &
m(0, 2));
731 for (
auto &bd : *brokenSideDataPtr) {
733 auto t_w = this->getFTensor0IntegrationWeight();
734 auto t_lagrange = getFTensor1FromMat<3>(*lagrangePtr);
735 auto t_adjoint_lambda = getFTensor2FromMat<3, 3>(bd.getVarFlux());
736 auto t_t1 = get_ftensor1(*this->tangent1DiffPtr);
737 auto t_t2 = get_ftensor1(*this->tangent2DiffPtr);
740 auto nb_base_functions = row_data.
getN().size2();
748 (
t_kd(
i,
m) * t_diff_base(N1));
750 (
t_kd(
k,
m) * t_diff_base(N0)) * t_t2(
i);
751 auto t_vec = getFTensor1FromPtr<3>(&*OP::locF.data().begin());
753 for (; rr != row_data.
getIndices().size() / 3; ++rr) {
754 t_vec(
m) += (0.5 * t_w) * (t_adjoint_lambda(
i,
J) * t_lagrange(
i)) *
759 for (; rr < nb_base_functions; ++rr)
770 OP::locF *= *scalarBetaPtr;
775template <
typename OpBase>
782 const std::string row_field,
783 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
784 boost::shared_ptr<MatrixDouble> adjoint_hybrid_ptr,
785 boost::shared_ptr<MatrixDouble> tangent1_ptr,
786 boost::shared_ptr<MatrixDouble> tangent2_ptr,
788 Tag th, boost::shared_ptr<Range> ents_ptr =
nullptr)
789 :
OP(row_field, tangent1_ptr, tangent2_ptr, assemble_vec,
th, ents_ptr),
790 brokenSideDataPtr(broken_side_data_ptr),
791 adjointHybridPtr(adjoint_hybrid_ptr), scalarBetaPtr(beta_ptr) {}
800template <
typename OpBase>
817 OP::locF.resize(row_data.
getIndices().size(),
false);
822 &
m(0, 0), &
m(0, 1), &
m(0, 2));
825 for (
auto &bd : *brokenSideDataPtr) {
826 auto t_w = this->getFTensor0IntegrationWeight();
827 auto t_t1 = get_ftensor1(*this->tangent1DiffPtr);
828 auto t_t2 = get_ftensor1(*this->tangent2DiffPtr);
829 auto t_flux = getFTensor2FromMat<3, 3>(bd.getFlux());
830 auto t_adjoint_lambda = getFTensor1FromMat<3>(*adjointHybridPtr);
832 auto nb_base_functions = row_data.
getN().size2();
833 auto sense = bd.getSense();
838 auto t_vec = getFTensor1FromPtr<3>(&*OP::locF.data().begin());
840 for (; rr != row_data.
getIndices().size() / 3; ++rr) {
842 t_normal_diff(
j,
m) =
844 (
t_kd(
i,
m) * t_diff_base(N1));
845 t_normal_diff(
j,
m) +=
849 t_vec(
m) += (sense * 0.5 * t_w) * (t_adjoint_lambda(
i) * t_flux(
i,
J)) *
855 for (; rr < nb_base_functions; ++rr)
866 OP::locF *= *scalarBetaPtr;
#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
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.
const VectorFieldEntities & getFieldEntities() const
Get field entities (const version)
VectorInt iNdices
Global indices on entity.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
virtual int getSense() const
Get entity sense for conforming approximation fields.
VectorFieldEntities fieldEntities
Field entities.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
std::array< std::array< boost::shared_ptr< MatrixDouble >, LASTBASE >, LastDerivative > baseFunctionsAndBaseDerivatives
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 > tangent2DiffPtr
boost::shared_ptr< MatrixDouble > tangent1DiffPtr
MoFEMErrorCode aSsemble(EntitiesFieldData::EntData &row_data)
SmartPetscObj< Vec > assembleVec
OpBrokenTopoBase(const std::string row_field, boost::shared_ptr< MatrixDouble > tangent1_diff_ptr, boost::shared_ptr< MatrixDouble > tangent2_diff_ptr, SmartPetscObj< Vec > assemble_vec, Tag th, boost::shared_ptr< Range > ents_ptr=nullptr)
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideData
OpGetBrokenBaseSideData(const std::string field_name, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data)
Element used to execute operators on side of the element.
const std::string sideFEName
boost::shared_ptr< E > sideFEPtr
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
OpSetFlux(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< MatrixDouble > flux_ptr)
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideData
boost::shared_ptr< MatrixDouble > fluxPtr
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
boost::shared_ptr< MatrixDouble > fluxVarPtr
OpSetVarFlux(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< MatrixDouble > flux_var_ptr)
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideData
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenSideDataPtr
OpTopoDerivativeBrokenSpaceConstrainDFluxImpl(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_side_data_ptr, boost::shared_ptr< MatrixDouble > lamgrange_ptr, boost::shared_ptr< MatrixDouble > tangent1_ptr, boost::shared_ptr< MatrixDouble > tangent2_ptr, boost::shared_ptr< double > beta_ptr, SmartPetscObj< Vec > assemble_vec, Tag th, boost::shared_ptr< Range > ents_ptr=nullptr)
boost::shared_ptr< double > scalarBetaPtr
boost::shared_ptr< MatrixDouble > lagrangePtr
boost::shared_ptr< double > scalarBetaPtr
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenSideDataPtr
boost::shared_ptr< MatrixDouble > adjointHybridPtr
OpTopoDerivativeBrokenSpaceConstrainDHybridImpl(const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_side_data_ptr, boost::shared_ptr< MatrixDouble > adjoint_hybrid_ptr, boost::shared_ptr< MatrixDouble > tangent1_ptr, boost::shared_ptr< MatrixDouble > tangent2_ptr, boost::shared_ptr< double > beta_ptr, SmartPetscObj< Vec > assemble_vec, Tag th, boost::shared_ptr< Range > ents_ptr=nullptr)
intrusive_ptr for managing petsc objects