10#ifndef __FORMSBROKENSPACECONSTRAINTIMPL_HPP__
11#define __FORMSBROKENSPACECONSTRAINTIMPL_HPP__
27 prev_side_fe_ptr->numeredEntFiniteElementPtr->getFEUId();
40 prev_side_fe_ptr->numeredEntFiniteElementPtr;
45 "sideFEName is different");
74 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data);
83template <
typename OpBase>
86 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data)
88 brokenBaseSideData(broken_base_side_data) {}
90template <
typename OpBase>
99 brokenBaseSideData->resize(OP::getLoopSize());
101 const auto n_in_the_loop = OP::getNinTheLoop();
102 const auto face_sense = OP::getSkeletonSense();
105 if (face_sense != -1 && face_sense != 1)
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;
120 side_data.getData().fieldData = row_data.
fieldData;
123 auto set_base = [&](
auto &side_data) {
124 auto base = side_data.getData().getBase();
125 for (
auto dd = 0; dd != BaseDerivatives::LastDerivative; ++dd) {
127 side_data.getData().baseFunctionsAndBaseDerivatives[dd][base] =
128 boost::make_shared<MatrixDouble>(*base_ptr);
133 set_data((*brokenBaseSideData)[n_in_the_loop]);
134 set_base((*brokenBaseSideData)[n_in_the_loop]);
144 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
145 boost::shared_ptr<MatrixDouble> flux_ptr);
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),
162template <
typename OpBase>
166 auto swap_flux = [&](
auto &side_data) { side_data.getFlux().swap(*fluxPtr); };
167 swap_flux((*brokenBaseSideData)[OP::getNinTheLoop()]);
176 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
177 boost::shared_ptr<Range> ents_ptr =
nullptr)
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),
202template <
typename OpBase>
209 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
214 if (!brokenBaseSideData) {
219 auto do_work_rhs = [
this](
int row_side, EntityType row_type,
227 OP::nbIntegrationPts = OP::getGaussPts().size2();
229 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
231 OP::locF.resize(OP::nbRows,
false);
234 CHKERR this->iNtegrate(row_data);
237 CHKERR this->aSsemble(row_data);
241 auto do_work_lhs = [
this](
int row_side,
int col_side, EntityType row_type,
247 auto check_if_assemble_transpose = [&] {
249 if (OP::rowSide != OP::colSide || OP::rowType != OP::colType)
253 }
else if (OP::assembleTranspose) {
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);
266 CHKERR this->iNtegrate(row_data, col_data);
268 CHKERR this->aSsemble(row_data, col_data, check_if_assemble_transpose());
272 switch (OP::opType) {
278 OP::nbIntegrationPts = OP::getGaussPts().size2();
279 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
284 for (
auto &bd : *brokenBaseSideData) {
287 if (bd.getData().getSpace() !=
HDIV && bd.getData().getSpace() !=
HCURL) {
289 (std::string(
"Expect Hdiv or Hcurl space but received ") +
293 if (!bd.getData().getNSharedPtr(bd.getData().getBase())) {
295 "base functions not set");
302 row_side, bd.getSide(),
305 row_type, bd.getType(),
308 row_data, bd.getData(),
318 for (
auto &bd : *brokenBaseSideData) {
319 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData(),
325 (std::string(
"wrong op type ") +
333template <
int FIELD_DIM, IntegrationType I,
typename OpBrokenBase>
336template <
int FIELD_DIM,
typename OpBrokenBase>
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,
349 scalarBetaPtr(beta_ptr) {}
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)
357 boost::make_shared<
double>(beta),
358 assmb_transpose, only_transpose, ents_ptr) {}
366template <
int FIELD_DIM,
typename OpBase>
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)
380 auto t_w = this->getFTensor0IntegrationWeight();
381 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
382 size_t nb_base_functions = col_data.
getN().size2() / 3;
387 "number of dofs not divisible by field dimension");
392 "number of dofs exceeds number of base functions");
401 for (; cc != nb_col_dofs /
FIELD_DIM; cc++) {
403 for (
auto rr = 0; rr != nb_row_dofs /
FIELD_DIM; ++rr) {
405 (t_w * t_row_base) * (t_normal(
J) * t_col_base(
J));
412 for (; cc < nb_base_functions; ++cc)
419 for (
auto rr = 0; rr != nb_row_dofs /
FIELD_DIM; ++rr) {
420 for (
auto cc = 0; cc != nb_col_dofs /
FIELD_DIM; ++cc) {
429 OP::locMat *= *scalarBetaPtr;
434template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
437template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
440template <
int FIELD_DIM,
typename OpBrokenBase>
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) {}
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)
459 boost::make_shared<
double>(beta),
469template <
int FIELD_DIM,
typename OpBase>
478 auto t_w = this->getFTensor0IntegrationWeight();
479 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
480 auto t_lagrange = getFTensor1FromMat<FIELD_DIM>(*lagrangePtr);
483 auto nb_base_functions = row_data.
getN().size2() / 3;
486 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
489 t_vec(
i) += (t_w * (t_row_base(
J) * t_normal(
J))) * t_lagrange(
i);
493 for (; rr < nb_base_functions; ++rr)
501 OP::locF *= *scalarBetaPtr;
506template <
int FIELD_DIM,
typename OpBase>
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) {}
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)
525 boost::make_shared<
double>(beta),
534template <
int FIELD_DIM,
typename OpBase>
543 OP::locF.resize(row_data.
getIndices().size(),
false);
546 for (
auto &bd : *brokenSideDataPtr) {
547 auto t_w = this->getFTensor0IntegrationWeight();
548 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
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();
554 auto t_vec = getFTensor1FromPtr<FIELD_DIM>(&*OP::locF.data().begin());
557 t_vec(
i) += (sense * t_w) * t_row_base * t_normal(
J) * t_flux(
i,
J);
561 for (; rr < nb_base_functions; ++rr)
570 OP::locF *= *scalarBetaPtr;
#define FTENSOR_INDEX(DIM, I)
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HCURL
field with continuous tangents
@ HDIV
field with continuous normal traction
static const char *const FieldSpaceNames[]
#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
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
constexpr auto field_name
OpBaseImpl< PETSC, EdgeEleOp > OpBase
EntitiesFieldData::EntData entData
Data on single entity (This is passed as argument to DataOperator::doWork)
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, 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
VectorDofs dOfs
DoFs on entity.
const VectorInt & getIndices() const
Get global indices of dofs 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
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.
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)
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