v0.8.23
Public Member Functions | Public Attributes | List of all members
FieldApproximationH1::OpApproxVolume< FUNEVAL > Struct Template Reference

Gauss point operators to calculate matrices and vectors. More...

#include <users_modules/basic_finite_elements/src/FieldApproximation.hpp>

Inheritance diagram for FieldApproximationH1::OpApproxVolume< FUNEVAL >:
[legend]
Collaboration diagram for FieldApproximationH1::OpApproxVolume< FUNEVAL >:
[legend]

Public Member Functions

 OpApproxVolume (const std::string &field_name, Mat _A, std::vector< Vec > &vec_F, FUNEVAL &function_evaluator)
 
virtual ~OpApproxVolume ()
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 calculate matrix More...
 
MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 calculate vector More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space)
 
 UserDataOperator (const std::string &field_name, const char type)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
double getMeasure ()
 get measure of element More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
MatrixDoublegetHoCoordsAtGaussPts ()
 coordinate at Gauss points (if hierarchical approximation of element geometry) More...
 
MatrixDoublegetHoGaussPtsJac ()
 
MatrixDoublegetHoGaussPtsInvJac ()
 
VectorDoublegetHoGaussPtsDetJac ()
 
auto getFTenosr0HoMeasure ()
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
auto getFTensor1HoCoordsAtGaussPts ()
 Get coordinates at integration points taking geometry from field. More...
 
auto getFTensor2HoGaussPtsJac ()
 
auto getFTensor2HoGaussPtsInvJac ()
 
const VolumeElementForcesAndSourcesCoregetVolumeFE ()
 return pointer to Generic Volume Finite Element object More...
 
MoFEMErrorCode getDivergenceOfHDivBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, VectorDouble &div)
 Get divergence of base functions at integration point. More...
 
MoFEMErrorCode getCurlOfHCurlBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, MatrixDouble &curl)
 Get curl of base functions at integration point. More...
 
DEPRECATED auto getTensor1CoordsAtGaussPts ()
 
DEPRECATED auto getTensor1HoCoordsAtGaussPts ()
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPLAST, const bool symm=true)
 
 UserDataOperator (const std::string &field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
virtual ~UserDataOperator ()
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
const std::string & getFEName () const
 Get name of the element. More...
 
Vec getSnesF () const
 
Vec getSnesX () const
 
Mat getSnesA () const
 
Mat getSnesB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTSa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
DEPRECATED MoFEMErrorCode getPorblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
DEPRECATED MoFEMErrorCode getPorblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
- 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 left 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 doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
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

Mat A
 
std::vector< Vec > & vecF
 
FUNEVAL & functionEvaluator
 
MatrixDouble NN
 
MatrixDouble transNN
 
std::vector< VectorDouble > Nf
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- 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
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType { OPROW = 1 << 0, OPCOL = 1 << 1, OPROWCOL = 1 << 2, OPLAST = 1 << 3 }
 Controls loop over entities on element. More...
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

template<typename FUNEVAL>
struct FieldApproximationH1::OpApproxVolume< FUNEVAL >

Gauss point operators to calculate matrices and vectors.

Function work on volumes (Terahedrons & Bricks)

Definition at line 44 of file FieldApproximation.hpp.

Constructor & Destructor Documentation

◆ OpApproxVolume()

template<typename FUNEVAL >
FieldApproximationH1::OpApproxVolume< FUNEVAL >::OpApproxVolume ( const std::string &  field_name,
Mat  _A,
std::vector< Vec > &  vec_F,
FUNEVAL &  function_evaluator 
)

Definition at line 51 of file FieldApproximation.hpp.

54  field_name, UserDataOperator::OPROW | UserDataOperator::OPROWCOL),
55  A(_A), vecF(vec_F), functionEvaluator(function_evaluator) {}
ForcesAndSourcesCore::UserDataOperator UserDataOperator

◆ ~OpApproxVolume()

template<typename FUNEVAL >
virtual FieldApproximationH1::OpApproxVolume< FUNEVAL >::~OpApproxVolume ( )
virtual

Definition at line 56 of file FieldApproximation.hpp.

56 {}

Member Function Documentation

◆ doWork() [1/2]

template<typename FUNEVAL >
MoFEMErrorCode FieldApproximationH1::OpApproxVolume< FUNEVAL >::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
DataForcesAndSourcesCore::EntData &  row_data,
DataForcesAndSourcesCore::EntData &  col_data 
)

calculate matrix

Definition at line 64 of file FieldApproximation.hpp.

67  {
69 
70  if (A == PETSC_NULL)
72  if (row_data.getIndices().size() == 0)
74  if (col_data.getIndices().size() == 0)
76 
77  const auto &dof_ptr = row_data.getFieldDofs()[0];
78  int rank = dof_ptr->getNbOfCoeffs();
79 
80  int nb_row_dofs = row_data.getIndices().size() / rank;
81  int nb_col_dofs = col_data.getIndices().size() / rank;
82 
83  NN.resize(nb_row_dofs, nb_col_dofs, false);
84  NN.clear();
85 
86  unsigned int nb_gauss_pts = row_data.getN().size1();
87  for (unsigned int gg = 0; gg != nb_gauss_pts; gg++) {
88  double w = getVolume() * getGaussPts()(3, gg);
89  if (getHoCoordsAtGaussPts().size1() == nb_gauss_pts) {
90  w *= getHoGaussPtsDetJac()[gg];
91  }
92  // noalias(NN) += w*outer_prod(row_data.getN(gg),col_data.getN(gg));
93  cblas_dger(CblasRowMajor, nb_row_dofs, nb_col_dofs, w,
94  &row_data.getN()(gg, 0), 1, &col_data.getN()(gg, 0), 1,
95  &*NN.data().begin(), nb_col_dofs);
96  }
97 
98  if ((row_type != col_type) || (row_side != col_side)) {
99  transNN.resize(nb_col_dofs, nb_row_dofs, false);
100  ublas::noalias(transNN) = trans(NN);
101  }
102 
103  double *data = &*NN.data().begin();
104  double *trans_data = &*transNN.data().begin();
105  VectorInt row_indices, col_indices;
106  row_indices.resize(nb_row_dofs);
107  col_indices.resize(nb_col_dofs);
108 
109  for (int rr = 0; rr < rank; rr++) {
110 
111  if ((row_data.getIndices().size() % rank) != 0) {
112  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
113  "data inconsistency");
114  }
115 
116  if ((col_data.getIndices().size() % rank) != 0) {
117  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
118  "data inconsistency");
119  }
120 
121  unsigned int nb_rows;
122  unsigned int nb_cols;
123  int *rows;
124  int *cols;
125 
126  if (rank > 1) {
127 
128  ublas::noalias(row_indices) = ublas::vector_slice<VectorInt>(
129  row_data.getIndices(),
130  ublas::slice(rr, rank, row_data.getIndices().size() / rank));
131  ublas::noalias(col_indices) = ublas::vector_slice<VectorInt>(
132  col_data.getIndices(),
133  ublas::slice(rr, rank, col_data.getIndices().size() / rank));
134 
135  nb_rows = row_indices.size();
136  nb_cols = col_indices.size();
137  rows = &*row_indices.data().begin();
138  cols = &*col_indices.data().begin();
139 
140  } else {
141 
142  nb_rows = row_data.getIndices().size();
143  nb_cols = col_data.getIndices().size();
144  rows = &*row_data.getIndices().data().begin();
145  cols = &*col_data.getIndices().data().begin();
146  }
147 
148  if (nb_rows != NN.size1()) {
149  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
150  "data inconsistency");
151  }
152  if (nb_cols != NN.size2()) {
153  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
154  "data inconsistency");
155  }
156 
157  CHKERR MatSetValues(A, nb_rows, rows, nb_cols, cols, data, ADD_VALUES);
158  if ((row_type != col_type) || (row_side != col_side)) {
159  if (nb_rows != transNN.size2()) {
160  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
161  "data inconsistency");
162  }
163  if (nb_cols != transNN.size1()) {
164  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
165  "data inconsistency");
166  }
167  CHKERR MatSetValues(A, nb_cols, cols, nb_rows, rows, trans_data,
168  ADD_VALUES);
169  }
170  }
171 
173  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
void cblas_dger(const enum CBLAS_ORDER order, const int M, const int N, const double alpha, const double *X, const int incX, const double *Y, const int incY, double *A, const int lda)
Definition: cblas_dger.c:12
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
MoFEMErrorCode MatSetValues(Mat M, const DataForcesAndSourcesCore::EntData &row_data, const DataForcesAndSourcesCore::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
#define CHKERR
Inline error check.
Definition: definitions.h:595
ublas::vector< int, IntAllocator > VectorInt
Definition: Types.hpp:75
MatrixDouble & getHoCoordsAtGaussPts()
coordinate at Gauss points (if hierarchical approximation of element geometry)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ doWork() [2/2]

template<typename FUNEVAL >
MoFEMErrorCode FieldApproximationH1::OpApproxVolume< FUNEVAL >::doWork ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData &  data 
)

calculate vector

Definition at line 177 of file FieldApproximation.hpp.

178  {
180 
181  if (data.getIndices().size() == 0)
183 
184  // PetscAttachDebugger();
185 
186  const auto &dof_ptr = data.getFieldDofs()[0];
187  unsigned int rank = dof_ptr->getNbOfCoeffs();
188 
189  int nb_row_dofs = data.getIndices().size() / rank;
190 
191  if (getCoordsAtGaussPts().size1() != data.getN().size1()) {
192  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
193  "data inconsistency");
194  }
195  if (getCoordsAtGaussPts().size2() != 3) {
196  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
197  "data inconsistency");
198  }
199 
200  // itegration
201  unsigned int nb_gauss_pts = data.getN().size1();
202  for (unsigned int gg = 0; gg != nb_gauss_pts; gg++) {
203 
204  double x, y, z, w;
205  w = getVolume() * getGaussPts()(3, gg);
206  if (getHoCoordsAtGaussPts().size1() == nb_gauss_pts) {
207  // intergation points global positions if higher order geometry is
208  // given
209  x = getHoCoordsAtGaussPts()(gg, 0);
210  y = getHoCoordsAtGaussPts()(gg, 1);
211  z = getHoCoordsAtGaussPts()(gg, 2);
212  // correction of jacobian for higher order geometry
213  w *= getHoGaussPtsDetJac()[gg];
214  } else {
215  // intergartion point global positions for linear tetrahedral element
216  x = getCoordsAtGaussPts()(gg, 0);
217  y = getCoordsAtGaussPts()(gg, 1);
218  z = getCoordsAtGaussPts()(gg, 2);
219  }
220 
221  std::vector<VectorDouble> fun_val;
222 
223  fun_val = functionEvaluator(x, y, z);
224 
225  if (fun_val.size() != vecF.size()) {
226  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
227  "data inconsistency");
228  }
229 
230  Nf.resize(fun_val.size());
231  for (unsigned int lhs = 0; lhs != fun_val.size(); lhs++) {
232 
233  if (!gg) {
234  Nf[lhs].resize(data.getIndices().size());
235  Nf[lhs].clear();
236  }
237 
238  if (fun_val[lhs].size() != rank) {
239  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
240  "data inconsistency");
241  }
242 
243  for (unsigned int rr = 0; rr != rank; rr++) {
244  cblas_daxpy(nb_row_dofs, w * (fun_val[lhs])[rr],
245  &data.getN()(gg, 0), 1, &(Nf[lhs])[rr], rank);
246  }
247  }
248  }
249 
250  for (unsigned int lhs = 0; lhs != vecF.size(); lhs++) {
251 
252  CHKERR VecSetValues(vecF[lhs], data.getIndices().size(),
253  &data.getIndices()[0], &(Nf[lhs])[0], ADD_VALUES);
254  }
255 
257  }
void cblas_daxpy(const int N, const double alpha, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_daxpy.c:11
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
MoFEMErrorCode VecSetValues(Vec V, const DataForcesAndSourcesCore::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
#define CHKERR
Inline error check.
Definition: definitions.h:595
MatrixDouble & getHoCoordsAtGaussPts()
coordinate at Gauss points (if hierarchical approximation of element geometry)
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

Member Data Documentation

◆ A

template<typename FUNEVAL >
Mat FieldApproximationH1::OpApproxVolume< FUNEVAL >::A

Definition at line 47 of file FieldApproximation.hpp.

◆ functionEvaluator

template<typename FUNEVAL >
FUNEVAL& FieldApproximationH1::OpApproxVolume< FUNEVAL >::functionEvaluator

Definition at line 49 of file FieldApproximation.hpp.

◆ Nf

template<typename FUNEVAL >
std::vector<VectorDouble> FieldApproximationH1::OpApproxVolume< FUNEVAL >::Nf

Definition at line 60 of file FieldApproximation.hpp.

◆ NN

template<typename FUNEVAL >
MatrixDouble FieldApproximationH1::OpApproxVolume< FUNEVAL >::NN

Definition at line 58 of file FieldApproximation.hpp.

◆ transNN

template<typename FUNEVAL >
MatrixDouble FieldApproximationH1::OpApproxVolume< FUNEVAL >::transNN

Definition at line 59 of file FieldApproximation.hpp.

◆ vecF

template<typename FUNEVAL >
std::vector<Vec>& FieldApproximationH1::OpApproxVolume< FUNEVAL >::vecF

Definition at line 48 of file FieldApproximation.hpp.


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