13#ifndef __BILINEAR_FORMS_INTEGRATORS_IMPL_HPP__
14#define __BILINEAR_FORMS_INTEGRATORS_IMPL_HPP__
19 std::map<std::pair<int, int>, boost::shared_ptr<MatrixDouble>>;
25template <
int SPACE_DIM,
typename OpBase>
29 const std::string col_field_name,
31 boost::shared_ptr<Range> ents_ptr =
nullptr)
32 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
34 if (row_field_name == col_field_name)
44template <
int FIELD_DIM,
int SPACE_DIM,
typename OpBase>
48 const std::string col_field_name,
50 boost::shared_ptr<Range> ents_ptr =
nullptr)
51 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
53 if (row_field_name == col_field_name)
63template <
int BASE_DIM,
int FIELD_DIM, IntegrationType I,
typename OpBase>
66template <
typename OpBase>
70 const std::string row_field_name,
const std::string col_field_name,
71 ScalarFun beta = [](
double,
double,
double)
constexpr {
return 1; },
72 boost::shared_ptr<Range> ents_ptr =
nullptr,
73 boost::shared_ptr<MatrixDouble> cache_mat =
nullptr)
74 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
76 if (row_field_name == col_field_name)
89 return integrateImpl(row_data, col_data, OpBase::getMeasure());
94template <
int FIELD_DIM,
typename OpBase>
106 return integrateImpl(row_data, col_data, OpBase::getMeasure());
113template <
int FIELD_DIM,
typename OpBase>
116 OpMassImpl(
const std::string row_field_name,
const std::string col_field_name,
118 boost::shared_ptr<Range> ents_ptr =
nullptr)
119 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
121 if (row_field_name == col_field_name)
131template <
typename OpBase>
133 OpMassImpl(
const std::string row_field_name,
const std::string col_field_name,
135 boost::shared_ptr<Range> ents_ptr =
nullptr)
136 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
138 if (row_field_name == col_field_name)
148template <
typename OpBase>
150 OpMassImpl(
const std::string row_field_name,
const std::string col_field_name,
152 boost::shared_ptr<Range> ents_ptr =
nullptr)
153 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
155 if (row_field_name == col_field_name)
165template <
int BASE_DIM,
int FIELD_DIM, IntegrationType I,
typename OpBase>
168template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
170 :
public OpMassImpl<1, FIELD_DIM, I, OpBase> {
172 const std::string col_field_name,
const double beta,
173 boost::shared_ptr<Range> ents_ptr =
nullptr)
175 row_field_name, col_field_name,
177 cacheLocMats(cache), scalarBeta(beta) {}
180 EntitiesFieldData::EntData &col_data);
191template <
int SPACE_DIM,
int S,
typename OpBase>
196 const std::string col_field_name,
197 boost::shared_ptr<MatrixDouble> mat_D,
198 boost::shared_ptr<Range> ents_ptr =
nullptr,
200 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
201 matD(mat_D), betaCoeff(beta) {
202 if (row_field_name == col_field_name)
207 boost::shared_ptr<MatrixDouble>
matD;
217template <
int SPACE_DIM,
int S,
typename OpBase>
222 const std::string col_field_name,
223 boost::shared_ptr<MatrixDouble> mat_D,
225 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL), matD(mat_D),
229 boost::shared_ptr<MatrixDouble>
matD;
239template <
int SPACE_DIM,
int S,
typename OpBase>
244 const std::string col_field_name,
245 boost::shared_ptr<MatrixDouble> mat_D,
246 boost::shared_ptr<Range> ents_ptr =
nullptr)
247 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
249 if (row_field_name == col_field_name)
254 boost::shared_ptr<MatrixDouble>
matD;
259template <
int SPACE_DIM, IntegrationType I,
typename OpBase>
262template <
int SPACE_DIM,
typename OpBase>
265 const std::string col_field_name,
267 const bool assemble_transpose =
false,
268 const bool only_transpose =
false)
269 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
270 alphaConstant(alpha_fun) {
271 this->assembleTranspose = assemble_transpose;
272 this->onlyTranspose = only_transpose;
286template <
int SPACE_DIM,
typename OpBase, CoordinateTypes CoordSys>
291 const std::string col_field_name,
ConstantFun alpha_fun,
292 const bool assemble_transpose,
293 const bool only_transpose =
false)
294 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
295 alphaConstant(alpha_fun) {
296 this->assembleTranspose = assemble_transpose;
297 this->onlyTranspose = only_transpose;
301 const std::string col_field_name,
ConstantFun alpha_fun,
302 ScalarFun beta_fun,
const bool assemble_transpose,
303 const bool only_transpose =
false)
304 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
305 alphaConstant(alpha_fun), betaCoeff(beta_fun) {
306 this->assembleTranspose = assemble_transpose;
307 this->onlyTranspose = only_transpose;
324template <
int SPACE_DIM,
typename OpBase, CoordinateTypes COORDINATE_SYSTEM>
328 const std::string row_field_name,
const std::string col_field_name,
329 ScalarFun alpha_fun = [](
double,
double,
double)
constexpr {
return 1; },
330 const bool assemble_transpose =
false,
const bool only_transpose =
false)
331 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
332 alphaConstant(alpha_fun) {
333 this->assembleTranspose = assemble_transpose;
334 this->onlyTranspose = only_transpose;
349template <
int SPACE_DIM,
typename OpBase>
353 const std::string col_field_name,
355 const bool assemble_transpose =
false,
356 const bool only_transpose =
false)
357 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
358 alphaConstant(alpha_fun) {
359 this->assembleTranspose = assemble_transpose;
360 this->onlyTranspose = only_transpose;
370template <
int SPACE_DIM,
typename OpBase>
374 const std::string col_field_name,
376 const bool assemble_transpose =
false,
377 const bool only_transpose =
false)
378 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
379 alphaConstant(alpha_fun) {
380 this->assembleTranspose = assemble_transpose;
381 this->onlyTranspose = only_transpose;
391template <
int SPACE_DIM, IntegrationType I,
typename OpBase>
394template <
int SPACE_DIM,
typename OpBase>
398 const std::string col_field_name,
399 ConstantFun alpha_fun,
const bool assemble_transpose,
400 const bool only_transpose =
false)
401 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
402 alphaConstant(alpha_fun) {
403 this->assembleTranspose = assemble_transpose;
404 this->onlyTranspose = only_transpose;
408 const std::string col_field_name,
410 const bool assemble_transpose,
411 const bool only_transpose =
false)
412 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
413 alphaConstant(alpha_fun), betaCoeff(beta_coeff) {
414 this->assembleTranspose = assemble_transpose;
415 this->onlyTranspose = only_transpose;
427template <
int SPACE_DIM,
typename OpBase>
434 const double vol = OpBase::getMeasure();
436 auto t_w = OpBase::getFTensor0IntegrationWeight();
440 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
444 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
446 const double alpha = t_w * beta;
471template <
int FIELD_DIM,
int SPACE_DIM,
typename OpBase>
479 const double vol = OpBase::getMeasure();
481 auto t_w = OpBase::getFTensor0IntegrationWeight();
485 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
490 auto get_t_vec = [&](
const int rr) {
491 std::array<double *, FIELD_DIM> ptrs;
500 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
502 const double alpha = t_w * beta;
513 t_vec(
j) += alpha * (t_row_grad(
i) * t_col_grad(
i));
530template <
typename OpBase>
537 auto log_error = [&]() {
547 "Wrong number of base functions on rows %lu < %d",
553 "Wrong number of base functions on cols %lu < %d",
559 "Wrong number of integration points on rows %lu != %d",
565 "Wrong number of integration points on cols %lu != %d",
571 auto t_w = OpBase::getFTensor0IntegrationWeight();
575 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
578 const double beta = betaCoeff(t_coords(0), t_coords(1), t_coords(2));
580 const double alpha = t_w * beta;
590 *a_mat_ptr += alpha * (t_row_base * t_col_base);
607template <
int FIELD_DIM,
typename OpBase>
613 auto integrate = [&](
auto &mat) {
616 auto t_w = OpBase::getFTensor0IntegrationWeight();
620 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
623 auto get_t_vec = [&](
const int rr) {
624 std::array<double *, FIELD_DIM> ptrs;
626 ptrs[
i] = &mat(rr +
i,
i);
633 const double beta = betaCoeff(t_coords(0), t_coords(1), t_coords(2));
635 const double alpha = t_w * beta;
646 t_vec(
i) += alpha * (t_row_base * t_col_base);
666template <
int FIELD_DIM,
typename OpBase>
672 size_t nb_base_functions = row_data.
getN().size2() / 3;
674 const double vol = OpBase::getMeasure();
676 auto t_w = OpBase::getFTensor0IntegrationWeight();
680 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
683 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
685 const double alpha = t_w * beta;
695 (*a_mat_ptr) += alpha * (t_row_base(
i) * t_col_base(
i));
701 for (; rr < nb_base_functions; ++rr)
709template <
typename OpBase>
716 auto get_t_vec = [&](
const int rr) {
720 size_t nb_base_functions = row_data.
getN().size2() / 3;
722 const double vol = OpBase::getMeasure();
724 auto t_w = OpBase::getFTensor0IntegrationWeight();
728 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
731 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
733 const double alpha = t_w * beta;
739 auto t_vec = get_t_vec(2 * rr);
743 t_vec(
I) += alpha * (t_row_base(
k) * t_col_base(
k));
749 for (; rr < nb_base_functions; ++rr)
757template <
typename OpBase>
764 auto get_t_vec = [&](
const int rr) {
769 size_t nb_base_functions = row_data.
getN().size2() / 3;
771 const double vol = OpBase::getMeasure();
773 auto t_w = OpBase::getFTensor0IntegrationWeight();
777 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
780 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
782 const double alpha = t_w * beta;
788 auto t_vec = get_t_vec(3 * rr);
792 t_vec(
i) += alpha * (t_row_base(
k) * t_col_base(
k));
798 for (; rr < nb_base_functions; ++rr)
806template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
812 const auto vol = this->getMeasure();
813 const auto row_type = this->rowType;
814 const auto col_type = this->colType;
815 auto &loc_mat = this->locMat;
817 auto p = std::make_pair(row_type, col_type);
819 if (cacheLocMats[p]) {
820 if (cacheLocMats[p]->size1() != loc_mat.size1() &&
821 cacheLocMats[p]->size2() != loc_mat.size2()) {
822 cacheLocMats[p]->resize(loc_mat.size1(), loc_mat.size2());
823 CHKERR this->integrateImpl(row_data, col_data, 1);
824 *(cacheLocMats[p]) = loc_mat;
826 loc_mat = *(cacheLocMats[p]);
828 loc_mat *= scalarBeta * this->getMeasure();
830 CHKERR this->integrateImpl(row_data, col_data,
831 scalarBeta * this->getMeasure());
836template <
int SPACE_DIM,
int S,
typename OpBase>
843 const size_t nb_row_dofs = row_data.
getIndices().size();
844 const size_t nb_col_dofs = col_data.
getIndices().size();
846 if (nb_row_dofs && nb_col_dofs) {
857 double vol = OpBase::getMeasure();
860 auto t_w = OpBase::getFTensor0IntegrationWeight();
867 auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, S>(*matD);
870 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
876 double a = t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
883 auto t_m = OpBase::template getLocMat<SPACE_DIM>(
SPACE_DIM * rr);
889 t_rowD(
l,
j,
k) = t_D(
i,
j,
k,
l) * (
a * t_row_diff_base(
i));
896 for (
int cc = 0; cc <= nb_cols; ++cc) {
899 t_m(
i,
j) += t_rowD(
i,
j,
k) * t_col_diff_base(
k);
924 auto t_m_rr = getFTensor2FromArray<SPACE_DIM, SPACE_DIM, SPACE_DIM>(
926 auto t_m_cc = getFTensor2FromArray<SPACE_DIM, SPACE_DIM>(
930 t_m_rr(
i,
j) = t_m_cc(
j,
i);
941template <
int SPACE_DIM,
int S,
typename OpBase>
955 auto &row_hessian = row_data.
getN(BaseDerivatives::SecondDerivative);
956 auto &col_hessian = col_data.
getN(BaseDerivatives::SecondDerivative);
961 "Wrong number of integration pts (%d != %d)",
964 if (row_hessian.size2() !=
967 "Wrong number of base functions (%d != %d)",
973 "Wrong number of base functions (%d < %d)", row_hessian.size2(),
978 "Wrong number of integration pts (%d != %d)",
983 "Wrong number of base functions (%d < %d)", col_hessian.size2(),
989 double vol = OpBase::getMeasure();
992 auto t_w = OpBase::getFTensor0IntegrationWeight();
994 auto t_row_diff2 = getFTensor2SymmetricLowerFromPtr<SPACE_DIM>(
995 &*row_hessian.data().begin());
999 auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, S>(*matD);
1005 double a = t_w * vol;
1015 t_rowD(
k,
l) = t_D(
i,
j,
k,
l) * (
a * t_row_diff2(
i,
j));
1019 getFTensor2SymmetricLowerFromPtr<SPACE_DIM>(&col_hessian(gg, 0));
1025 *m_ptr += t_rowD(
i,
j) * t_col_diff2(
i,
j);
1050template <
int SPACE_DIM,
int S,
typename OpBase>
1057 const size_t nb_row_dofs = row_data.
getIndices().size();
1058 const size_t nb_col_dofs = col_data.
getIndices().size();
1060 if (nb_row_dofs && nb_col_dofs) {
1068 double vol = OpBase::getMeasure();
1071 auto t_w = OpBase::getFTensor0IntegrationWeight();
1078 getFTensor4FromMat<SPACE_DIM, SPACE_DIM, SPACE_DIM, SPACE_DIM, S>(
1082 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1088 double a = t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
1095 auto t_m = OpBase::template getLocMat<SPACE_DIM>(
SPACE_DIM * rr);
1099 t_row(
i,
k,
l) = t_D(
i,
j,
k,
l) * (
a * t_row_diff_base(
j));
1108 t_m(
i,
k) += t_row(
i,
k,
l) * t_col_diff_base(
l);
1134template <
int SPACE_DIM,
typename OpBase>
1140 auto t_w = this->getFTensor0IntegrationWeight();
1142 size_t nb_base_functions = row_data.
getN().size2() / 3;
1145 const double alpha_constant = alphaConstant();
1148 const double alpha = alpha_constant * this->getMeasure() * t_w;
1152 const double t_row_div_base = t_row_diff_base(
i,
i);
1155 this->locMat(rr, cc) += alpha * t_row_div_base * t_col_base;
1160 for (; rr < nb_base_functions; ++rr)
1169template <
int SPACE_DIM,
typename OpBase, CoordinateTypes CoordSys>
1176 auto t_w = this->getFTensor0IntegrationWeight();
1177 auto t_coords = this->getFTensor1CoordsAtGaussPts();
1179 size_t nb_base_functions = row_data.
getN().size2() / 3;
1182 const double alpha_constant = alphaConstant() * this->getMeasure();
1185 const double alpha =
1186 alpha_constant * t_w * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
1190 auto t_mat_diag = getFTensor1FromArrayDiag<SPACE_DIM, SPACE_DIM>(
1192 const double t_row_div_base = t_row_diff_base(
i,
i);
1196 t_mat_diag(
i) += alpha * t_row_div_base * t_col_base;
1198 t_mat_diag(
i) += t_row_base(0) * (alpha / t_coords(0)) * t_col_base;
1206 for (; rr < nb_base_functions; ++rr) {
1218template <
int SPACE_DIM,
typename OpBase, CoordinateTypes COORDINATE_SYSTEM>
1228 "Number of rows in matrix should be multiple of space dimensions");
1232 if constexpr (COORDINATE_SYSTEM ==
POLAR || COORDINATE_SYSTEM ==
SPHERICAL)
1234 "%s coordiante not implemented",
1237 auto t_w = this->getFTensor0IntegrationWeight();
1238 auto t_coords = this->getFTensor1CoordsAtGaussPts();
1239 size_t nb_base_functions_row = row_data.
getN().size2();
1241 const double vol = this->getMeasure();
1244 const double alpha =
1245 alphaConstant(t_coords(0), t_coords(1), t_coords(2)) * t_w * vol;
1248 auto t_m = getFTensor1FromPtr<SPACE_DIM>(
OpBase::locMat.data().data());
1251 if constexpr (COORDINATE_SYSTEM ==
CARTESIAN) {
1253 const double r_val = alpha * t_row_base;
1256 t_m(
i) += r_val * t_col_diff_base(
i);
1267 const double r_val = alpha * t_row_base;
1271 t_m(
i) += r_val * t_col_diff_base(
i);
1272 t_m(0) += (r_val / t_coords(0)) * t_col_base;
1281 for (; rr < nb_base_functions_row; ++rr)
1291template <
int SPACE_DIM,
typename OpBase>
1298 auto t_w = this->getFTensor0IntegrationWeight();
1300 size_t nb_base_functions = row_data.
getN().size2() / 3;
1302 auto &mat = this->locMat;
1303 const double alpha_constant = alphaConstant();
1306 const double alpha = alpha_constant * this->getMeasure() * t_w;
1312 mat(rr, cc) += alpha * t_row_base(
i) * t_col_diff_base(
i);
1317 for (; rr < nb_base_functions; ++rr)
1326template <
int SPACE_DIM,
typename OpBase>
1333 auto t_w = this->getFTensor0IntegrationWeight();
1335 size_t nb_base_functions = row_data.
getN().size2();
1338 auto get_t_vec = [&](
const int rr) {
1339 std::array<double *, SPACE_DIM> ptrs;
1345 const double alpha_constant = alphaConstant();
1348 const double alpha = alpha_constant * this->getMeasure() * t_w;
1355 t_vec(
i) += alpha * t_row_base * t_col_diff_base(
i);
1361 for (; rr < nb_base_functions; ++rr)
1370template <
int SPACE_DIM,
typename OpBase>
1376 auto t_w = this->getFTensor0IntegrationWeight();
1377 auto t_coords = this->getFTensor1CoordsAtGaussPts();
1379 size_t nb_base_functions = row_data.
getN().size2() / 3;
1381 const double alpha_constant = alphaConstant() * this->getMeasure();
1384 const double alpha =
1385 alpha_constant * betaCoeff(t_coords(0), t_coords(1), t_coords(2)) * t_w;
1389 auto t_mat_diag = getFTensor1FromArrayDiag<SPACE_DIM, SPACE_DIM>(
1394 t_mat_diag(
i) += alpha * t_row_base(
j) * t_col_diff_base(
j);
1401 for (; rr < nb_base_functions; ++rr)
constexpr int SPACE_DIM
[Define dimension]
static const char *const CoordinateTypesNames[]
Coordinate system names.
#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 ...
@ MOFEM_DATA_INCONSISTENCY
CoordinateTypes
Coodinate system.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MOFEM_LOG(channel, severity)
Log.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
std::map< std::pair< int, int >, boost::shared_ptr< MatrixDouble > > CacheMatsTypeType
boost::function< double()> ConstantFun
Constant function type.
constexpr IntegrationType I
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
const VectorFieldEntities & getFieldEntities() const
Get field entities (const version)
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
int rowSide
row side number
int colSide
column side number
int nbRows
number of dofs on rows
EntityType colType
column type
int nbIntegrationPts
number of integration points
MatrixDouble locMat
local entity block matrix
int nbCols
number if dof on column
int nbRowBaseFunctions
number or row base functions
EntityType rowType
row type
OpGradGradImpl(const std::string row_field_name, const std::string col_field_name, ScalarFun beta=scalar_fun_one, boost::shared_ptr< Range > ents_ptr=nullptr)
FTensor::Index< 'i', SPACE_DIM > i
summit Index
OpGradGradImpl(const std::string row_field_name, const std::string col_field_name, ScalarFun beta=scalar_fun_one, boost::shared_ptr< Range > ents_ptr=nullptr)
FTensor::Index< 'i', SPACE_DIM > i
summit Index
boost::shared_ptr< MatrixDouble > matD
OpGradGradSymTensorGradGradImpl(const std::string row_field_name, const std::string col_field_name, boost::shared_ptr< MatrixDouble > mat_D, boost::shared_ptr< Range > ents_ptr=nullptr)
FTensor::Index< 'i', SPACE_DIM > i
summit Index
FTensor::Index< 'i', SPACE_DIM > i
summit Index
boost::shared_ptr< MatrixDouble > matD
OpGradSymTensorGradImpl(const std::string row_field_name, const std::string col_field_name, boost::shared_ptr< MatrixDouble > mat_D, boost::shared_ptr< Range > ents_ptr=nullptr, ScalarFun beta=scalar_fun_one)
boost::shared_ptr< MatrixDouble > matD
FTensor::Index< 'i', SPACE_DIM > i
summit Index
OpGradTensorGradImpl(const std::string row_field_name, const std::string col_field_name, boost::shared_ptr< MatrixDouble > mat_D, ScalarFun beta=scalar_fun_one)
CacheMatsTypeType cacheLocMats
OpMassCacheImpl(CacheMatsTypeType cache, const std::string row_field_name, const std::string col_field_name, const double beta, boost::shared_ptr< Range > ents_ptr=nullptr)
OpMassImpl(const std::string row_field_name, const std::string col_field_name, ScalarFun beta=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr, boost::shared_ptr< MatrixDouble > cache_mat=nullptr)
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Integrate grad-grad operator.
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &data)
Class dedicated to integrate operator.
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Integrate grad-grad operator.
OpMassImpl(const std::string row_field_name, const std::string col_field_name, ScalarFun beta=scalar_fun_one, boost::shared_ptr< Range > ents_ptr=nullptr)
OpMassImpl(const std::string row_field_name, const std::string col_field_name, ScalarFun beta=scalar_fun_one, boost::shared_ptr< Range > ents_ptr=nullptr)
OpMassImpl(const std::string row_field_name, const std::string col_field_name, ScalarFun beta=scalar_fun_one, boost::shared_ptr< Range > ents_ptr=nullptr)
ConstantFun alphaConstant
FTensor::Index< 'i', SPACE_DIM > i
summit Index
OpMixDivTimesScalarImpl(const std::string row_field_name, const std::string col_field_name, ConstantFun alpha_fun, const bool assemble_transpose=false, const bool only_transpose=false)
OpMixDivTimesVecImpl(const std::string row_field_name, const std::string col_field_name, ConstantFun alpha_fun, const bool assemble_transpose, const bool only_transpose=false)
FTensor::Index< 'i', SPACE_DIM > i
summit Index
OpMixDivTimesVecImpl(const std::string row_field_name, const std::string col_field_name, ConstantFun alpha_fun, ScalarFun beta_fun, const bool assemble_transpose, const bool only_transpose=false)
OpMixScalarTimesDivImpl(const std::string row_field_name, const std::string col_field_name, ScalarFun alpha_fun=[](double, double, double) constexpr { return 1;}, const bool assemble_transpose=false, const bool only_transpose=false)
FTensor::Index< 'i', SPACE_DIM > i
summit Index
FTensor::Index< 'i', SPACE_DIM > i
summit Index
OpMixTensorTimesGradImpl(const std::string row_field_name, const std::string col_field_name, ConstantFun alpha_fun, const bool assemble_transpose, const bool only_transpose=false)
OpMixTensorTimesGradImpl(const std::string row_field_name, const std::string col_field_name, ConstantFun alpha_fun, ScalarFun beta_coeff, const bool assemble_transpose, const bool only_transpose=false)
FTensor::Index< 'j', SPACE_DIM > j
summit Index
ConstantFun alphaConstant
FTensor::Index< 'i', SPACE_DIM > i
summit Index
OpMixVectorTimesGradImpl(const std::string row_field_name, const std::string col_field_name, ConstantFun alpha_fun, const bool assemble_transpose=false, const bool only_transpose=false)
FTensor::Index< 'i', SPACE_DIM > i
summit Index
ConstantFun alphaConstant
OpMixVectorTimesGradImpl(const std::string row_field_name, const std::string col_field_name, ConstantFun alpha_fun, const bool assemble_transpose=false, const bool only_transpose=false)