v0.8.19
Public Member Functions | Public Attributes | List of all members
MoFEM::OpGetDataAndGradient< RANK, DIM > Struct Template Reference

Get field values and gradients at Gauss points. More...

#include <src/finite_elements/DataOperators.hpp>

Inheritance diagram for MoFEM::OpGetDataAndGradient< RANK, DIM >:
[legend]
Collaboration diagram for MoFEM::OpGetDataAndGradient< RANK, DIM >:
[legend]

Public Member Functions

 OpGetDataAndGradient (MatrixDouble &data_at_gauss_pt, MatrixDouble &data_grad_at_gauss_pt)
 
template<int R>
FTensor::Tensor1< double *, R > getValAtGaussPtsTensor (MatrixDouble &data)
 
template<int R, int D>
FTensor::Tensor2< double *, R, D > getGradAtGaussPtsTensor (MatrixDouble &data)
 
MoFEMErrorCode calculateValAndGrad (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Calculate gradient and values at integration points. More...
 
MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on left hand side. More...
 
template<>
FTensor::Tensor1< double *, 3 > getValAtGaussPtsTensor (MatrixDouble &data)
 Specialization for field with 3 coefficients in 3 dimension. More...
 
template<>
FTensor::Tensor2< double *, 3, 3 > getGradAtGaussPtsTensor (MatrixDouble &data)
 Specialization for field with 3 coefficients in 3 dimension. More...
 
template<>
MoFEMErrorCode calculateValAndGrad (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Specialization for field with 3 coefficients in 3 dimension. More...
 
template<>
MoFEMErrorCode calculateValAndGrad (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Specialization for field with for scalar field in 3 dimension. More...
 
template<>
FTensor::Tensor1< double *, 3 > getValAtGaussPtsTensor (MatrixDouble &data)
 
template<>
FTensor::Tensor2< double *, 3, 3 > getGradAtGaussPtsTensor (MatrixDouble &data)
 
template<>
MoFEMErrorCode calculateValAndGrad (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 
template<>
MoFEMErrorCode calculateValAndGrad (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true, const bool do_vertices=true, const bool do_edges=true, const bool do_quads=true, const bool do_tris=true, const bool do_tets=true, const bool do_prisms=true)
 
virtual ~DataOperator ()
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data, bool symm=true)
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data)
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool error_if_no_base=true)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Public Attributes

MatrixDoubledataAtGaussPts
 
MatrixDoubledataGradAtGaussPts
 
- Public Attributes inherited from MoFEM::DataOperator
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
bool doVertices
 If false skip vertices. More...
 
bool doEdges
 If false skip edges. More...
 
bool doQuads
 
bool doTris
 
bool doTets
 
bool doPrisms
 

Detailed Description

template<int RANK, int DIM>
struct MoFEM::OpGetDataAndGradient< RANK, DIM >

Get field values and gradients at Gauss points.

Definition at line 464 of file DataOperators.hpp.

Constructor & Destructor Documentation

◆ OpGetDataAndGradient()

template<int RANK, int DIM>
MoFEM::OpGetDataAndGradient< RANK, DIM >::OpGetDataAndGradient ( MatrixDouble data_at_gauss_pt,
MatrixDouble data_grad_at_gauss_pt 
)

Definition at line 469 of file DataOperators.hpp.

471  : dataAtGaussPts(data_at_gauss_pt),
472  dataGradAtGaussPts(data_grad_at_gauss_pt) {}

Member Function Documentation

◆ calculateValAndGrad() [1/5]

template<int RANK, int DIM>
MoFEMErrorCode MoFEM::OpGetDataAndGradient< RANK, DIM >::calculateValAndGrad ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)

Calculate gradient and values at integration points.

Parameters
sideside of entity on element
typetype of entity
datadata stored on entity (dofs values, dofs indices, etc.)
Returns
error code

Definition at line 497 of file DataOperators.hpp.

498  {
500  const int nb_base_functions = data.getN().size2();
501  bool constant_diff = false;
502  if (type == MBVERTEX && data.getDiffN().size1() * data.getDiffN().size2() ==
503  DIM * nb_base_functions) {
504  constant_diff = true;
505  }
506  const int nb_dofs = data.getFieldData().size();
507  for (unsigned int gg = 0; gg < data.getN().size1(); gg++) {
508  double *data_ptr, *n_ptr, *diff_n_ptr;
509  n_ptr = &data.getN()(gg, 0);
510  if (constant_diff) {
511  diff_n_ptr = &data.getDiffN()(0, 0);
512  } else {
513  diff_n_ptr = &data.getDiffN()(gg, 0);
514  }
515  data_ptr = &*data.getFieldData().data().begin();
516  for (int rr = 0; rr < RANK; rr++, data_ptr++) {
517  dataAtGaussPts(gg, rr) +=
518  cblas_ddot(nb_dofs / RANK, n_ptr, 1, data_ptr, RANK);
519  double *diff_n_ptr2 = diff_n_ptr;
520  for (unsigned int dd = 0; dd < DIM; dd++, diff_n_ptr2++) {
521  dataGradAtGaussPts(gg, DIM * rr + dd) +=
522  cblas_ddot(nb_dofs / RANK, diff_n_ptr2, DIM, data_ptr, RANK);
523  }
524  }
525  }
527  }
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T *> &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:499
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506
double cblas_ddot(const int N, const double *X, const int incX, const double *Y, const int incY)
Definition: cblas_ddot.c:12

◆ calculateValAndGrad() [2/5]

template<>
MoFEMErrorCode MoFEM::OpGetDataAndGradient< 3, 3 >::calculateValAndGrad ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)

Specialization for field with 3 coefficients in 3 dimension.

◆ calculateValAndGrad() [3/5]

template<>
MoFEMErrorCode MoFEM::OpGetDataAndGradient< 1, 3 >::calculateValAndGrad ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)

Specialization for field with for scalar field in 3 dimension.

◆ calculateValAndGrad() [4/5]

template<>
MoFEMErrorCode MoFEM::OpGetDataAndGradient< 3, 3 >::calculateValAndGrad ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)

Definition at line 1683 of file DataOperators.cpp.

1684  {
1686  if (data.getBase() == NOBASE)
1688  const unsigned int nb_gauss_pts = data.getN().size1();
1689  const unsigned int nb_base_functions = data.getN().size2();
1690  const unsigned int nb_dofs = data.getFieldData().size();
1691  if (!nb_dofs)
1693  auto t_n = data.getFTensor0N();
1694  auto t_val = getValAtGaussPtsTensor<3>(dataAtGaussPts);
1695  auto t_grad = getGradAtGaussPtsTensor<3, 3>(dataGradAtGaussPts);
1698  if (type == MBVERTEX &&
1699  data.getDiffN().data().size() == 3 * nb_base_functions) {
1700  for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
1701  auto t_data = data.getFTensor1FieldData<3>();
1702  auto t_diff_n = data.getFTensor1DiffN<3>();
1703  unsigned int bb = 0;
1704  for (; bb != nb_dofs / 3; ++bb) {
1705  t_val(i) += t_data(i) * t_n;
1706  t_grad(i, j) += t_data(i) * t_diff_n(j);
1707  ++t_n;
1708  ++t_diff_n;
1709  ++t_data;
1710  }
1711  ++t_val;
1712  ++t_grad;
1713  for (; bb != nb_base_functions; ++bb) {
1714  ++t_n;
1715  }
1716  }
1717  } else {
1718  auto t_diff_n = data.getFTensor1DiffN<3>();
1719  for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
1720  auto t_data = data.getFTensor1FieldData<3>();
1721  unsigned int bb = 0;
1722  for (; bb != nb_dofs / 3; ++bb) {
1723  t_val(i) += t_data(i) * t_n;
1724  t_grad(i, j) += t_data(i) * t_diff_n(j);
1725  ++t_n;
1726  ++t_diff_n;
1727  ++t_data;
1728  }
1729  ++t_val;
1730  ++t_grad;
1731  for (; bb != nb_base_functions; ++bb) {
1732  ++t_n;
1733  ++t_diff_n;
1734  }
1735  }
1736  }
1738 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:499
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506

◆ calculateValAndGrad() [5/5]

template<>
MoFEMErrorCode MoFEM::OpGetDataAndGradient< 1, 3 >::calculateValAndGrad ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)

Definition at line 1741 of file DataOperators.cpp.

1742  {
1744  const unsigned int nb_gauss_pts = data.getN().size1();
1745  const unsigned int nb_base_functions = data.getN().size2();
1746  // bool constant_diff = false;
1747  const unsigned int nb_dofs = data.getFieldData().size();
1748  auto t_n = data.getFTensor0N();
1750  FTensor::Tensor0<double *>(&*dataAtGaussPts.data().begin(), 1);
1751  double *ptr = &*dataGradAtGaussPts.data().begin();
1752  FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> t_grad(ptr, &ptr[1],
1753  &ptr[2]);
1755  if (type == MBVERTEX &&
1756  data.getDiffN().data().size() == 3 * nb_base_functions) {
1757  for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
1758  auto t_data = data.getFTensor0FieldData();
1759  auto t_diff_n = data.getFTensor1DiffN<3>();
1760  unsigned int bb = 0;
1761  for (; bb != nb_dofs / 3; ++bb) {
1762  t_val += t_data * t_n;
1763  t_grad(i) += t_data * t_diff_n(i);
1764  ++t_n;
1765  ++t_diff_n;
1766  ++t_data;
1767  }
1768  ++t_val;
1769  ++t_grad;
1770  for (; bb != nb_base_functions; ++bb) {
1771  ++t_n;
1772  }
1773  }
1774  } else {
1775  auto t_diff_n = data.getFTensor1DiffN<3>();
1776  for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
1777  auto t_data = data.getFTensor0FieldData();
1778  unsigned int bb = 0;
1779  for (; bb != nb_dofs / 3; ++bb) {
1780  t_val = t_data * t_n;
1781  t_grad(i) += t_data * t_diff_n(i);
1782  ++t_n;
1783  ++t_diff_n;
1784  ++t_data;
1785  }
1786  ++t_val;
1787  ++t_grad;
1788  for (; bb != nb_base_functions; ++bb) {
1789  ++t_n;
1790  ++t_diff_n;
1791  }
1792  }
1793  }
1795 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:499
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506

◆ doWork()

template<int RANK, int DIM>
MoFEMErrorCode MoFEM::OpGetDataAndGradient< RANK, DIM >::doWork ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)
virtual

Operator for linear form, usually to calculate values on left hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 529 of file DataOperators.hpp.

530  {
531 
533 
534  if (data.getFieldData().size() == 0) {
536  }
537 
538  unsigned int nb_dofs = data.getFieldData().size();
539  if (nb_dofs == 0)
541 
542  if (nb_dofs % RANK != 0) {
543  SETERRQ4(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
544  "data inconsistency, type %d, side %d, nb_dofs %d, rank %d",
545  type, side, nb_dofs, RANK);
546  }
547  if (nb_dofs / RANK > data.getN().size2()) {
548  std::cerr << side << " " << type << " "
549  << ApproximationBaseNames[data.getBase()] << std::endl;
550  std::cerr << data.getN() << std::endl;
551  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
552  "data inconsistency nb_dofs >= data.N.size2(), i.e. %u >= %u",
553  nb_dofs, data.getN().size2());
554  }
555 
556  if (type == MBVERTEX) {
557  dataAtGaussPts.resize(data.getN().size1(), RANK, false);
558  dataGradAtGaussPts.resize(data.getN().size1(), RANK * DIM, false);
559  dataAtGaussPts.clear();
560  dataGradAtGaussPts.clear();
561  }
562 
563  CHKERR calculateValAndGrad(side, type, data);
564 
566  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:475
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506
static const char *const ApproximationBaseNames[]
Definition: definitions.h:156
#define CHKERR
Inline error check.
Definition: definitions.h:594
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405
MoFEMErrorCode calculateValAndGrad(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Calculate gradient and values at integration points.

◆ getGradAtGaussPtsTensor() [1/3]

template<>
FTensor::Tensor2< double *, 3, 3 > MoFEM::OpGetDataAndGradient< 3, 3 >::getGradAtGaussPtsTensor< 3, 3 > ( MatrixDouble data)

Specialization for field with 3 coefficients in 3 dimension.

◆ getGradAtGaussPtsTensor() [2/3]

template<int RANK, int DIM>
template<int R, int D>
FTensor::Tensor2<double *, R, D> MoFEM::OpGetDataAndGradient< RANK, DIM >::getGradAtGaussPtsTensor ( MatrixDouble data)

Return tensor associated with matrix storing gradient values

Definition at line 486 of file DataOperators.hpp.

486  {
487  THROW_MESSAGE("Not implemented");
488  }
#define THROW_MESSAGE(a)
Throw MoFEM exception.
Definition: definitions.h:618

◆ getGradAtGaussPtsTensor() [3/3]

template<>
FTensor::Tensor2< double *, 3, 3 > MoFEM::OpGetDataAndGradient< 3, 3 >::getGradAtGaussPtsTensor< 3, 3 > ( MatrixDouble data)

Definition at line 1675 of file DataOperators.cpp.

1675  {
1676  double *ptr = &*data.data().begin();
1677  return FTensor::Tensor2<double *, 3, 3>(ptr, &ptr[1], &ptr[2], &ptr[3],
1678  &ptr[4], &ptr[5], &ptr[6], &ptr[7],
1679  &ptr[8], 9);
1680 }

◆ getValAtGaussPtsTensor() [1/3]

template<>
FTensor::Tensor1< double *, 3 > MoFEM::OpGetDataAndGradient< 3, 3 >::getValAtGaussPtsTensor< 3 > ( MatrixDouble data)

Specialization for field with 3 coefficients in 3 dimension.

◆ getValAtGaussPtsTensor() [2/3]

template<int RANK, int DIM>
template<int R>
FTensor::Tensor1<double *, R> MoFEM::OpGetDataAndGradient< RANK, DIM >::getValAtGaussPtsTensor ( MatrixDouble data)

Return tensor associated with matrix storing values

Definition at line 478 of file DataOperators.hpp.

478  {
479  THROW_MESSAGE("Not implemented");
480  }
#define THROW_MESSAGE(a)
Throw MoFEM exception.
Definition: definitions.h:618

◆ getValAtGaussPtsTensor() [3/3]

template<>
FTensor::Tensor1< double *, 3 > MoFEM::OpGetDataAndGradient< 3, 3 >::getValAtGaussPtsTensor< 3 > ( MatrixDouble data)

Definition at line 1667 of file DataOperators.cpp.

1667  {
1668  double *ptr = &*data.data().begin();
1669  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2], 3);
1670 }

Member Data Documentation

◆ dataAtGaussPts

template<int RANK, int DIM>
MatrixDouble& MoFEM::OpGetDataAndGradient< RANK, DIM >::dataAtGaussPts

Definition at line 466 of file DataOperators.hpp.

◆ dataGradAtGaussPts

template<int RANK, int DIM>
MatrixDouble& MoFEM::OpGetDataAndGradient< RANK, DIM >::dataGradAtGaussPts

Definition at line 467 of file DataOperators.hpp.


The documentation for this struct was generated from the following file: