v0.8.4
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 468 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 473 of file DataOperators.hpp.

476  :
477  dataAtGaussPts(data_at_gauss_pt),
478  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 503 of file DataOperators.hpp.

505  {
507  const int nb_base_functions = data.getN().size2();
508  bool constant_diff = false;
509  if(
510  type == MBVERTEX&&
511  data.getDiffN().size1()*data.getDiffN().size2()==DIM*nb_base_functions
512  ) {
513  constant_diff = true;
514  }
515  const int nb_dofs = data.getFieldData().size();
516  for(unsigned int gg = 0;gg<data.getN().size1();gg++) {
517  double *data_ptr,*n_ptr,*diff_n_ptr;
518  n_ptr = &data.getN()(gg,0);
519  if(constant_diff) {
520  diff_n_ptr = &data.getDiffN()(0,0);
521  } else {
522  diff_n_ptr = &data.getDiffN()(gg,0);
523  }
524  data_ptr = &*data.getFieldData().data().begin();
525  for(int rr = 0;rr<RANK;rr++,data_ptr++) {
526  dataAtGaussPts(gg,rr) += cblas_ddot(nb_dofs/RANK,n_ptr,1,data_ptr,RANK);
527  double *diff_n_ptr2 = diff_n_ptr;
528  for(unsigned int dd = 0;dd<DIM;dd++,diff_n_ptr2++) {
529  dataGradAtGaussPts(gg,DIM*rr+dd) += cblas_ddot(nb_dofs/RANK,diff_n_ptr2,DIM,data_ptr,RANK);
530  }
531  }
532  }
534  }
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:522
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
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 1672 of file DataOperators.cpp.

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

◆ calculateValAndGrad() [5/5]

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

Definition at line 1730 of file DataOperators.cpp.

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

◆ 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 536 of file DataOperators.hpp.

540  {
541 
543 
544  try {
545 
546  if(data.getFieldData().size() == 0) {
548  }
549 
550  unsigned int nb_dofs = data.getFieldData().size();
551  if(nb_dofs == 0) MoFEMFunctionReturnHot(0);
552 
553  if(nb_dofs % RANK != 0) {
554  SETERRQ4(PETSC_COMM_SELF,MOFEM_DATA_INCONSISTENCY,
555  "data inconsistency, type %d, side %d, nb_dofs %d, rank %d",
556  type,side,nb_dofs,RANK
557  );
558  }
559  if(nb_dofs/RANK > data.getN().size2()) {
560  std::cerr << side << " " << type << " " << ApproximationBaseNames[data.getBase()] << std::endl;
561  std::cerr << data.getN() << std::endl;
562  SETERRQ2(PETSC_COMM_SELF,MOFEM_DATA_INCONSISTENCY,
563  "data inconsistency nb_dofs >= data.N.size2(), i.e. %u >= %u",nb_dofs,data.getN().size2()
564  );
565  }
566 
567  if(type == MBVERTEX) {
568  dataAtGaussPts.resize(data.getN().size1(),RANK,false);
569  dataGradAtGaussPts.resize(data.getN().size1(),RANK*DIM,false);
570  dataAtGaussPts.clear();
571  dataGradAtGaussPts.clear();
572  }
573 
574  ierr = calculateValAndGrad(side,type,data); CHKERRG(ierr);
575 
576  } catch (std::exception& ex) {
577  std::ostringstream ss;
578  ss << "throw in method: " << ex.what() << " at line " << __LINE__ << " in file " << __FILE__;
579  SETERRQ(PETSC_COMM_SELF,MOFEM_STD_EXCEPTION_THROW,ss.str().c_str());
580  }
581 
583  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:522
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:565
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
static const char *const ApproximationBaseNames[]
Definition: definitions.h:153
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
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 492 of file DataOperators.hpp.

492  {
493  THROW_MESSAGE("Not implemented");
494  }
#define THROW_MESSAGE(a)
Throw MoFEM exception.
Definition: definitions.h:641

◆ getGradAtGaussPtsTensor() [3/3]

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

Definition at line 1664 of file DataOperators.cpp.

1664  {
1665  double *ptr = &*data.data().begin();
1666  return FTensor::Tensor2<double *, 3, 3>(ptr, &ptr[1], &ptr[2], &ptr[3],
1667  &ptr[4], &ptr[5], &ptr[6], &ptr[7],
1668  &ptr[8], 9);
1669 }

◆ 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 484 of file DataOperators.hpp.

484  {
485  THROW_MESSAGE("Not implemented");
486  }
#define THROW_MESSAGE(a)
Throw MoFEM exception.
Definition: definitions.h:641

◆ getValAtGaussPtsTensor() [3/3]

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

Definition at line 1656 of file DataOperators.cpp.

1656  {
1657  double *ptr = &*data.data().begin();
1658  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2], 3);
1659 }

Member Data Documentation

◆ dataAtGaussPts

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

Definition at line 470 of file DataOperators.hpp.

◆ dataGradAtGaussPts

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

Definition at line 471 of file DataOperators.hpp.


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