|
| v0.14.0
|
Go to the documentation of this file.
6 #ifndef __HO_DATA_OPERATORS_HPP__
7 #define __HO_DATA_OPERATORS_HPP__
29 boost::shared_ptr<MatrixDouble>
jacPtr;
40 template <
int FIELD_DIM = 3>
72 template <
int DIM = 3>
74 using OpSetHOInvJacToScalarBasesImpl::OpSetHOInvJacToScalarBasesImpl;
94 boost::shared_ptr<MatrixDouble> inv_jac_ptr)
99 "Pointer for invJacPtr not allocated");
144 template <
int SPACE_DIM>
166 "Pointer for detPtr not allocated");
184 boost::shared_ptr<VectorDouble> det_ptr,
185 boost::shared_ptr<MatrixDouble> jac_ptr)
210 boost::shared_ptr<MatrixDouble> jac_inv_ptr)
218 "Pointer for jacPtr not allocated");
280 template <
int FIELD_DIM = 3>
303 boost::shared_ptr<MatrixDouble> normals_at_gauss_pts =
nullptr)
322 boost::shared_ptr<MatrixDouble> tangent_at_pts =
nullptr)
341 boost::shared_ptr<MatrixDouble> normals_at_pts =
nullptr,
342 boost::shared_ptr<MatrixDouble> tangent1_at_pts =
nullptr,
343 boost::shared_ptr<MatrixDouble> tangent2_at_pts =
nullptr)
347 if (normals_at_pts || tangent1_at_pts || tangent2_at_pts)
348 if (normals_at_pts && tangent1_at_pts && tangent2_at_pts)
350 "All elements in constructor have to set pointer");
368 template <
int FIELD_DIM = 3>
374 boost::shared_ptr<MatrixDouble> tangents_at_pts =
nullptr)
396 const FieldSpace space, boost::shared_ptr<VectorDouble> det_jac_ptr);
413 template <
int FE_DIM,
int PROBLEM_DIM,
int SPACE_DIM>
struct AddHOOps;
418 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
419 std::vector<FieldSpace> spaces, std::string geom_field_name =
"",
420 boost::shared_ptr<MatrixDouble> jac =
nullptr,
421 boost::shared_ptr<VectorDouble> det =
nullptr,
422 boost::shared_ptr<MatrixDouble> inv_jac =
nullptr);
428 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
429 std::vector<FieldSpace> spaces, std::string geom_field_name =
"",
430 boost::shared_ptr<MatrixDouble> jac =
nullptr,
431 boost::shared_ptr<VectorDouble> det =
nullptr,
432 boost::shared_ptr<MatrixDouble> inv_jac =
nullptr);
438 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
439 std::vector<FieldSpace> spaces, std::string geom_field_name =
"");
444 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
445 std::vector<FieldSpace> space, std::string geom_field_name =
"");
451 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
452 std::vector<FieldSpace> space, std::string geom_field_name =
"");
458 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
459 std::vector<FieldSpace> space, std::string geom_field_name =
"",
460 boost::shared_ptr<MatrixDouble> jac =
nullptr,
461 boost::shared_ptr<VectorDouble> det =
nullptr,
462 boost::shared_ptr<MatrixDouble> inv_jac =
nullptr);
465 template <
int FIELD_DIM>
473 if (
type == MBVERTEX)
474 getCoordsAtGaussPts().clear();
476 auto t_coords = getFTensor1CoordsAtGaussPts();
477 const auto nb_integration_pts = data.
getN().size1();
478 const auto nb_base_functions = data.
getN().size2();
479 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
482 for (; bb != nb_dofs; ++bb) {
483 t_coords(
i) += t_base * t_dof(
i);
487 for (; bb < nb_base_functions; ++bb)
495 template <
int FIELD_DIM>
506 &
m(0, 0), &
m(0, 1), &
m(0, 2));
512 int nb_gauss_pts = data.
getN().size1();
513 auto &tangent1_at_gauss_pts = getTangent1AtGaussPts();
514 auto &tangent2_at_gauss_pts = getTangent2AtGaussPts();
515 tangent1_at_gauss_pts.resize(nb_gauss_pts, 3,
false);
516 tangent2_at_gauss_pts.resize(nb_gauss_pts, 3,
false);
520 tangent1_at_gauss_pts.clear();
521 tangent2_at_gauss_pts.clear();
528 if (2 * data.
getN().size2() != data.
getDiffN().size2()) {
530 "expected two direcatives in local element coordinates");
534 "expected that number of dofs is multiplicative of field "
541 for (; nn != nb_dofs; nn++) {
547 "Data inconsistency for base %s",
555 const int nb_base_functions = data.
getN().size2();
558 auto t_t1 = get_ftensor1(tangent1_at_gauss_pts);
559 auto t_t2 = get_ftensor1(tangent2_at_gauss_pts);
560 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
563 for (; bb != nb_dofs /
FIELD_DIM; ++bb) {
565 t_t1(
i) += t_data(
i) * t_diff_base(N0);
566 t_t2(
i) += t_data(
i) * t_diff_base(N1);
571 for (; bb != nb_base_functions; ++bb) {
584 if (moab::CN::Dimension(
type) == 2) {
586 auto &normal_at_gauss_pts = getNormalsAtGaussPts();
587 auto &tangent1_at_gauss_pts = getTangent1AtGaussPts();
588 auto &tangent2_at_gauss_pts = getTangent2AtGaussPts();
590 const auto nb_gauss_pts = tangent1_at_gauss_pts.size1();
591 normal_at_gauss_pts.resize(nb_gauss_pts, 3,
false);
593 auto t_normal = get_ftensor1(normal_at_gauss_pts);
594 auto t_t1 = get_ftensor1(tangent1_at_gauss_pts);
595 auto t_t2 = get_ftensor1(tangent2_at_gauss_pts);
596 for (
unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
610 template <
int FIELD_DIM>
618 return *tangentsAtPts;
620 return getTangentAtGaussPts();
623 auto &tangent = get_tangent();
624 int nb_gauss_pts = getGaussPts().size2();
626 if (
type == MBVERTEX) {
627 tangent.resize(nb_gauss_pts, 3,
false);
633 const int nb_base_functions = data.
getN().size2();
634 double *diff_base_function = &*data.
getDiffN().data().begin();
635 auto tangent_at_gauss_pts =
636 getFTensor1FromPtr<FIELD_DIM, 3>(&*tangent.data().begin());
642 "Inconsistent number of dofs and filed dimension");
644 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
647 for (; bb < size; ++bb) {
648 tangent_at_gauss_pts(
i) += field_data(
i) * (*diff_base_function);
650 ++diff_base_function;
652 for (; bb != nb_base_functions; ++bb) {
653 ++diff_base_function;
655 ++tangent_at_gauss_pts;
673 template <
typename E>
675 bool hcurl,
bool hdiv,
bool l2) {
676 std::vector<FieldSpace> spaces;
678 spaces.push_back(
H1);
680 spaces.push_back(
HCURL);
682 spaces.push_back(
HDIV);
684 spaces.push_back(
L2);
698 template <
typename E>
700 bool hcurl,
bool hdiv) {
701 std::vector<FieldSpace> spaces;
703 spaces.push_back(
HCURL);
705 spaces.push_back(
HDIV);
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
Data on single entity (This is passed as argument to DataOperator::doWork)
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 > detPtr
Calculate tangent vector on edge form HO geometry approximation.
@ OPSPACE
operator do Work is execute on space data
OpSetHOWeights(boost::shared_ptr< VectorDouble > det_ptr)
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
const static char *const ApproximationBaseNames[]
@ L2
field with C-1 continuity
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
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 > tangentsAtPts
const VectorDouble & getFieldData() const
get dofs values
OpCalculateHOJacForVolume(boost::shared_ptr< MatrixDouble > jac_ptr)
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
boost::shared_ptr< VectorDouble > detJacPtr
Scale base functions by inverses of measure of element.
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
bool & doVertices
\deprectaed If false skip vertices
Calculate HO coordinates at gauss points.
Calculate jacobian for face element.
OpSetInvJacToScalarBasesBasic(FieldSpace space, boost::shared_ptr< MatrixDouble > inv_jac_ptr)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
FieldSpace
approximation spaces
boost::shared_ptr< MatrixDouble > jacPtr
implementation of Data Operators for Forces and Sources
Calculate jacobian on Hex or other volume which is not simplex.
default operator for TRI element
Set inverse jacobian to base functions.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
FieldApproximationBase & getBase()
Get approximation base.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Modify integration weights on face to take in account higher-order geometry.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpSetHOInvJacVectorBase(const FieldSpace space, boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Calculate normals at Gauss points of triangle element.
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
friend class UserDataOperator
FTensor::Index< 'i', SPACE_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< MatrixDouble > invJacPtr
constexpr auto field_name
Add operators pushing bases from local to physical configuration.
default operator for EDGE element
structure to get information form mofem into EntitiesFieldData
boost::shared_ptr< MatrixDouble > jacPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
bool & doEdges
\deprectaed If false skip edges
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1FieldData()
Return FTensor of rank 1, i.e. vector from filed data coefficients.
OpGetHOTangentsOnEdge(std::string field_name, boost::shared_ptr< MatrixDouble > tangents_at_pts=nullptr)
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
MatrixDouble diffHdivInvJac
FTensor::Index< 'j', 3 > j
Set inverse jacobian to base functions.
@ HCURL
field with continuous tangents
@ MOFEM_DATA_INCONSISTENCY
OpGetHONormalsOnFace(std::string field_name)
FTensor::Index< 'm', 3 > m
transform local reference derivatives of shape function to global derivatives if higher order geometr...
Modify integration weights on face to take in account higher-order geometry.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
OpCalculateHOCoords(const std::string field_name)
FTensor::Index< 'k', 3 > k
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HDIV
field with continuous normal traction
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Transform local reference derivatives of shape functions to global derivatives.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpScaleBaseBySpaceInverseOfMeasure(const FieldSpace space, boost::shared_ptr< VectorDouble > det_jac_ptr)
@ OPROW
operator doWork function is executed on FE rows