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)
87 return integrateImpl(row_data, col_data, OpBase::getMeasure());
91template <
int FIELD_DIM,
typename OpBase>
103 return integrateImpl(row_data, col_data, OpBase::getMeasure());
107template <
int FIELD_DIM,
typename OpBase>
110 OpMassImpl(
const std::string row_field_name,
const std::string col_field_name,
112 boost::shared_ptr<Range> ents_ptr =
nullptr)
113 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
115 if (row_field_name == col_field_name)
125template <
typename OpBase>
127 OpMassImpl(
const std::string row_field_name,
const std::string col_field_name,
129 boost::shared_ptr<Range> ents_ptr =
nullptr)
130 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
132 if (row_field_name == col_field_name)
142template <
typename OpBase>
144 OpMassImpl(
const std::string row_field_name,
const std::string col_field_name,
146 boost::shared_ptr<Range> ents_ptr =
nullptr)
147 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
149 if (row_field_name == col_field_name)
159template <
int BASE_DIM,
int FIELD_DIM, IntegrationType I,
typename OpBase>
162template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
164 :
public OpMassImpl<1, FIELD_DIM, I, OpBase> {
166 const std::string col_field_name,
const double beta,
167 boost::shared_ptr<Range> ents_ptr =
nullptr)
169 row_field_name, col_field_name,
171 cacheLocMats(cache), scalarBeta(beta) {}
174 EntitiesFieldData::EntData &col_data);
185template <
int SPACE_DIM,
int S,
typename OpBase>
190 const std::string col_field_name,
191 boost::shared_ptr<MatrixDouble> mat_D,
192 boost::shared_ptr<Range> ents_ptr =
nullptr,
194 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
195 matD(mat_D), betaCoeff(beta) {
196 if (row_field_name == col_field_name)
201 boost::shared_ptr<MatrixDouble>
matD;
211template <
int SPACE_DIM,
int S,
typename OpBase>
216 const std::string col_field_name,
217 boost::shared_ptr<MatrixDouble> mat_D,
219 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL), matD(mat_D),
223 boost::shared_ptr<MatrixDouble>
matD;
233template <
int SPACE_DIM,
int S,
typename OpBase>
238 const std::string col_field_name,
239 boost::shared_ptr<MatrixDouble> mat_D,
240 boost::shared_ptr<Range> ents_ptr =
nullptr)
241 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL, ents_ptr),
243 if (row_field_name == col_field_name)
248 boost::shared_ptr<MatrixDouble>
matD;
253template <
int SPACE_DIM, IntegrationType I,
typename OpBase>
256template <
int SPACE_DIM,
typename OpBase>
259 const std::string col_field_name,
261 const bool assemble_transpose =
false,
262 const bool only_transpose =
false)
263 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
264 alphaConstant(alpha_fun) {
265 this->assembleTranspose = assemble_transpose;
266 this->onlyTranspose = only_transpose;
280template <
int SPACE_DIM,
typename OpBase, CoordinateTypes CoordSys>
285 const std::string col_field_name,
ConstantFun alpha_fun,
286 const bool assemble_transpose,
287 const bool only_transpose =
false)
288 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
289 alphaConstant(alpha_fun) {
290 this->assembleTranspose = assemble_transpose;
291 this->onlyTranspose = only_transpose;
295 const std::string col_field_name,
ConstantFun alpha_fun,
296 ScalarFun beta_fun,
const bool assemble_transpose,
297 const bool only_transpose =
false)
298 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
299 alphaConstant(alpha_fun), betaCoeff(beta_fun) {
300 this->assembleTranspose = assemble_transpose;
301 this->onlyTranspose = only_transpose;
318template <
int SPACE_DIM,
typename OpBase, CoordinateTypes COORDINATE_SYSTEM>
322 const std::string row_field_name,
const std::string col_field_name,
323 ScalarFun alpha_fun = [](
double,
double,
double)
constexpr {
return 1; },
324 const bool assemble_transpose =
false,
const bool only_transpose =
false)
325 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
326 alphaConstant(alpha_fun) {
327 this->assembleTranspose = assemble_transpose;
328 this->onlyTranspose = only_transpose;
343template <
int SPACE_DIM,
typename OpBase>
347 const std::string col_field_name,
349 const bool assemble_transpose =
false,
350 const bool only_transpose =
false)
351 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
352 alphaConstant(alpha_fun) {
353 this->assembleTranspose = assemble_transpose;
354 this->onlyTranspose = only_transpose;
364template <
int SPACE_DIM,
typename OpBase>
368 const std::string col_field_name,
370 const bool assemble_transpose =
false,
371 const bool only_transpose =
false)
372 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
373 alphaConstant(alpha_fun) {
374 this->assembleTranspose = assemble_transpose;
375 this->onlyTranspose = only_transpose;
385template <
int SPACE_DIM, IntegrationType I,
typename OpBase>
388template <
int SPACE_DIM,
typename OpBase>
392 const std::string col_field_name,
393 ConstantFun alpha_fun,
const bool assemble_transpose,
394 const bool only_transpose =
false)
395 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
396 alphaConstant(alpha_fun) {
397 this->assembleTranspose = assemble_transpose;
398 this->onlyTranspose = only_transpose;
402 const std::string col_field_name,
404 const bool assemble_transpose,
405 const bool only_transpose =
false)
406 :
OpBase(row_field_name, col_field_name,
OpBase::OPROWCOL),
407 alphaConstant(alpha_fun), betaCoeff(beta_coeff) {
408 this->assembleTranspose = assemble_transpose;
409 this->onlyTranspose = only_transpose;
421template <
int SPACE_DIM,
typename OpBase>
428 const double vol = OpBase::getMeasure();
430 auto t_w = OpBase::getFTensor0IntegrationWeight();
434 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
438 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
440 const double alpha = t_w * beta;
465template <
int FIELD_DIM,
int SPACE_DIM,
typename OpBase>
473 const double vol = OpBase::getMeasure();
475 auto t_w = OpBase::getFTensor0IntegrationWeight();
479 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
484 auto get_t_vec = [&](
const int rr) {
485 std::array<double *, FIELD_DIM> ptrs;
494 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
496 const double alpha = t_w * beta;
507 t_vec(
j) += alpha * (t_row_grad(
i) * t_col_grad(
i));
524template <
typename OpBase>
531 auto log_error = [&]() {
541 "Wrong number of base functions on rows %lu < %d",
547 "Wrong number of base functions on cols %lu < %d",
553 "Wrong number of integration points on rows %lu != %d",
559 "Wrong number of integration points on cols %lu != %d",
565 auto t_w = OpBase::getFTensor0IntegrationWeight();
569 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
572 const double beta = betaCoeff(t_coords(0), t_coords(1), t_coords(2));
574 const double alpha = t_w * beta;
584 *a_mat_ptr += alpha * (t_row_base * t_col_base);
601template <
int FIELD_DIM,
typename OpBase>
607 auto integrate = [&](
auto &mat) {
610 auto t_w = OpBase::getFTensor0IntegrationWeight();
614 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
617 auto get_t_vec = [&](
const int rr) {
618 std::array<double *, FIELD_DIM> ptrs;
620 ptrs[
i] = &mat(rr +
i,
i);
627 const double beta = betaCoeff(t_coords(0), t_coords(1), t_coords(2));
629 const double alpha = t_w * beta;
640 t_vec(
i) += alpha * (t_row_base * t_col_base);
660template <
int FIELD_DIM,
typename OpBase>
666 size_t nb_base_functions = row_data.
getN().size2() / 3;
668 const double vol = OpBase::getMeasure();
670 auto t_w = OpBase::getFTensor0IntegrationWeight();
674 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
677 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
679 const double alpha = t_w * beta;
689 (*a_mat_ptr) += alpha * (t_row_base(
i) * t_col_base(
i));
695 for (; rr < nb_base_functions; ++rr)
703template <
typename OpBase>
710 auto get_t_vec = [&](
const int rr) {
714 size_t nb_base_functions = row_data.
getN().size2() / 3;
716 const double vol = OpBase::getMeasure();
718 auto t_w = OpBase::getFTensor0IntegrationWeight();
722 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
725 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
727 const double alpha = t_w * beta;
733 auto t_vec = get_t_vec(2 * rr);
737 t_vec(
I) += alpha * (t_row_base(
k) * t_col_base(
k));
743 for (; rr < nb_base_functions; ++rr)
751template <
typename OpBase>
758 auto get_t_vec = [&](
const int rr) {
763 size_t nb_base_functions = row_data.
getN().size2() / 3;
765 const double vol = OpBase::getMeasure();
767 auto t_w = OpBase::getFTensor0IntegrationWeight();
771 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
774 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
776 const double alpha = t_w * beta;
782 auto t_vec = get_t_vec(3 * rr);
786 t_vec(
i) += alpha * (t_row_base(
k) * t_col_base(
k));
792 for (; rr < nb_base_functions; ++rr)
800template <
int FIELD_DIM, IntegrationType I,
typename OpBase>
806 const auto vol = this->getMeasure();
807 const auto row_type = this->rowType;
808 const auto col_type = this->colType;
809 auto &loc_mat = this->locMat;
811 auto p = std::make_pair(row_type, col_type);
813 if (cacheLocMats[p]) {
814 if (cacheLocMats[p]->size1() != loc_mat.size1() &&
815 cacheLocMats[p]->size2() != loc_mat.size2()) {
816 cacheLocMats[p]->resize(loc_mat.size1(), loc_mat.size2());
817 CHKERR this->integrateImpl(row_data, col_data, 1);
818 *(cacheLocMats[p]) = loc_mat;
820 loc_mat = *(cacheLocMats[p]);
822 loc_mat *= scalarBeta * this->getMeasure();
824 CHKERR this->integrateImpl(row_data, col_data,
825 scalarBeta * this->getMeasure());
830template <
int SPACE_DIM,
int S,
typename OpBase>
837 const size_t nb_row_dofs = row_data.
getIndices().size();
838 const size_t nb_col_dofs = col_data.
getIndices().size();
840 if (nb_row_dofs && nb_col_dofs) {
851 double vol = OpBase::getMeasure();
854 auto t_w = OpBase::getFTensor0IntegrationWeight();
861 auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, S>(*matD);
864 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
870 double a = t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
877 auto t_m = OpBase::template getLocMat<SPACE_DIM>(
SPACE_DIM * rr);
883 t_rowD(
l,
j,
k) = t_D(
i,
j,
k,
l) * (
a * t_row_diff_base(
i));
890 for (
int cc = 0; cc <= nb_cols; ++cc) {
893 t_m(
i,
j) += t_rowD(
i,
j,
k) * t_col_diff_base(
k);
918 auto t_m_rr = getFTensor2FromArray<SPACE_DIM, SPACE_DIM, SPACE_DIM>(
920 auto t_m_cc = getFTensor2FromArray<SPACE_DIM, SPACE_DIM>(
924 t_m_rr(
i,
j) = t_m_cc(
j,
i);
935template <
int SPACE_DIM,
int S,
typename OpBase>
949 auto &row_hessian = row_data.
getN(BaseDerivatives::SecondDerivative);
950 auto &col_hessian = col_data.
getN(BaseDerivatives::SecondDerivative);
955 "Wrong number of integration pts (%d != %d)",
958 if (row_hessian.size2() !=
961 "Wrong number of base functions (%d != %d)",
967 "Wrong number of base functions (%d < %d)", row_hessian.size2(),
972 "Wrong number of integration pts (%d != %d)",
977 "Wrong number of base functions (%d < %d)", col_hessian.size2(),
983 double vol = OpBase::getMeasure();
986 auto t_w = OpBase::getFTensor0IntegrationWeight();
988 auto t_row_diff2 = getFTensor2SymmetricLowerFromPtr<SPACE_DIM>(
989 &*row_hessian.data().begin());
993 auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, S>(*matD);
999 double a = t_w * vol;
1009 t_rowD(
k,
l) = t_D(
i,
j,
k,
l) * (
a * t_row_diff2(
i,
j));
1013 getFTensor2SymmetricLowerFromPtr<SPACE_DIM>(&col_hessian(gg, 0));
1019 *m_ptr += t_rowD(
i,
j) * t_col_diff2(
i,
j);
1044template <
int SPACE_DIM,
int S,
typename OpBase>
1051 const size_t nb_row_dofs = row_data.
getIndices().size();
1052 const size_t nb_col_dofs = col_data.
getIndices().size();
1054 if (nb_row_dofs && nb_col_dofs) {
1062 double vol = OpBase::getMeasure();
1065 auto t_w = OpBase::getFTensor0IntegrationWeight();
1072 getFTensor4FromMat<SPACE_DIM, SPACE_DIM, SPACE_DIM, SPACE_DIM, S>(
1076 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1082 double a = t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
1089 auto t_m = OpBase::template getLocMat<SPACE_DIM>(
SPACE_DIM * rr);
1093 t_row(
i,
k,
l) = t_D(
i,
j,
k,
l) * (
a * t_row_diff_base(
j));
1102 t_m(
i,
k) += t_row(
i,
k,
l) * t_col_diff_base(
l);
1128template <
int SPACE_DIM,
typename OpBase>
1134 auto t_w = this->getFTensor0IntegrationWeight();
1136 size_t nb_base_functions = row_data.
getN().size2() / 3;
1139 const double alpha_constant = alphaConstant();
1142 const double alpha = alpha_constant * this->getMeasure() * t_w;
1146 const double t_row_div_base = t_row_diff_base(
i,
i);
1149 this->locMat(rr, cc) += alpha * t_row_div_base * t_col_base;
1154 for (; rr < nb_base_functions; ++rr)
1163template <
int SPACE_DIM,
typename OpBase, CoordinateTypes CoordSys>
1170 auto t_w = this->getFTensor0IntegrationWeight();
1171 auto t_coords = this->getFTensor1CoordsAtGaussPts();
1173 size_t nb_base_functions = row_data.
getN().size2() / 3;
1176 const double alpha_constant = alphaConstant() * this->getMeasure();
1179 const double alpha =
1180 alpha_constant * t_w * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
1184 auto t_mat_diag = getFTensor1FromArrayDiag<SPACE_DIM, SPACE_DIM>(
1186 const double t_row_div_base = t_row_diff_base(
i,
i);
1190 t_mat_diag(
i) += alpha * t_row_div_base * t_col_base;
1192 t_mat_diag(
i) += t_row_base(0) * (alpha / t_coords(0)) * t_col_base;
1200 for (; rr < nb_base_functions; ++rr) {
1212template <
int SPACE_DIM,
typename OpBase, CoordinateTypes COORDINATE_SYSTEM>
1222 "Number of rows in matrix should be multiple of space dimensions");
1226 if constexpr (COORDINATE_SYSTEM ==
POLAR || COORDINATE_SYSTEM ==
SPHERICAL)
1228 "%s coordiante not implemented",
1231 auto t_w = this->getFTensor0IntegrationWeight();
1232 auto t_coords = this->getFTensor1CoordsAtGaussPts();
1233 size_t nb_base_functions_row = row_data.
getN().size2();
1235 const double vol = this->getMeasure();
1238 const double alpha =
1239 alphaConstant(t_coords(0), t_coords(1), t_coords(2)) * t_w * vol;
1242 auto t_m = getFTensor1FromPtr<SPACE_DIM>(
OpBase::locMat.data().data());
1245 if constexpr (COORDINATE_SYSTEM ==
CARTESIAN) {
1247 const double r_val = alpha * t_row_base;
1250 t_m(
i) += r_val * t_col_diff_base(
i);
1261 const double r_val = alpha * t_row_base;
1265 t_m(
i) += r_val * t_col_diff_base(
i);
1266 t_m(0) += (r_val / t_coords(0)) * t_col_base;
1275 for (; rr < nb_base_functions_row; ++rr)
1285template <
int SPACE_DIM,
typename OpBase>
1292 auto t_w = this->getFTensor0IntegrationWeight();
1294 size_t nb_base_functions = row_data.
getN().size2() / 3;
1296 auto &mat = this->locMat;
1297 const double alpha_constant = alphaConstant();
1300 const double alpha = alpha_constant * this->getMeasure() * t_w;
1306 mat(rr, cc) += alpha * t_row_base(
i) * t_col_diff_base(
i);
1311 for (; rr < nb_base_functions; ++rr)
1320template <
int SPACE_DIM,
typename OpBase>
1327 auto t_w = this->getFTensor0IntegrationWeight();
1329 size_t nb_base_functions = row_data.
getN().size2();
1332 auto get_t_vec = [&](
const int rr) {
1333 std::array<double *, SPACE_DIM> ptrs;
1339 const double alpha_constant = alphaConstant();
1342 const double alpha = alpha_constant * this->getMeasure() * t_w;
1349 t_vec(
i) += alpha * t_row_base * t_col_diff_base(
i);
1355 for (; rr < nb_base_functions; ++rr)
1364template <
int SPACE_DIM,
typename OpBase>
1370 auto t_w = this->getFTensor0IntegrationWeight();
1371 auto t_coords = this->getFTensor1CoordsAtGaussPts();
1373 size_t nb_base_functions = row_data.
getN().size2() / 3;
1375 const double alpha_constant = alphaConstant() * this->getMeasure();
1378 const double alpha =
1379 alpha_constant * betaCoeff(t_coords(0), t_coords(1), t_coords(2)) * t_w;
1383 auto t_mat_diag = getFTensor1FromArrayDiag<SPACE_DIM, SPACE_DIM>(
1388 t_mat_diag(
i) += alpha * t_row_base(
j) * t_col_diff_base(
j);
1395 for (; rr < nb_base_functions; ++rr)
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
double scalar_fun_one(const double, const double, const double)
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
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 dofs 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 &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)