6#ifndef __USER_DATA_OPERATORS_HPP__
7 #define __USER_DATA_OPERATORS_HPP__
22template <
class T,
class A>
28 boost::shared_ptr<ublas::vector<T, A>> data_ptr,
29 const EntityType zero_type = MBVERTEX)
38 boost::shared_ptr<ublas::vector<T, A>> data_ptr,
57 boost::shared_ptr<ublas::vector<T, A>>
dataPtr;
67template <
class T,
class A>
100 if (type ==
zeroType || vec.size() != nb_gauss_pts) {
101 vec.resize(nb_gauss_pts,
false);
114 for (
int i = 0;
i != local_indices.size(); ++
i)
115 if (local_indices[
i] != -1)
123 const size_t nb_base_functions = data.
getN().size2();
126 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
129 for (; bb != nb_dofs; ++bb) {
130 values_at_gauss_pts += field_data * base_function;
135 for (; bb < nb_base_functions; ++bb)
137 ++values_at_gauss_pts;
155template <PetscData::DataContext CTX>
160 const std::string
field_name, boost::shared_ptr<VectorDouble> data_ptr,
161 const EntityType zero_at_type = MBVERTEX)
176 if (type ==
zeroAtType || vec.size() != nb_gauss_pts) {
177 vec.resize(nb_gauss_pts,
false);
182 const size_t nb_dofs = local_indices.size();
187 auto get_array = [&](
const auto ctx,
auto vec) {
193 <<
"In this case field degrees of freedom are read from vector. "
194 "That usually happens when time solver is used, and acces to "
195 "first or second rates is needed. You probably not set ts_u, "
196 "ts_u_t, or ts_u_tt and associated data structure, i.e. "
197 "data_ctx to CTX_SET_X, CTX_SET_X_T, or CTX_SET_X_TT "
202 CHKERR VecGetArrayRead(vec, &array);
206 auto restore_array = [&](
auto vec) {
207 return VecRestoreArrayRead(vec, &array);
222 "That case is not implemented");
225 std::array<double, MAX_DOFS_ON_ENTITY> dot_dofs_vector;
226 for (
int i = 0;
i != local_indices.size(); ++
i)
227 if (local_indices[
i] != -1)
228 dot_dofs_vector[
i] = array[local_indices[
i]];
230 dot_dofs_vector[
i] = 0;
244 "That case is not implemented");
247 const size_t nb_base_functions = data.
getN().size2();
251 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
253 for (; bb != nb_dofs; ++bb) {
254 values_at_gauss_pts += dot_dofs_vector[bb] * base_function;
259 for (; bb < nb_base_functions; ++bb)
261 ++values_at_gauss_pts;
292template <
int Tensor_Dim,
class T,
class L,
class A>
298 boost::shared_ptr<ublas::matrix<T, L, A>> data_ptr,
299 const EntityType zero_type = MBVERTEX)
318 boost::shared_ptr<ublas::matrix<T, L, A>>
dataPtr;
322template <
int Tensor_Dim,
class T,
class L,
class A>
328 "Not implemented for T = %s and dim = %d",
338template <
int Tensor_Dim>
344 boost::shared_ptr<MatrixDouble> data_ptr,
345 const EntityType zero_type = MBVERTEX)
354 boost::shared_ptr<MatrixDouble> data_ptr,
356 const EntityType zero_type = MBVERTEX)
378template <
int Tensor_Dim>
380 Tensor_Dim,
double, ublas::row_major,
385 const size_t nb_gauss_pts = getGaussPts().size2();
386 auto &mat = *dataPtr;
387 if (type == zeroType || mat.size2() != nb_gauss_pts) {
388 mat.resize(Tensor_Dim, nb_gauss_pts,
false);
395 if (dataVec.use_count()) {
396 dotVector.resize(nb_dofs,
false);
398 CHKERR VecGetArrayRead(dataVec, &array);
400 for (
int i = 0;
i != local_indices.size(); ++
i)
401 if (local_indices[
i] != -1)
402 dotVector[
i] = array[local_indices[
i]];
405 CHKERR VecRestoreArrayRead(dataVec, &array);
410 const size_t nb_base_functions = data.
getN().size2();
412 auto values_at_gauss_pts = getFTensor1FromMat<Tensor_Dim>(mat);
414 const size_t size = nb_dofs / Tensor_Dim;
415 if (nb_dofs % Tensor_Dim) {
417 "Nb. of DOFs is inconsistent with Tensor_Dim");
419 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
423 if (field_data.l2() != field_data.l2()) {
424 MOFEM_LOG(
"SELF", Sev::error) <<
"field data: " << field_data;
426 "Wrong number in coefficients");
431 for (; bb != size; ++bb) {
435 if (base_function != base_function) {
436 MOFEM_LOG(
"SELF", Sev::error) <<
"base function: " << base_function;
438 "Wrong number number in base functions");
443 values_at_gauss_pts(
I) += field_data(
I) * base_function;
449 for (; bb < nb_base_functions; ++bb)
451 ++values_at_gauss_pts;
455 if (dataVec.use_count()) {
467template <
int Tensor_Dim>
470 Tensor_Dim, double, ublas::row_major, DoubleAllocator> {
473 Tensor_Dim,
double, ublas::row_major,
487template <
int Tensor_Dim, CoordinateTypes COORDINATE_SYSTEM = CARTESIAN>
492 const std::string
field_name, boost::shared_ptr<VectorDouble> data_ptr,
493 const EntityType zero_type = MBVERTEX)
506 if constexpr (COORDINATE_SYSTEM ==
POLAR || COORDINATE_SYSTEM ==
SPHERICAL)
508 "%s coordiante not implemented",
514 vec.resize(nb_gauss_pts,
false);
522 const size_t nb_base_functions = data.
getN().size2();
525 const size_t size = nb_dofs / Tensor_Dim;
527 if (nb_dofs % Tensor_Dim) {
529 "Number of dofs should multiple of dimensions");
534 if constexpr (COORDINATE_SYSTEM ==
CARTESIAN) {
536 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
539 for (; bb != size; ++bb) {
540 values_at_gauss_pts += field_data(
I) * diff_base_function(
I);
542 ++diff_base_function;
546 for (; bb < nb_base_functions; ++bb)
547 ++diff_base_function;
548 ++values_at_gauss_pts;
558 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
561 for (; bb != size; ++bb) {
562 values_at_gauss_pts += field_data(
I) * diff_base_function(
I);
563 values_at_gauss_pts +=
564 base_function * (field_data(0) / t_coords(0));
567 ++diff_base_function;
571 for (; bb < nb_base_functions; ++bb) {
573 ++diff_base_function;
575 ++values_at_gauss_pts;
596template <
int Tensor_Dim, PetscData::DataContext CTX>
601 const std::string
field_name, boost::shared_ptr<MatrixDouble> data_ptr,
602 const EntityType zero_at_type = MBVERTEX,
bool throw_error =
true)
615 const size_t nb_dofs = local_indices.size();
620 mat.resize(Tensor_Dim, nb_gauss_pts,
false);
634 auto get_array = [&](
const auto ctx,
auto vec) {
640 <<
"In this case field degrees of freedom are read from vector. "
641 "That usually happens when time solver is used, and access to "
642 "first or second rates is needed. You probably not set ts_u, "
643 "ts_u_t, or ts_u_tt and associated data structure, i.e. "
644 "data_ctx to CTX_SET_X, CTX_SET_DX, CTX_SET_X_T, or "
645 "CTX_SET_X_TT respectively";
649 CHKERR VecGetArrayRead(vec, &array);
653 auto restore_array = [&](
auto vec) {
654 return VecRestoreArrayRead(vec, &array);
672 "That case is not implemented");
676 for (
int i = 0;
i != local_indices.size(); ++
i)
677 if (local_indices[
i] != -1)
697 "That case is not implemented");
700 const size_t nb_base_functions = data.
getN().size2();
702 auto values_at_gauss_pts = getFTensor1FromMat<Tensor_Dim>(mat);
705 const size_t size = nb_dofs / Tensor_Dim;
706 if (nb_dofs % Tensor_Dim) {
709 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
710 auto field_data = getFTensor1FromArray<Tensor_Dim, Tensor_Dim>(
dotVector);
712 for (; bb != size; ++bb) {
713 values_at_gauss_pts(
I) += field_data(
I) * base_function;
719 for (; bb < nb_base_functions; ++bb)
721 ++values_at_gauss_pts;
738template <
int Tensor_Dim>
748template <
int Tensor_Dim>
758template <
int Tensor_Dim>
772template <
int Tensor_Dim0,
int Tensor_Dim1,
class T,
class L,
class A>
778 boost::shared_ptr<ublas::matrix<T, L, A>> data_ptr,
779 const EntityType zero_type = MBVERTEX)
791 boost::shared_ptr<ublas::matrix<T, L, A>>
dataPtr;
795template <
int Tensor_Dim0,
int Tensor_Dim1,
class T,
class L,
class A>
801 "Not implemented for T = %s, dim0 = %d and dim1 = %d",
803 Tensor_Dim0, Tensor_Dim1);
807template <
int Tensor_Dim0,
int Tensor_Dim1>
815 ublas::matrix<double, ublas::row_major, DoubleAllocator>>
817 const EntityType zero_type = MBVERTEX)
828 ublas::matrix<double, ublas::row_major, DoubleAllocator>>
842 boost::shared_ptr<ublas::matrix<double, ublas::row_major, DoubleAllocator>>
849template <
int Tensor_Dim0,
int Tensor_Dim1>
851 Tensor_Dim0, Tensor_Dim1,
double, ublas::row_major,
857 const size_t nb_gauss_pts = data.
getN().size1();
858 if (type == zeroType) {
859 mat.resize(Tensor_Dim0 * Tensor_Dim1, nb_gauss_pts,
false);
865 if (dataVec.use_count()) {
866 dotVector.resize(nb_dofs,
false);
868 CHKERR VecGetArrayRead(dataVec, &array);
870 for (
int i = 0;
i != local_indices.size(); ++
i)
871 if (local_indices[
i] != -1)
872 dotVector[
i] = array[local_indices[
i]];
875 CHKERR VecRestoreArrayRead(dataVec, &array);
880 const size_t nb_base_functions = data.
getN().size2();
882 auto values_at_gauss_pts =
883 getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(mat);
886 const size_t size = nb_dofs / (Tensor_Dim0 * Tensor_Dim1);
887 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
890 for (; bb != size; ++bb) {
891 values_at_gauss_pts(
i,
j) += field_data(
i,
j) * base_function;
895 for (; bb < nb_base_functions; ++bb)
897 ++values_at_gauss_pts;
900 if (dataVec.use_count()) {
912template <
int Tensor_Dim0,
int Tensor_Dim1>
915 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
918 Tensor_Dim0, Tensor_Dim1,
double, ublas::row_major,
927template <
int Tensor_Dim0,
int Tensor_Dim1>
932 boost::shared_ptr<MatrixDouble> data_ptr,
933 const EntityType zero_at_type = MBVERTEX)
948 mat.resize(Tensor_Dim0 * Tensor_Dim1, nb_gauss_pts,
false);
952 const size_t nb_dofs = local_indices.size();
957 for (
size_t i = 0;
i != local_indices.size(); ++
i)
958 if (local_indices[
i] != -1)
964 const size_t nb_base_functions = data.
getN().size2();
967 auto values_at_gauss_pts =
968 getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(mat);
971 const size_t size = nb_dofs / (Tensor_Dim0 * Tensor_Dim1);
972 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
974 getFTensor2FromPtr<Tensor_Dim0, Tensor_Dim1>(&*
dotVector.begin());
976 for (; bb != size; ++bb) {
977 values_at_gauss_pts(
i,
j) += field_data(
i,
j) * base_function;
981 for (; bb < nb_base_functions; ++bb)
983 ++values_at_gauss_pts;
1002template <
int Tensor_Dim>
1007 const std::string
field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1008 const EntityType zero_type = MBEDGE,
const int zero_side = 0)
1017 const std::string
field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1019 const int zero_side = 0)
1034 mat.resize((Tensor_Dim * (Tensor_Dim + 1)) / 2, nb_gauss_pts,
false);
1043 const double *array;
1046 for (
int i = 0;
i != local_indices.size(); ++
i)
1047 if (local_indices[
i] != -1)
1055 const int nb_base_functions = data.
getN().size2();
1057 auto values_at_gauss_pts = getFTensor2SymmetricFromMat<Tensor_Dim>(mat);
1060 const int size = nb_dofs / ((Tensor_Dim * (Tensor_Dim + 1)) / 2);
1061 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
1064 for (; bb != size; ++bb) {
1065 values_at_gauss_pts(
i,
j) += field_data(
i,
j) * base_function;
1069 for (; bb < nb_base_functions; ++bb)
1071 ++values_at_gauss_pts;
1096template <
int Tensor_Dim>
1101 const std::string
field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1102 const EntityType zero_type = MBEDGE,
const int zero_side = 0)
1115 constexpr auto symm_size = (Tensor_Dim * (Tensor_Dim + 1)) / 2;
1117 mat.resize(symm_size, nb_gauss_pts,
false);
1121 const int nb_dofs = local_indices.size();
1129 <<
"In this case field degrees of freedom are read from vector. "
1130 "That usually happens when time solver is used, and acces to "
1131 "first rates is needed. You probably not set "
1132 "ts_u_t and associated data structure data_ctx to CTX_SET_X_T "
1139 const double *array;
1141 for (
int i = 0;
i != local_indices.size(); ++
i)
1142 if (local_indices[
i] != -1)
1148 const int nb_base_functions = data.
getN().size2();
1151 auto values_at_gauss_pts = getFTensor2SymmetricFromMat<Tensor_Dim>(mat);
1154 const int size = nb_dofs / symm_size;
1155 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
1156 auto field_data = getFTensorDotData<Tensor_Dim>();
1158 for (; bb != size; ++bb) {
1159 values_at_gauss_pts(
i,
j) += field_data(
i,
j) * base_function;
1163 for (; bb < nb_base_functions; ++bb)
1165 ++values_at_gauss_pts;
1178 static_assert(Dim || !Dim,
"not implemented");
1187 &dotVector[0], &dotVector[1], &dotVector[2], &dotVector[3], &dotVector[4],
1196 &dotVector[0], &dotVector[1], &dotVector[2]);
1210template <
int Tensor_Dim,
class T,
class L,
class A>
1215 Tensor_Dim, T,
L,
A>::OpCalculateVectorFieldValues_General;
1222template <
int Tensor_Dim>
1226 Tensor_Dim, double, ublas::row_major, DoubleAllocator> {
1229 Tensor_Dim,
double, ublas::row_major,
1248template <
int Tensor_Dim>
1250 Tensor_Dim,
double, ublas::row_major,
1255 const size_t nb_gauss_pts = this->getGaussPts().size2();
1256 auto &mat = *this->dataPtr;
1257 if (type == this->zeroType) {
1258 mat.resize(Tensor_Dim, nb_gauss_pts,
false);
1266 const int nb_base_functions = data.
getN().size2();
1268 auto gradients_at_gauss_pts = getFTensor1FromMat<Tensor_Dim>(mat);
1271 if (nb_dofs > nb_base_functions)
1273 "Number of base functions inconsistent with number of DOFs "
1275 nb_dofs, nb_base_functions);
1277 if (data.
getDiffN().size2() != nb_base_functions * Tensor_Dim)
1280 "Number of base functions inconsistent with number of derivatives "
1282 data.
getDiffN().size2(), nb_base_functions);
1284 if (data.
getDiffN().size1() != nb_gauss_pts)
1287 "Number of base functions inconsistent with number of integration "
1289 data.
getDiffN().size2(), nb_gauss_pts);
1294 for (
int gg = 0; gg < nb_gauss_pts; ++gg) {
1297 for (; bb != nb_dofs; ++bb) {
1298 gradients_at_gauss_pts(
I) += field_data * diff_base_function(
I);
1300 ++diff_base_function;
1303 for (; bb < nb_base_functions; ++bb)
1304 ++diff_base_function;
1305 ++gradients_at_gauss_pts;
1318template <
int Tensor_Dim>
1321 Tensor_Dim, double, ublas::row_major, DoubleAllocator> {
1323 Tensor_Dim,
double, ublas::row_major,
1331template <
int Tensor_Dim>
1334 Tensor_Dim, double, ublas::row_major, DoubleAllocator> {
1337 Tensor_Dim,
double, ublas::row_major,
1351template <
int Tensor_Dim>
1356 const size_t nb_gauss_pts = this->getGaussPts().size2();
1357 auto &mat = *this->dataPtr;
1358 if (type == this->zeroType) {
1359 mat.resize(Tensor_Dim * Tensor_Dim, nb_gauss_pts,
false);
1367 const int nb_base_functions = data.
getN().size2();
1369 auto &hessian_base = data.
getN(BaseDerivatives::SecondDerivative);
1371 if (hessian_base.size1() != nb_gauss_pts) {
1373 "Wrong number of integration pts (%ld != %d)",
1374 hessian_base.size1(), nb_gauss_pts);
1376 if (hessian_base.size2() != nb_base_functions * Tensor_Dim * Tensor_Dim) {
1378 "Wrong number of base functions (%ld != %d)",
1379 hessian_base.size2() / (Tensor_Dim * Tensor_Dim),
1382 if (hessian_base.size2() < nb_dofs * Tensor_Dim * Tensor_Dim) {
1384 "Wrong number of base functions (%ld < %d)",
1385 hessian_base.size2(), nb_dofs * Tensor_Dim * Tensor_Dim);
1389 auto t_diff2_base_function = getFTensor2FromPtr<Tensor_Dim, Tensor_Dim>(
1390 &*hessian_base.data().begin());
1392 auto hessian_at_gauss_pts =
1393 getFTensor2FromMat<Tensor_Dim, Tensor_Dim>(mat);
1397 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
1400 for (; bb != nb_dofs; ++bb) {
1401 hessian_at_gauss_pts(
I,
J) +=
1402 field_data * t_diff2_base_function(
I,
J);
1404 ++t_diff2_base_function;
1407 for (; bb < nb_base_functions; ++bb) {
1408 ++t_diff2_base_function;
1411 ++hessian_at_gauss_pts;
1430template <
int Tensor_Dim0,
int Tensor_Dim1,
class T,
class L,
class A>
1436 Tensor_Dim0, Tensor_Dim1, T,
L,
A>::OpCalculateTensor2FieldValues_General;
1439template <
int Tensor_Dim0,
int Tensor_Dim1>
1443 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
1446 Tensor_Dim0, Tensor_Dim1,
double, ublas::row_major,
1465template <
int Tensor_Dim0,
int Tensor_Dim1>
1467 Tensor_Dim0, Tensor_Dim1,
double, ublas::row_major,
1473 "Data pointer not allocated");
1475 const size_t nb_gauss_pts = this->getGaussPts().size2();
1476 auto &mat = *this->dataPtr;
1477 if (type == this->zeroType) {
1478 mat.resize(Tensor_Dim0 * Tensor_Dim1, nb_gauss_pts,
false);
1487 if (this->dataVec.use_count()) {
1488 this->dotVector.resize(nb_dofs,
false);
1489 const double *array;
1490 CHKERR VecGetArrayRead(this->dataVec, &array);
1492 for (
int i = 0;
i != local_indices.size(); ++
i)
1493 if (local_indices[
i] != -1)
1494 this->dotVector[
i] = array[local_indices[
i]];
1496 this->dotVector[
i] = 0;
1497 CHKERR VecRestoreArrayRead(this->dataVec, &array);
1501 const int nb_base_functions = data.
getN().size2();
1503 auto gradients_at_gauss_pts =
1504 getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(mat);
1507 int size = nb_dofs / Tensor_Dim0;
1508 if (nb_dofs % Tensor_Dim0) {
1510 "Data inconsistency");
1512 for (
int gg = 0; gg < nb_gauss_pts; ++gg) {
1516 if (field_data.l2() != field_data.l2()) {
1517 MOFEM_LOG(
"SELF", Sev::error) <<
"field data " << field_data;
1519 "Wrong number in coefficients");
1524 for (; bb < size; ++bb) {
1527 if (diff_base_function.l2() != diff_base_function.l2()) {
1529 <<
"diff_base_function: " << diff_base_function;
1531 "Wrong number number in base functions");
1536 gradients_at_gauss_pts(
I,
J) += field_data(
I) * diff_base_function(
J);
1538 ++diff_base_function;
1542 for (; bb != nb_base_functions; ++bb)
1543 ++diff_base_function;
1544 ++gradients_at_gauss_pts;
1547 if (this->dataVec.use_count()) {
1560template <
int Tensor_Dim0,
int Tensor_Dim1>
1563 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
1566 Tensor_Dim0, Tensor_Dim1,
double, ublas::row_major,
1575template <
int Tensor_Dim0,
int Tensor_Dim1>
1580 boost::shared_ptr<MatrixDouble> data_ptr,
1581 const EntityType zero_at_type = MBVERTEX)
1594 const int nb_dofs = local_indices.size();
1595 const int nb_gauss_pts = this->
getGaussPts().size2();
1599 mat.resize(Tensor_Dim0 * Tensor_Dim1, nb_gauss_pts,
false);
1606 const double *array;
1608 for (
int i = 0;
i != local_indices.size(); ++
i)
1609 if (local_indices[
i] != -1)
1615 const int nb_base_functions = data.
getN().size2();
1617 auto gradients_at_gauss_pts =
1618 getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(mat);
1621 int size = nb_dofs / Tensor_Dim0;
1622 if (nb_dofs % Tensor_Dim0) {
1626 for (
int gg = 0; gg < nb_gauss_pts; ++gg) {
1627 auto field_data = getFTensor1FromPtr<Tensor_Dim0>(&*
dotVector.begin());
1629 for (; bb < size; ++bb) {
1630 gradients_at_gauss_pts(
I,
J) += field_data(
I) * diff_base_function(
J);
1632 ++diff_base_function;
1636 for (; bb != nb_base_functions; ++bb)
1637 ++diff_base_function;
1638 ++gradients_at_gauss_pts;
1654template <
int Tensor_Dim0,
int Tensor_Dim1,
class T,
class L,
class A>
1660 const std::string
field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1661 const EntityType zero_type = MBVERTEX)
1667template <
int Tensor_Dim0,
int Tensor_Dim1>
1671 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
1674 const std::string
field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1675 const EntityType zero_type = MBVERTEX)
1697template <
int Tensor_Dim0,
int Tensor_Dim1>
1699 Tensor_Dim0, Tensor_Dim1,
double, ublas::row_major,
1705 "Data pointer not allocated");
1707 const size_t nb_gauss_pts = this->getGaussPts().size2();
1708 constexpr size_t msize = (Tensor_Dim0 * (Tensor_Dim0 + 1)) / 2;
1709 auto &mat = *this->dataPtr;
1710 if (type == this->zeroType) {
1711 mat.resize(msize * Tensor_Dim1, nb_gauss_pts,
false);
1720 const int nb_base_functions = data.
getN().size2();
1722 auto gradients_at_gauss_pts =
1723 getFTensor3DgFromMat<Tensor_Dim0, Tensor_Dim1>(mat);
1727 int size = nb_dofs / msize;
1728 if (nb_dofs % msize) {
1730 "Data inconsistency");
1732 for (
int gg = 0; gg < nb_gauss_pts; ++gg) {
1735 for (; bb < size; ++bb) {
1736 gradients_at_gauss_pts(
I,
J,
K) +=
1737 field_data(
I,
J) * diff_base_function(
K);
1739 ++diff_base_function;
1743 for (; bb != nb_base_functions; ++bb)
1744 ++diff_base_function;
1745 ++gradients_at_gauss_pts;
1757template <
int Tensor_Dim0,
int Tensor_Dim1>
1760 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
1763 const std::string
field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1764 const EntityType zero_type = MBVERTEX)
1766 Tensor_Dim0, Tensor_Dim1,
double, ublas::row_major,
1770template <
int Tensor_Dim0,
int Tensor_Dim1>
1773 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
1776 Tensor_Dim0, Tensor_Dim1,
double, ublas::row_major,
1795template <
int Tensor_Dim0,
int Tensor_Dim1>
1801 "Data pointer not allocated");
1803 const size_t nb_gauss_pts = this->getGaussPts().size2();
1804 auto &mat = *this->dataPtr;
1805 if (type == this->zeroType) {
1806 mat.resize(Tensor_Dim0 * Tensor_Dim1 * Tensor_Dim1, nb_gauss_pts,
false);
1815 if (this->dataVec.use_count()) {
1816 this->dotVector.resize(nb_dofs,
false);
1817 const double *array;
1818 CHKERR VecGetArrayRead(this->dataVec, &array);
1820 for (
int i = 0;
i != local_indices.size(); ++
i)
1821 if (local_indices[
i] != -1)
1822 this->dotVector[
i] = array[local_indices[
i]];
1824 this->dotVector[
i] = 0;
1825 CHKERR VecRestoreArrayRead(this->dataVec, &array);
1829 const int nb_base_functions = data.
getN().size2();
1831 auto &hessian_base = data.
getN(BaseDerivatives::SecondDerivative);
1833 if (hessian_base.size1() != nb_gauss_pts) {
1835 "Wrong number of integration pts (%d != %d)",
1836 hessian_base.size1(), nb_gauss_pts);
1838 if (hessian_base.size2() !=
1839 nb_base_functions * Tensor_Dim1 * Tensor_Dim1) {
1841 "Wrong number of base functions (%d != %d)",
1842 hessian_base.size2() / (Tensor_Dim1 * Tensor_Dim1),
1845 if (hessian_base.size2() < nb_dofs * Tensor_Dim1 * Tensor_Dim1) {
1847 "Wrong number of base functions (%d < %d)",
1848 hessian_base.size2(), nb_dofs * Tensor_Dim1 * Tensor_Dim1);
1852 auto t_diff2_base_function = getFTensor2FromPtr<Tensor_Dim1, Tensor_Dim1>(
1853 &*hessian_base.data().begin());
1855 auto t_hessian_at_gauss_pts =
1856 getFTensor3FromMat<Tensor_Dim0, Tensor_Dim1, Tensor_Dim1>(mat);
1862 int size = nb_dofs / Tensor_Dim0;
1864 if (nb_dofs % Tensor_Dim0) {
1866 "Data inconsistency");
1870 for (
int gg = 0; gg < nb_gauss_pts; ++gg) {
1873 for (; bb < size; ++bb) {
1874 t_hessian_at_gauss_pts(
I,
J,
K) +=
1875 field_data(
I) * t_diff2_base_function(
J,
K);
1877 ++t_diff2_base_function;
1881 for (; bb != nb_base_functions; ++bb)
1882 ++t_diff2_base_function;
1883 ++t_hessian_at_gauss_pts;
1886 if (this->dataVec.use_count()) {
1908template <
int DIM_01,
int DIM_23,
int S = 0>
1920 boost::shared_ptr<MatrixDouble> in_mat,
1921 boost::shared_ptr<MatrixDouble> out_mat,
1922 boost::shared_ptr<MatrixDouble> d_mat)
1935 boost::shared_ptr<MatrixDouble> out_mat,
1936 boost::shared_ptr<MatrixDouble> d_mat)
1949 const size_t nb_gauss_pts =
getGaussPts().size2();
1950 auto t_D = getFTensor4DdgFromMat<DIM_01, DIM_23, S>(*(
dMat));
1951 auto t_in = getFTensor2SymmetricFromMat<DIM_01>(*(
inMat));
1952 outMat->resize((DIM_23 * (DIM_23 + 1)) / 2, nb_gauss_pts,
false);
1953 auto t_out = getFTensor2SymmetricFromMat<DIM_23>(*(
outMat));
1954 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1955 t_out(
i,
j) = t_D(
i,
j,
k,
l) * t_in(
k,
l);
1970 boost::shared_ptr<MatrixDouble>
dMat;
1983 boost::shared_ptr<MatrixDouble> in_mat,
1984 boost::shared_ptr<MatrixDouble> out_mat)
1995 boost::shared_ptr<MatrixDouble> out_mat)
2006 const size_t nb_gauss_pts =
getGaussPts().size2();
2007 auto t_in = getFTensor2FromMat<DIM, DIM>(*(
inMat));
2008 outMat->resize((DIM * (DIM + 1)) / 2, nb_gauss_pts,
false);
2009 auto t_out = getFTensor2SymmetricFromMat<DIM>(*(
outMat));
2010 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
2011 t_out(
i,
j) = (t_in(
i,
j) || t_in(
j,
i)) / 2;
2034 boost::shared_ptr<MatrixDouble> in_mat,
2035 boost::shared_ptr<MatrixDouble> out_mat)
2047 boost::shared_ptr<MatrixDouble> in_mat,
2048 boost::shared_ptr<MatrixDouble> out_mat)
2062 noalias(*
outMat) = (*scalePtr) * (*inMat);
2081template <
int Base_Dim,
int Field_Dim,
class T,
class L,
class A>
2087template <
int Field_Dim>
2093 boost::shared_ptr<MatrixDouble> data_ptr,
2095 const EntityType zero_type = MBEDGE,
2096 const int zero_side = 0)
2099 dataPtr(data_ptr), dataVec(data_vec), zeroType(zero_type),
2100 zeroSide(zero_side) {
2106 boost::shared_ptr<MatrixDouble> data_ptr,
2107 const EntityType zero_type = MBEDGE,
2108 const int zero_side = 0)
2130template <
int Field_Dim>
2132 3, Field_Dim,
double, ublas::row_major,
2136 const size_t nb_integration_points = this->getGaussPts().size2();
2137 if (type == zeroType && side == zeroSide) {
2138 dataPtr->resize(Field_Dim, nb_integration_points,
false);
2145 if (dataVec.use_count()) {
2146 dotVector.resize(nb_dofs,
false);
2147 const double *array;
2148 CHKERR VecGetArrayRead(dataVec, &array);
2150 for (
int i = 0;
i != local_indices.size(); ++
i)
2151 if (local_indices[
i] != -1)
2152 dotVector[
i] = array[local_indices[
i]];
2155 CHKERR VecRestoreArrayRead(dataVec, &array);
2159 const size_t nb_base_functions = data.
getN().size2() / 3;
2162 auto t_data = getFTensor1FromMat<Field_Dim>(*dataPtr);
2163 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
2166 for (; bb != nb_dofs; ++bb) {
2167 t_data(
i) += t_n_hdiv(
i) * t_dof;
2171 for (; bb != nb_base_functions; ++bb)
2176 if (dataVec.use_count()) {
2186template <
int Base_Dim,
int Field_Dim = Base_Dim>
2189 Base_Dim, Field_Dim, double, ublas::row_major, DoubleAllocator> {
2191 Base_Dim, Field_Dim,
double, ublas::row_major,
2198template <
int Base_Dim,
int Field_Dim = Base_Dim>
2201template <
int Field_Dim>
2206 boost::shared_ptr<MatrixDouble> data_ptr,
2207 const EntityType zero_type = MBEDGE,
2208 const int zero_side = 0)
2211 dataPtr(data_ptr), zeroType(zero_type), zeroSide(zero_side) {
2232template <
int Field_Dim>
2237 const size_t nb_integration_points = this->getGaussPts().size2();
2238 if (type == zeroType && side == zeroSide) {
2239 dataPtr->resize(Field_Dim, nb_integration_points,
false);
2244 const size_t nb_dofs = local_indices.size();
2247 std::array<double, MAX_DOFS_ON_ENTITY> dot_dofs_vector;
2248 const double *array;
2249 CHKERR VecGetArrayRead(getFEMethod()->ts_u_t, &array);
2250 for (
size_t i = 0;
i != nb_dofs; ++
i)
2251 if (local_indices[
i] != -1)
2252 dot_dofs_vector[
i] = array[local_indices[
i]];
2254 dot_dofs_vector[
i] = 0;
2255 CHKERR VecRestoreArrayRead(getFEMethod()->ts_u_t, &array);
2257 const size_t nb_base_functions = data.
getN().size2() / 3;
2260 auto t_data = getFTensor1FromMat<Field_Dim>(*dataPtr);
2261 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
2263 for (; bb != nb_dofs; ++bb) {
2264 t_data(
i) += t_n_hdiv(
i) * dot_dofs_vector[bb];
2267 for (; bb != nb_base_functions; ++bb)
2283template <
int BASE_DIM,
int SPACE_DIM>
2288 boost::shared_ptr<VectorDouble> data_ptr,
2289 const EntityType zero_type = MBEDGE,
2290 const int zero_side = 0)
2301 const size_t nb_integration_points =
getGaussPts().size2();
2303 dataPtr->resize(nb_integration_points,
false);
2309 const size_t nb_base_functions = data.
getN().size2() /
BASE_DIM;
2314 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
2317 for (; bb != nb_dofs; ++bb) {
2318 t_data += t_dof * t_n_diff_hdiv(
j,
j);
2322 for (; bb != nb_base_functions; ++bb)
2342template <
int BASE_DIM,
int SPACE_DIM>
2347 boost::shared_ptr<MatrixDouble> data_ptr,
2348 const EntityType zero_type = MBEDGE,
2349 const int zero_side = 0)
2360 const size_t nb_integration_points =
getGaussPts().size2();
2368 const size_t nb_base_functions = data.
getN().size2() /
BASE_DIM;
2372 auto t_data = getFTensor2FromMat<BASE_DIM, SPACE_DIM>(*
dataPtr);
2373 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
2376 for (; bb != nb_dofs; ++bb) {
2377 t_data(
i,
j) += t_dof * t_base_diff(
i,
j);
2381 for (; bb != nb_base_functions; ++bb)
2401template <
int BASE_DIM,
int SPACE_DIM>
2406 boost::shared_ptr<MatrixDouble> data_ptr,
2407 const EntityType zero_type = MBEDGE,
2408 const int zero_side = 0)
2419 const size_t nb_integration_points =
getGaussPts().size2();
2429 const int nb_base_functions = data.
getN().size2() /
BASE_DIM;
2432 auto &hessian_base = data.
getN(BaseDerivatives::SecondDerivative);
2433 if (hessian_base.size1() != nb_integration_points) {
2435 "Wrong number of integration pts (%d != %d)",
2436 hessian_base.size1(), nb_integration_points);
2438 if (hessian_base.size2() !=
2441 "Wrong number of base functions (%d != %d)",
2447 "Wrong number of base functions (%d < %d)", hessian_base.size2(),
2458 auto t_data = getFTensor3FromMat<BASE_DIM, SPACE_DIM, SPACE_DIM>(*
dataPtr);
2459 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
2462 for (; bb != nb_dofs; ++bb) {
2463 t_data(
i,
j,
k) += t_dof * t_base_diff2(
i,
j,
k);
2468 for (; bb != nb_base_functions; ++bb)
2487template <
int Tensor_Dim1,
int Tensor_Dim2>
2492 boost::shared_ptr<VectorDouble> data_ptr,
2493 const EntityType zero_type = MBEDGE,
2494 const int zero_side = 0)
2505 const size_t nb_integration_points =
getGaussPts().size2();
2507 dataPtr->resize(nb_integration_points,
false);
2512 const int nb_dofs = local_indices.size();
2515 std::array<double, MAX_DOFS_ON_ENTITY> dot_dofs_vector;
2516 const double *array;
2518 for (
size_t i = 0;
i != local_indices.size(); ++
i)
2519 if (local_indices[
i] != -1)
2520 dot_dofs_vector[
i] = array[local_indices[
i]];
2522 dot_dofs_vector[
i] = 0;
2525 const size_t nb_base_functions = data.
getN().size2() / Tensor_Dim1;
2529 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
2531 for (; bb != nb_dofs; ++bb) {
2533 for (
auto ii = 0; ii != Tensor_Dim2; ++ii)
2534 div += t_n_diff_hdiv(ii, ii);
2535 t_data += dot_dofs_vector[bb] * div;
2538 for (; bb != nb_base_functions; ++bb)
2574 boost::shared_ptr<MatrixDouble> data_ptr,
2575 const EntityType zero_type = MBEDGE,
2576 const int zero_side = 0);
2598 boost::shared_ptr<MatrixDouble> data_ptr,
2599 const EntityType zero_type = MBVERTEX,
2600 const int zero_side = 0);
2623 boost::shared_ptr<MatrixDouble> data_ptr,
2624 const EntityType zero_type = MBVERTEX,
2625 const int zero_side = 0);
2643template <
int Tensor_Dim0,
int Tensor_Dim1>
2648 boost::shared_ptr<MatrixDouble> data_ptr,
2649 boost::shared_ptr<double> scale_ptr,
2651 const EntityType zero_type = MBEDGE,
2652 const int zero_side = 0)
2662 boost::shared_ptr<MatrixDouble> data_ptr,
2663 const EntityType zero_type = MBEDGE,
2664 const int zero_side = 0)
2672 const size_t nb_integration_points =
getGaussPts().size2();
2674 dataPtr->resize(Tensor_Dim0 * Tensor_Dim1, nb_integration_points,
false);
2682 const double *array;
2685 for (
int i = 0;
i != local_indices.size(); ++
i)
2686 if (local_indices[
i] != -1)
2695 const size_t nb_base_functions = data.
getN().size2() / 3;
2699 auto t_data = getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(*
dataPtr);
2700 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
2703 for (; bb != nb_dofs / Tensor_Dim0; ++bb) {
2704 t_data(
i,
j) += (
scale * t_dof(
i)) * t_n_hvec(
j);
2708 for (; bb < nb_base_functions; ++bb)
2734template <
int Tensor_Dim0,
int Tensor_Dim1>
2741 const std::string
field_name, boost::shared_ptr<MatrixDouble> data_ptr,
2743 boost::shared_ptr<Range> broken_range_ptr =
nullptr,
2744 boost::shared_ptr<double> scale_ptr =
nullptr,
2745 const EntityType zero_type = MBEDGE,
const int zero_side = 0)
2773template <
int Tensor_Dim0,
int Tensor_Dim1>
2778 const size_t nb_integration_points = OP::getGaussPts().size2();
2779 if (type == zeroType) {
2780 dataPtr->resize(Tensor_Dim0 * Tensor_Dim1, nb_integration_points,
false);
2787 if (dataVec.use_count()) {
2788 dotVector.resize(nb_dofs,
false);
2789 const double *array;
2790 CHKERR VecGetArrayRead(dataVec, &array);
2792 for (
int i = 0;
i != local_indices.size(); ++
i)
2793 if (local_indices[
i] != -1)
2794 dotVector[
i] = array[local_indices[
i]];
2797 CHKERR VecRestoreArrayRead(dataVec, &array);
2808 auto get_get_side_face_dofs = [&]() {
2809 auto fe_type = OP::getFEType();
2813 std::vector<int> side_face_dofs;
2814 side_face_dofs.reserve(data.
getIndices().size() / Tensor_Dim0);
2818 auto it = side_dof_map.get<1>().begin();
2819 it != side_dof_map.get<1>().end(); ++it
2822 if ((Tensor_Dim0 * it->dof) >= data.
getIndices().size()) {
2825 if (it->type == brokenType) {
2826 if (brokenRangePtr) {
2827 auto ent = OP::getSideEntity(it->side, brokenType);
2828 if (brokenRangePtr->find(ent) != brokenRangePtr->end()) {
2829 side_face_dofs.push_back(it->dof);
2832 side_face_dofs.push_back(it->dof);
2837 return side_face_dofs;
2840 auto side_face_dofs = get_get_side_face_dofs();
2844 auto t_data = getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(*dataPtr);
2845 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
2846 for (
auto b : side_face_dofs) {
2848 auto t_dof = getFTensor1FromPtr<Tensor_Dim0>(data.
getFieldData().data() +
2850 t_data(
i,
j) += t_dof(
i) * t_row_base(
j);
2854 *dataPtr *= (scalePtr) ? *scalePtr : 1.0;
2856 if (dataVec.use_count()) {
2870template <
int Tensor_Dim0,
int Tensor_Dim1>
2875 const std::string
field_name, boost::shared_ptr<MatrixDouble> data_ptr,
2876 boost::shared_ptr<double> scale_ptr,
2878 const EntityType zero_type = MBEDGE,
const int zero_side = 0)
2888 boost::shared_ptr<MatrixDouble> data_ptr,
2889 const EntityType zero_type = MBEDGE,
2890 const int zero_side = 0)
2898 const size_t nb_integration_points =
getGaussPts().size2();
2900 dataPtr->resize(Tensor_Dim0 * Tensor_Dim1, nb_integration_points,
false);
2909 const double *array;
2912 for (
int i = 0;
i != local_indices.size(); ++
i)
2913 if (local_indices[
i] != -1)
2922 const size_t nb_base_functions =
2923 data.
getN().size2() / (Tensor_Dim0 * Tensor_Dim1);
2926 auto t_n_hten = data.
getFTensor2N<Tensor_Dim0, Tensor_Dim1>();
2927 auto t_data = getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(*
dataPtr);
2928 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
2931 for (; bb != nb_dofs; ++bb) {
2932 t_data(
i,
j) += (
scale * t_dof) * t_n_hten(
i,
j);
2936 for (; bb < nb_base_functions; ++bb)
2964template <
int Tensor_Dim0,
int Tensor_Dim1,
2970 boost::shared_ptr<MatrixDouble> data_ptr,
2972 const EntityType zero_type = MBEDGE,
2973 const int zero_side = 0)
2983 boost::shared_ptr<MatrixDouble> data_ptr,
2984 const EntityType zero_type = MBEDGE,
2985 const int zero_side = 0)
2992 const size_t nb_integration_points =
getGaussPts().size2();
2994 dataPtr->resize(Tensor_Dim0, nb_integration_points,
false);
3002 const double *array;
3005 for (
int i = 0;
i != local_indices.size(); ++
i)
3006 if (local_indices[
i] != -1)
3014 const size_t nb_base_functions = data.
getN().size2() / 3;
3018 auto t_data = getFTensor1FromMat<Tensor_Dim0>(*
dataPtr);
3021 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
3024 for (; bb != nb_dofs / Tensor_Dim0; ++bb) {
3025 double div = t_n_diff_hvec(
j,
j);
3026 t_data(
i) += t_dof(
i) * div;
3028 t_data(
i) += t_base(0) * (t_dof(
i) / t_coords(0));
3034 for (; bb < nb_base_functions; ++bb) {
3067template <
int Tensor_Dim0,
int Tensor_Dim1,
3075 const std::string
field_name, boost::shared_ptr<MatrixDouble> data_ptr,
3077 boost::shared_ptr<Range> broken_range_ptr =
nullptr,
3078 boost::shared_ptr<double> scale_ptr =
nullptr,
3079 const EntityType zero_type = MBEDGE)
3091 const size_t nb_integration_points =
getGaussPts().size2();
3092 if (type ==
zeroType && side == 0) {
3093 dataPtr->resize(Tensor_Dim0, nb_integration_points,
false);
3102 const double *array;
3105 for (
int i = 0;
i != local_indices.size(); ++
i)
3106 if (local_indices[
i] != -1)
3121 auto get_get_side_face_dofs = [&]() {
3126 std::vector<int> side_face_dofs;
3127 side_face_dofs.reserve(data.
getIndices().size() / Tensor_Dim0);
3131 auto it = side_dof_map.get<1>().begin();
3132 it != side_dof_map.get<1>().end(); ++it
3135 if ((Tensor_Dim0 * it->dof) >= data.
getIndices().size()) {
3142 side_face_dofs.push_back(it->dof);
3145 side_face_dofs.push_back(it->dof);
3150 return side_face_dofs;
3153 auto side_face_dofs = get_get_side_face_dofs();
3157 auto t_data = getFTensor1FromMat<Tensor_Dim0>(*
dataPtr);
3159 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
3160 for (
auto b : side_face_dofs) {
3161 auto t_dof = getFTensor1FromPtr<Tensor_Dim0>(
3165 double div = t_diff_base(
j,
j);
3166 t_data(
i) += t_dof(
i) * div;
3168 t_data(
i) += t_base(0) * (t_dof(
i) / t_coords(0));
3199template <
int Tensor_Dim,
typename OpBase>
3203 boost::shared_ptr<MatrixDouble> data_ptr,
3204 boost::shared_ptr<double> scale_ptr,
3205 const EntityType zero_type = MBEDGE,
3206 const int zero_side = 0)
3214 boost::shared_ptr<MatrixDouble> data_ptr,
3215 const EntityType zero_type = MBEDGE,
3216 const int zero_side = 0)
3223 const size_t nb_integration_points = OpBase::getGaussPts().size2();
3224 if (type ==
zeroType && side == 0) {
3225 dataPtr->resize(Tensor_Dim, nb_integration_points,
false);
3231 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
3232 const size_t nb_base_functions = data.
getN().size2() / 3;
3234 auto t_data = getFTensor1FromMat<Tensor_Dim>(*
dataPtr);
3235 for (
size_t gg = 0; gg != nb_integration_points; ++gg) {
3237 t_normalized_normal(
j) = t_normal(
j);
3241 for (; bb != nb_dofs / Tensor_Dim; ++bb) {
3243 (scale_val * t_dof(
i)) * (t_base(
j) * t_normalized_normal(
j));
3247 for (; bb < nb_base_functions; ++bb) {
3290 boost::shared_ptr<MatrixDouble> inv_jac_ptr);
3293 template <
int D1,
int D2,
int J1,
int J2>
3297 static_assert(D2 == J2,
"Dimension of jacobian and dimension of <out> "
3298 "directive does not match");
3300 size_t nb_functions = diff_n.size2() / D1;
3302 size_t nb_gauss_pts = diff_n.size1();
3307 "Wrong number of Gauss Pts");
3308 if (diff_n.size2() % D1)
3310 "Number of directives of base functions and D1 dimension does "
3314 diffNinvJac.resize(diff_n.size1(), D2 * nb_functions,
false);
3318 auto t_diff_n = getFTensor1FromPtr<D2>(&*
diffNinvJac.data().begin());
3319 auto t_diff_n_ref = getFTensor1FromPtr<D1>(&*diff_n.data().begin());
3320 auto t_inv_jac = getFTensor2FromMat<J1, J2>(*
invJacPtr);
3321 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg, ++t_inv_jac) {
3322 for (
size_t dd = 0; dd != nb_functions; ++dd) {
3323 t_diff_n(
i) = t_inv_jac(
k,
i) * t_diff_n_ref(
k);
3368template <
int DERIVARIVE = 1>
3375template <
int DERIVARIVE = 1>
3382template <
int DERIVARIVE = 1>
3386 boost::shared_ptr<MatrixDouble> inv_jac_ptr)
3390template <
int DERIVARIVE = 1>
3394 boost::shared_ptr<MatrixDouble> inv_jac_ptr)
3412 invJacPtr(inv_jac_ptr) {}
3466 boost::shared_ptr<MatrixDouble> inv_jac_ptr);
3503 boost::shared_ptr<MatrixDouble> jac_ptr)
3520 2>::OpSetContravariantPiolaTransformOnFace2DImpl;
3703 boost::shared_ptr<VectorDouble> det_ptr,
3704 boost::shared_ptr<MatrixDouble> out_ptr)
3724 "Pointer for inPtr matrix not allocated");
3727 "Pointer for detPtr matrix not allocated");
3729 const auto nb_rows = inPtr->size1();
3730 const auto nb_integration_pts = inPtr->size2();
3734 detPtr->resize(nb_integration_pts,
false);
3735 auto t_in = getFTensor2FromMat<DIM, DIM>(*inPtr);
3737 for (
size_t gg = 0; gg != nb_integration_pts; ++gg) {
3746 outPtr->resize(nb_rows, nb_integration_pts,
false);
3747 auto t_in = getFTensor2FromMat<DIM, DIM>(*inPtr);
3748 auto t_out = getFTensor2FromMat<DIM, DIM>(*outPtr);
3750 for (
size_t gg = 0; gg != nb_integration_pts; ++gg) {
3773 boost::shared_ptr<VectorDouble> out_ptr)
3794 "Pointer for inPtr matrix not allocated");
3796 const auto nb_integration_pts = inPtr->size2();
3799 outPtr->resize(nb_integration_pts,
false);
3800 auto t_in = getFTensor2FromMat<DIM, DIM>(*inPtr);
3803 for (
size_t gg = 0; gg != nb_integration_pts; ++gg) {
3826 boost::shared_ptr<VectorDouble> out_ptr)
3847 "Pointer for inPtr matrix not allocated");
3849 const auto nb_integration_pts = inPtr->size2();
3852 outPtr->resize(nb_integration_pts,
false);
3853 auto t_in = getFTensor2SymmetricFromMat<DIM>(*inPtr);
3856 for (
size_t gg = 0; gg != nb_integration_pts; ++gg) {
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Tensor1< T, Tensor_Dim > normalize()
static const char *const CoordinateTypesNames[]
Coordinate system names.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FieldSpace
approximation spaces
@ L2
field with C-1 continuity
@ HCURL
field with continuous tangents
#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 THROW_MESSAGE(msg)
Throw MoFEM exception.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'J', DIM1 > J
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VecAllocator< double > DoubleAllocator
implementation of Data Operators for Forces and Sources
static MoFEMErrorCode invertTensor(FTensor::Tensor2< T1, DIM, DIM > &t, T2 &det, FTensor::Tensor2< T3, DIM, DIM > &inv_t)
static auto determinantTensor(FTensor::Tensor2< T, DIM, DIM > &t)
Calculate the determinant of a tensor of rank DIM.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
OpCalculateScalarFieldValuesFromPetscVecImpl< PetscData::CTX_SET_X_T > OpCalculateScalarFieldValuesDot
constexpr IntegrationType I
constexpr auto field_name
multi_index_container< DofsSideMapData, indexed_by< ordered_non_unique< tag< TypeSide_mi_tag >, composite_key< DofsSideMapData, member< DofsSideMapData, EntityType, &DofsSideMapData::type >, member< DofsSideMapData, int, &DofsSideMapData::side > > >, ordered_unique< tag< EntDofIdx_mi_tag >, member< DofsSideMapData, int, &DofsSideMapData::dof > > > > DofsSideMap
Map entity stype and side to element/entity dof index.
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
default operator for EDGE element
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::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
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....
const VectorDouble & getFieldData() const
get dofs values
FTensor::Tensor3< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 *Tensor_Dim2 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3Diff2N(FieldApproximationBase base)
Get second derivatives of base functions for Hvec space.
const VectorInt & getLocalIndices() const
get local indices of dofs on entity
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1FieldData()
Return FTensor of rank 1, i.e. vector from field data coefficients.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData()
Resturn scalar files as a FTensor of rank 0.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
FTensor::Tensor2_symmetric< FTensor::PackPtr< double *,(Tensor_Dim *(Tensor_Dim+1))/2 >, Tensor_Dim > getFTensor2SymmetricFieldData()
Return symmetric FTensor rank 2, i.e. matrix from field data coefficients.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2FieldData()
Return FTensor rank 2, i.e. matrix from field data coefficients.
const VectorInt & getIndices() const
Get global indices of dofs on entity.
default operator for TRI element
default operator for Flat Prism element
default operator for Flat Prism element
FlatPrism finite element.
EntityType getFEType() const
Get dimension of finite element.
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
EntityHandle getSideEntity(const int side_number, const EntityType type)
Get the side entity.
@ OPROW
operator doWork function is executed on FE rows
@ OPSPACE
operator do Work is execute on space data
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
structure to get information form mofem into EntitiesFieldData
Calculate divergence of tonsorial field using vectorial base.
boost::shared_ptr< double > scalePtr
boost::shared_ptr< Range > brokenRangePtr
SmartPetscObj< Vec > dataVec
OpCalculateBrokenHVecTensorDivergence(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, EntityType broken_type, boost::shared_ptr< Range > broken_range_ptr=nullptr, boost::shared_ptr< double > scale_ptr=nullptr, const EntityType zero_type=MBEDGE)
boost::shared_ptr< MatrixDouble > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
const EntityHandle zeroType
Get tensor field for H-div approximation.
boost::shared_ptr< MatrixDouble > dataPtr
boost::shared_ptr< Range > brokenRangePtr
boost::shared_ptr< double > scalePtr
OpCalculateBrokenHVecTensorField(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, EntityType broken_type, boost::shared_ptr< Range > broken_range_ptr=nullptr, boost::shared_ptr< double > scale_ptr=nullptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Calculate values of vector field at integration points.
const EntityHandle zeroType
SmartPetscObj< Vec > dataVec
Calculate field values (template specialization) for tensor field rank 1, i.e. vector field.
boost::shared_ptr< VectorDouble > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateDivergenceVectorFieldValues(const std::string field_name, boost::shared_ptr< VectorDouble > data_ptr, const EntityType zero_type=MBVERTEX)
const EntityHandle zeroType
Calculate tenor field using tensor base, i.e. Hdiv/Hcurl.
VectorDouble dotVector
Keeps temporary values of time derivatives.
OpCalculateHTensorTensorField(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, boost::shared_ptr< double > scale_ptr, SmartPetscObj< Vec > data_vec=SmartPetscObj< Vec >(), const EntityType zero_type=MBEDGE, const int zero_side=0)
SmartPetscObj< Vec > dataVec
boost::shared_ptr< MatrixDouble > dataPtr
const EntityHandle zeroType
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateHTensorTensorField(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
boost::shared_ptr< double > scalePtr
Calculate divergence of tonsorial field using vectorial base.
boost::shared_ptr< MatrixDouble > dataPtr
SmartPetscObj< Vec > dataVec
VectorDouble dotVector
Keeps temporary values of time derivatives.
const EntityHandle zeroType
OpCalculateHVecTensorDivergence(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBEDGE, const int zero_side=0)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateHVecTensorDivergence(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
Calculate tenor field using vectorial base, i.e. Hdiv/Hcurl.
OpCalculateHVecTensorField(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, boost::shared_ptr< double > scale_ptr, SmartPetscObj< Vec > data_vec=SmartPetscObj< Vec >(), const EntityType zero_type=MBEDGE, const int zero_side=0)
const EntityHandle zeroType
boost::shared_ptr< MatrixDouble > dataPtr
OpCalculateHVecTensorField(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
SmartPetscObj< Vec > dataVec
VectorDouble dotVector
Keeps temporary values of time derivatives.
boost::shared_ptr< double > scalePtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Calculate trace of vector (Hdiv/Hcurl) space.
const EntityHandle zeroType
boost::shared_ptr< MatrixDouble > dataPtr
boost::shared_ptr< double > scalePtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
OpCalculateHVecTensorTrace(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
FTensor::Index< 'j', Tensor_Dim > j
FTensor::Index< 'i', Tensor_Dim > i
OpCalculateHVecTensorTrace(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, boost::shared_ptr< double > scale_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
const EntityHandle zeroType
boost::shared_ptr< MatrixDouble > dataPtr
OpCalculateHVecVectorFieldDot(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
Get vector field for H-div approximation.
OpCalculateHVecVectorField_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBEDGE, const int zero_side=0)
boost::shared_ptr< MatrixDouble > dataPtr
OpCalculateHVecVectorField_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
SmartPetscObj< Vec > dataVec
const EntityHandle zeroType
Get vector field for H-div approximation.
Get vector field for H-div approximation.
Calculate gradient of vector field.
boost::shared_ptr< MatrixDouble > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
const EntityHandle zeroType
OpCalculateHVecVectorGradient(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
Calculate gradient of vector field.
boost::shared_ptr< MatrixDouble > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
const EntityHandle zeroType
OpCalculateHVecVectorHessian(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
boost::shared_ptr< MatrixDouble > dataPtr
const EntityHandle zeroType
const EntityHandle zeroType
boost::shared_ptr< MatrixDouble > dataPtr
boost::shared_ptr< MatrixDouble > dataPtr
const EntityHandle zeroType
Calculate curl of vector field.
Calculate divergence of vector field dot.
boost::shared_ptr< VectorDouble > dataPtr
OpCalculateHdivVectorDivergenceDot(const std::string field_name, boost::shared_ptr< VectorDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
const EntityHandle zeroType
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Calculate divergence of vector field.
boost::shared_ptr< VectorDouble > dataPtr
OpCalculateHdivVectorDivergence(const std::string field_name, boost::shared_ptr< VectorDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
const EntityHandle zeroType
Calculate inverse of jacobian for face element.
OpCalculateInvJacForFatPrism(MatrixDouble &inv_jac)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateInvJacForFatPrism(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
const boost::shared_ptr< MatrixDouble > invJacPtr
Calculate inverse of jacobian for face element.
OpCalculateInvJacForFlatPrism(MatrixDouble &inv_jac_f3)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Evaluate field gradient values for scalar field, i.e. gradient is tensor rank 1 (vector)
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
Evaluate field gradient values for scalar field, i.e. gradient is tensor rank 1 (vector),...
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
calculate gradient values of scalar field at integration points
Get rate of scalar field at integration points.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< VectorDouble > dataPtr
OpCalculateScalarFieldValuesFromPetscVecImpl(const std::string field_name, boost::shared_ptr< VectorDouble > data_ptr, const EntityType zero_at_type=MBVERTEX)
const EntityHandle zeroAtType
Scalar field values at integration points.
OpCalculateScalarFieldValues_General(const std::string field_name, boost::shared_ptr< ublas::vector< T, A > > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBVERTEX)
boost::shared_ptr< ublas::vector< T, A > > dataPtr
SmartPetscObj< Vec > dataVec
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
calculate values of scalar field at integration points
const EntityHandle zeroType
OpCalculateScalarFieldValues_General(const std::string field_name, boost::shared_ptr< ublas::vector< T, A > > data_ptr, const EntityType zero_type=MBVERTEX)
Get value at integration points for scalar field.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
calculate values of scalar field at integration points
Get time direvarive values at integration pts for tensor field rank 2, i.e. matrix field.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateTensor2FieldValuesDot(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_at_type=MBVERTEX)
EntityType zeroAtType
Zero values at Gauss point at this type.
VectorDouble dotVector
Keeps temporary values of time derivatives.
boost::shared_ptr< MatrixDouble > dataPtr
Data computed into this matrix.
SmartPetscObj< Vec > dataVec
OpCalculateTensor2FieldValues_General(const std::string field_name, boost::shared_ptr< ublas::matrix< double, ublas::row_major, DoubleAllocator > > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBVERTEX)
OpCalculateTensor2FieldValues_General(const std::string field_name, boost::shared_ptr< ublas::matrix< double, ublas::row_major, DoubleAllocator > > data_ptr, const EntityType zero_type=MBVERTEX)
boost::shared_ptr< ublas::matrix< double, ublas::row_major, DoubleAllocator > > dataPtr
const EntityHandle zeroType
Calculate field values for tenor field rank 2.
const EntityHandle zeroType
OpCalculateTensor2FieldValues_General(const std::string field_name, boost::shared_ptr< ublas::matrix< T, L, A > > data_ptr, const EntityType zero_type=MBVERTEX)
boost::shared_ptr< ublas::matrix< T, L, A > > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Get values at integration pts for tensor field rank 2, i.e. matrix field.
OpCalculateTensor2SymmetricFieldGradient_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBVERTEX)
Evaluate field gradient values for symmetric 2nd order tensor field, i.e. gradient is tensor rank 3.
OpCalculateTensor2SymmetricFieldGradient_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBVERTEX)
Get field gradients at integration pts for symmetric tensorial field rank 2.
OpCalculateTensor2SymmetricFieldGradient(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBVERTEX)
Calculate symmetric tensor field rates ant integratio pts.
const EntityHandle zeroType
boost::shared_ptr< MatrixDouble > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateTensor2SymmetricFieldValuesDot(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
Calculate symmetric tensor field values at integration pts.
OpCalculateTensor2SymmetricFieldValues(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBEDGE, const int zero_side=0)
OpCalculateTensor2SymmetricFieldValues(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
SmartPetscObj< Vec > dataVec
boost::shared_ptr< MatrixDouble > dataPtr
const EntityHandle zeroType
Calculates the trace of an input matrix.
boost::shared_ptr< MatrixDouble > inPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
FTensor::Index< 'i', DIM > i
boost::shared_ptr< VectorDouble > outPtr
OpCalculateTraceFromMat(boost::shared_ptr< MatrixDouble > in_ptr, boost::shared_ptr< VectorDouble > out_ptr)
Calculates the trace of an input matrix.
FTensor::Index< 'i', DIM > i
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< VectorDouble > outPtr
boost::shared_ptr< MatrixDouble > inPtr
OpCalculateTraceFromSymmMat(boost::shared_ptr< MatrixDouble > in_ptr, boost::shared_ptr< VectorDouble > out_ptr)
Get field gradients time derivative at integration pts for scalar field rank 0, i....
boost::shared_ptr< MatrixDouble > dataPtr
Data computed into this matrix.
OpCalculateVectorFieldGradientDot(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_at_type=MBVERTEX)
EntityType zeroAtType
Zero values at Gauss point at this type.
VectorDouble dotVector
Keeps temporary values of time derivatives.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Evaluate field gradient values for vector field, i.e. gradient is tensor rank 2.
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
calculate values of vector field at integration points
Approximate field values for given petsc vector.
boost::shared_ptr< MatrixDouble > dataPtr
OpCalculateVectorFieldValuesFromPetscVecImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_at_type=MBVERTEX, bool throw_error=true)
const EntityHandle zeroAtType
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< MatrixDouble > dataPtr
OpCalculateVectorFieldValues_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBVERTEX)
OpCalculateVectorFieldValues_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBVERTEX)
SmartPetscObj< Vec > dataVec
const EntityHandle zeroType
Calculate field values for tenor field rank 1, i.e. vector field.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
calculate values of vector field at integration points
OpCalculateVectorFieldValues_General(const std::string field_name, boost::shared_ptr< ublas::matrix< T, L, A > > data_ptr, const EntityType zero_type=MBVERTEX)
const EntityHandle zeroType
boost::shared_ptr< ublas::matrix< T, L, A > > dataPtr
Get values at integration pts for tensor field rank 1, i.e. vector field.
boost::shared_ptr< VectorDouble > detPtr
OpInvertMatrix(boost::shared_ptr< MatrixDouble > in_ptr, boost::shared_ptr< VectorDouble > det_ptr, boost::shared_ptr< MatrixDouble > out_ptr)
boost::shared_ptr< MatrixDouble > outPtr
boost::shared_ptr< MatrixDouble > inPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Make Hdiv space from Hcurl space in 2d.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Operator for fat prism element updating integration weights in the volume.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpMultiplyDeterminantOfJacobianAndWeightsForFatPrisms()
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< MatrixDouble > outMat
DEPRECATED OpScaleMatrix(const std::string field_name, const double scale, boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat)
boost::shared_ptr< double > scalePtr
OpScaleMatrix(boost::shared_ptr< double > scale_ptr, boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat)
boost::shared_ptr< MatrixDouble > inMat
OpSetInvJacH1ForFaceEmbeddedIn3DSpace(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
OpSetInvJacH1ForFace(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Transform local reference derivatives of shape functions to global derivatives.
const boost::shared_ptr< MatrixDouble > invJacPtr
OpSetInvJacH1ForFatPrism(MatrixDouble &inv_jac)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpSetInvJacH1ForFatPrism(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Transform local reference derivatives of shape functions to global derivatives.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpSetInvJacH1ForFlatPrism(MatrixDouble &inv_jac_f3)
boost::shared_ptr< MatrixDouble > invJacPtr
MatrixDouble diffHcurlInvJac
OpSetInvJacHcurlFaceImpl(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Transform local reference derivatives of shape function to global derivatives for face.
OpSetInvJacL2ForFaceEmbeddedIn3DSpace(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
OpSetInvJacL2ForFace(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Transform local reference derivatives of shape functions to global derivatives.
MoFEMErrorCode applyTransform(MatrixDouble &diff_n)
OpSetInvJacToScalarBasesBasic(FieldSpace space, boost::shared_ptr< MatrixDouble > inv_jac_ptr)
boost::shared_ptr< MatrixDouble > invJacPtr
FTensor::Index< 'i', DIM > i
boost::shared_ptr< MatrixDouble > outMat
OpSymmetrizeTensor(boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
FTensor::Index< 'j', DIM > j
DEPRECATED OpSymmetrizeTensor(const std::string field_name, boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat)
boost::shared_ptr< MatrixDouble > inMat
FTensor::Index< 'i', DIM_01 > i
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
FTensor::Index< 'k', DIM_23 > k
boost::shared_ptr< MatrixDouble > dMat
OpTensorTimesSymmetricTensor(boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat, boost::shared_ptr< MatrixDouble > d_mat)
FTensor::Index< 'j', DIM_01 > j
boost::shared_ptr< MatrixDouble > inMat
FTensor::Index< 'l', DIM_23 > l
DEPRECATED OpTensorTimesSymmetricTensor(const std::string field_name, boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat, boost::shared_ptr< MatrixDouble > d_mat)
boost::shared_ptr< MatrixDouble > outMat
static constexpr Switches CtxSetX
static constexpr Switches CtxSetX_TT
std::bitset< 8 > Switches
static constexpr Switches CtxSetX_T
intrusive_ptr for managing petsc objects