|
| v0.14.0
|
Go to the documentation of this file.
14 #ifndef __DATAOPERATORS_HPP
15 #define __DATAOPERATORS_HPP
17 using namespace boost::numeric;
31 DataOperator *op_ptr,
int row_side,
int col_side, EntityType row_type,
46 CHKERR doWorkLhsHook(
this, row_side, col_side, row_type, col_type,
50 "doWork function is not implemented for this operator");
74 "doWork function is not implemented for this operator");
80 const bool error_if_no_base =
false);
84 std::array<bool, MBMAXTYPE>
120 template <
bool ErrorIfNoBase>
122 const std::array<bool, MBMAXTYPE> &do_entities);
138 : tInvJac(&inv_jac(0, 0), &inv_jac(0, 1), &inv_jac(0, 2), &inv_jac(1, 0),
139 &inv_jac(1, 1), &inv_jac(1, 2), &inv_jac(2, 0), &inv_jac(2, 1),
161 : tInvJac(&inv_jac(0, 0), &inv_jac(0, 1), &inv_jac(0, 2), &inv_jac(1, 0),
162 &inv_jac(1, 1), &inv_jac(1, 2), &inv_jac(2, 0), &inv_jac(2, 1),
195 tJac(&jac(0, 0), &jac(0, 1), &jac(0, 2), &jac(1, 0), &jac(1, 1),
196 &jac(1, 2), &jac(2, 0), &jac(2, 1), &jac(2, 2)) {}
226 : tInvJac(&inv_jac(0, 0), &inv_jac(0, 1), &inv_jac(0, 2), &inv_jac(1, 0),
227 &inv_jac(1, 1), &inv_jac(1, 2), &inv_jac(2, 0), &inv_jac(2, 1),
247 : dataAtGaussPts(data_at_gauss_pt),
248 dataGradAtGaussPts(data_grad_at_gauss_pt) {}
261 template <
int R,
int D>
276 const int nb_base_functions = data.
getN().size2();
277 bool constant_diff =
false;
279 DIM * nb_base_functions) {
280 constant_diff =
true;
283 for (
unsigned int gg = 0; gg < data.
getN().size1(); gg++) {
284 double *data_ptr, *n_ptr, *diff_n_ptr;
285 n_ptr = &data.
getN()(gg, 0);
287 diff_n_ptr = &data.
getDiffN()(0, 0);
289 diff_n_ptr = &data.
getDiffN()(gg, 0);
292 for (
int rr = 0; rr < RANK; rr++, data_ptr++) {
293 dataAtGaussPts(gg, rr) +=
294 cblas_ddot(nb_dofs / RANK, n_ptr, 1, data_ptr, RANK);
295 double *diff_n_ptr2 = diff_n_ptr;
296 for (
unsigned int dd = 0;
dd < DIM;
dd++, diff_n_ptr2++) {
297 dataGradAtGaussPts(gg, DIM * rr +
dd) +=
298 cblas_ddot(nb_dofs / RANK, diff_n_ptr2, DIM, data_ptr, RANK);
318 if (nb_dofs % RANK != 0) {
320 "data inconsistency, type %d, side %d, nb_dofs %d, rank %d",
321 type, side, nb_dofs, RANK);
323 if (nb_dofs / RANK > data.
getN().size2()) {
324 std::cerr << side <<
" " <<
type <<
" "
326 std::cerr << data.
getN() << std::endl;
328 "data inconsistency nb_dofs >= data.N.size2(), i.e. %u >= %u",
329 nb_dofs, data.
getN().size2());
332 if (
type == MBVERTEX) {
333 dataAtGaussPts.resize(data.
getN().size1(), RANK,
false);
334 dataGradAtGaussPts.resize(data.
getN().size1(), RANK * DIM,
false);
335 dataAtGaussPts.clear();
336 dataGradAtGaussPts.clear();
351 OpGetDataAndGradient<3, 3>::getValAtGaussPtsTensor<3>(
MatrixDouble &data);
359 OpGetDataAndGradient<3, 3>::getGradAtGaussPtsTensor<3, 3>(
MatrixDouble &data);
394 : cOords_at_GaussPtF3(coords_at_gaussptf3),
395 nOrmals_at_GaussPtF3(normals_at_gaussptf3),
396 tAngent1_at_GaussPtF3(tangent1_at_gaussptf3),
397 tAngent2_at_GaussPtF3(tangent2_at_gaussptf3),
398 cOords_at_GaussPtF4(coords_at_gaussptf4),
399 nOrmals_at_GaussPtF4(normals_at_gaussptf4),
400 tAngent1_at_GaussPtF4(tangent1_at_gaussptf4),
401 tAngent2_at_GaussPtF4(tangent2_at_gaussptf4) {}
435 const int normal_shift = 0)
436 : normalRawPtr(&normal), normalsAtGaussPtsRawPtr(&normals_at_pts),
437 normalShift(normal_shift) {}
442 const int normal_shift = 0)
443 : normalRawPtr(normal_raw_ptr),
444 normalsAtGaussPtsRawPtr(normals_at_pts_ptr), normalShift(normal_shift) {
472 : nOrmal(normal), normalsAtGaussPts(normals_at_pts), tAngent0(tangent0),
473 tangent0AtGaussPt(tangent0_at_pts), tAngent1(tangent1),
474 tangent1AtGaussPt(tangent1_at_pts) {}
503 : tAngent(tangent), tangentAtGaussPt(tangent_at_pts) {}
511 #endif //__DATAOPERATORS_HPP
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Data on single entity (This is passed as argument to DataOperator::doWork)
void setSymm()
set if operator is executed taking in account symmetry
FTensor::Index< 'i', 3 > i
base operator to do operations at Gauss Pt. level
MatrixDouble & nOrmals_at_GaussPtF3
MatrixDouble & tAngent1_at_GaussPtF3
OpGetCoordsAndNormalsOnPrism(MatrixDouble &coords_at_gaussptf3, MatrixDouble &normals_at_gaussptf3, MatrixDouble &tangent1_at_gaussptf3, MatrixDouble &tangent2_at_gaussptf3, MatrixDouble &coords_at_gaussptf4, MatrixDouble &normals_at_gaussptf4, MatrixDouble &tangent1_at_gaussptf4, MatrixDouble &tangent2_at_gaussptf4)
MatrixDouble & tAngent2_at_GaussPtF3
MatrixDouble & dataAtGaussPts
const static char *const ApproximationBaseNames[]
MatrixDouble & tAngent2_at_GaussPtF4
Transform local reference derivatives of shape function to global derivatives.
MoFEMErrorCode calculateValAndGrad(int side, EntityType type, EntitiesFieldData::EntData &data)
Calculate gradient and values at integration points.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
MatrixDouble & cOords_at_GaussPtF3
UBlasMatrix< double > MatrixDouble
const VectorDouble & getFieldData() const
get dofs values
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
bool & doPrisms
\deprectaed
OpGetHOTangentOnEdge(MatrixDouble &tangent)
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
virtual MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
virtual MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
bool & doVertices
\deprectaed If false skip vertices
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
MatrixDouble & cOords_at_GaussPtF4
FTensor::Tensor2< double *, 3, 3 > tInvJac
FTensor::Index< 'k', 3 > k
#define CHKERR
Inline error check.
implementation of Data Operators for Forces and Sources
OpSetInvJacHdivAndHcurl(MatrixDouble3by3 &inv_jac)
FTensor::Index< 'i', 3 > i
brief Transform local reference derivatives of shape function to global derivatives
MatrixDouble diffHdivInvJac
FieldApproximationBase & getBase()
Get approximation base.
FTensor::Index< 'j', 3 > j
MatrixDouble & nOrmals_at_GaussPtF4
MatrixDouble & dataGradAtGaussPts
calculate normals at Gauss points of triangle element
bool & doQuads
\deprectaed
boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)> DoWorkLhsHookFunType
DoWorkRhsHookFunType doWorkRhsHook
EntitiesFieldData::EntData EntData
boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)> DoWorkRhsHookFunType
FTensor::Tensor2< double *, R, D > getGradAtGaussPtsTensor(MatrixDouble &data)
OpGetDataAndGradient(MatrixDouble &data_at_gauss_pt, MatrixDouble &data_grad_at_gauss_pt)
bool & doEdges
\deprectaed If false skip edges
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
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)
bool getSymm() const
Get if operator uses symmetry of DOFs or not.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
MatrixDouble & tAngent1_at_GaussPtF4
OpSetInvJacH1(MatrixDouble3by3 &inv_jac)
@ MOFEM_DATA_INCONSISTENCY
DoWorkLhsHookFunType doWorkLhsHook
FTensor::Tensor1< double *, R > getValAtGaussPtsTensor(MatrixDouble &data)
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
FTensor::Tensor2< double *, 3, 3 > tInvJac
data structure for finite element entity
Calculate tangent vector on edge form HO geometry approximation.
bool sYmm
If true assume that matrix is symmetric structure.
void unSetSymm()
unset if operator is executed for non symmetric problem
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Get field values and gradients at Gauss points.
FTensor::Index< 'j', 3 > j
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...