v0.14.0
Public Types | Protected Attributes | Friends | List of all members
MoFEM::EntitiesFieldData::EntData Struct Reference

Data on single entity (This is passed as argument to DataOperator::doWork) More...

#include <src/finite_elements/EntitiesFieldData.hpp>

Inheritance diagram for MoFEM::EntitiesFieldData::EntData:
[legend]
Collaboration diagram for MoFEM::EntitiesFieldData::EntData:
[legend]

Public Types

enum  BaseDerivatives {
  ZeroDerivative = 0, FirstDerivative, SecondDerivative, ThirdDerivative,
  ForthDerivative, LastDerivative
}
 

Protected Attributes

int sEnse
 Entity sense (orientation) More...
 
ApproximationOrder oRder
 Entity order. More...
 
FieldSpace sPace
 Entity space. More...
 
FieldApproximationBase bAse
 Field approximation base. More...
 
VectorInt iNdices
 Global indices on entity. More...
 
VectorInt localIndices
 Local indices on entity. More...
 
VectorDofs dOfs
 DoFs on entity. More...
 
VectorFieldEntities fieldEntities
 Field entities. More...
 
VectorDouble fieldData
 Field data on entity. More...
 
std::vector< BitRefLevelentDataBitRefLevel
 Bit ref level in entity. More...
 
std::array< std::array< boost::shared_ptr< MatrixDouble >, LASTBASE >, LastDerivativebaseFunctionsAndBaseDerivatives
 
std::array< boost::shared_ptr< MatrixDouble >, LASTBASE > & N
 Base functions. More...
 
std::array< boost::shared_ptr< MatrixDouble >, LASTBASE > & diffN
 Derivatives of base functions. More...
 
std::string bbFieldName
 field name More...
 
VectorInt bbNodeOrder
 order of nodes More...
 
std::map< std::string, boost::shared_ptr< MatrixDouble > > bbN
 
std::map< std::string, boost::shared_ptr< MatrixDouble > > bbDiffN
 
std::map< std::string, boost::shared_ptr< MatrixInt > > bbAlphaIndices
 Indices for Bernstein-Bezier (BB) base. More...
 
std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrderbbNByOrder
 BB base functions by order. More...
 
std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrderbbDiffNByOrder
 BB base functions direvatives by order. More...
 
std::array< boost::shared_ptr< MatrixInt >, MaxBernsteinBezierOrderbbAlphaIndicesByOrder
 BB alpha indices by order. More...
 
boost::shared_ptr< MatrixDoubleswapBaseNPtr
 Used by Bernstein base to keep temporally pointer. More...
 
boost::shared_ptr< MatrixDoubleswapBaseDiffNPtr
 Used by Bernstein base to keep temporally pointer. More...
 

Friends

struct OpAddParentEntData
 

Specializations for tensor base function

template<>
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getFTensor1N (FieldApproximationBase base)
 
template<>
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getFTensor1N (FieldApproximationBase base, const int gg, const int bb)
 
template<>
FTensor::Tensor2< FTensor::PackPtr< double *, 9 >, 3, 3 > getFTensor2N (FieldApproximationBase base)
 

Specializations for direcatives of base functions

template<>
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getFTensor1DiffN (const FieldApproximationBase base)
 Get spatial derivative of base function tensor for dimension 3d. More...
 
template<>
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getFTensor1DiffN ()
 Get spatial derivative of base function tensor for dimension 3d. More...
 
template<>
FTensor::Tensor2< FTensor::PackPtr< double *, 6 >, 3, 2 > getFTensor2DiffN (FieldApproximationBase base)
 
template<>
FTensor::Tensor2< FTensor::PackPtr< double *, 6 >, 3, 2 > getFTensor2DiffN (FieldApproximationBase base, const int gg, const int bb)
 
template<>
FTensor::Tensor3< FTensor::PackPtr< double *, 12 >, 3, 2, 2 > getFTensor3Diff2N (FieldApproximationBase base)
 

Specializations for field data

template<>
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getFTensor1FieldData ()
 
template<>
FTensor::Tensor2< FTensor::PackPtr< double *, 1 >, 1, 1 > getFTensor2FieldData ()
 
template<>
FTensor::Tensor2_symmetric< FTensor::PackPtr< double *, 6 >, 3 > getFTensor2SymmetricFieldData ()
 

Constructor and destructor

 EntData (const bool allocate_base_matrices=true)
 
virtual ~EntData ()=default
 

Sense, order and indices

virtual int getSense () const
 get entity sense, need to calculate base functions with conforming approximation fields More...
 
ApproximationOrder getOrder () const
 get approximation order More...
 
const VectorIntgetIndices () const
 Get global indices of dofs on entity. More...
 
const VectorIntAdaptor getIndicesUpToOrder (int order)
 get global indices of dofs on entity up to given order More...
 
const VectorIntgetLocalIndices () const
 get local indices of dofs on entity More...
 
const VectorIntAdaptor getLocalIndicesUpToOrder (int order)
 get local indices of dofs on entity up to given order More...
 
int & getSense ()
 
ApproximationOrdergetOrder ()
 
VectorIntgetIndices ()
 
VectorIntgetLocalIndices ()
 

Data on entity

const VectorDoublegetFieldData () const
 get dofs values More...
 
const VectorAdaptor getFieldDataUpToOrder (int order)
 get dofs values up to given order More...
 
const VectorDofsgetFieldDofs () const
 get dofs data stature FEDofEntity More...
 
VectorDoublegetFieldData ()
 get dofs data stature FEDofEntity More...
 
const VectorFieldEntitiesgetFieldEntities () const
 get field entities More...
 
VectorFieldEntitiesgetFieldEntities ()
 get field entities More...
 
virtual std::vector< BitRefLevel > & getEntDataBitRefLevel ()
 
template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1FieldData ()
 Return FTensor of rank 1, i.e. vector from filed data coeffinects. More...
 
template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2FieldData ()
 Return FTensor rank 2, i.e. matrix from filed data coeffinects. More...
 
template<int Tensor_Dim>
FTensor::Tensor2_symmetric< FTensor::PackPtr< double *,(Tensor_Dim *(Tensor_Dim+1))/2 >, Tensor_Dim > getFTensor2SymmetricFieldData ()
 Return symmetric FTensor rank 2, i.e. matrix from filed data coeffinects. More...
 
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData ()
 Resturn scalar files as a FTensor of rank 0. More...
 
VectorDofsgetFieldDofs ()
 

Base and space

FieldApproximationBasegetBase ()
 Get approximation base. More...
 
FieldSpacegetSpace ()
 Get field space. More...
 
virtual boost::shared_ptr< MatrixDouble > & getNSharedPtr (const FieldApproximationBase base, const BaseDerivatives direvatie)
 
virtual boost::shared_ptr< MatrixDouble > & getNSharedPtr (const FieldApproximationBase base)
 
virtual boost::shared_ptr< MatrixDouble > & getDiffNSharedPtr (const FieldApproximationBase base)
 

Get base functions for H1/L2

MatrixDoublegetN (const FieldApproximationBase base)
 get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb. of columns is equal to number of base functions on this entity. More...
 
MatrixDoublegetN (const std::string &field_name)
 get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb. of columns is equal to number of base functions on this entity. More...
 
MatrixDoublegetN ()
 get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb. of columns is equal to number of base functions on this entity. More...
 
MatrixDoublegetDiffN (const FieldApproximationBase base)
 get derivatives of base functions More...
 
MatrixDoublegetN (const FieldApproximationBase base, const BaseDerivatives derivative)
 Get base function derivative. More...
 
MatrixDoublegetDiffN (const std::string &field_name)
 get derivatives of base functions More...
 
MatrixDoublegetDiffN ()
 get derivatives of base functions More...
 
MatrixDoublegetN (const BaseDerivatives derivative)
 Get base function derivative. More...
 
const VectorAdaptor getN (const FieldApproximationBase base, const int gg)
 get base functions at Gauss pts More...
 
const VectorAdaptor getN (const int gg)
 get base functions at Gauss pts More...
 
const MatrixAdaptor getDiffN (const FieldApproximationBase base, const int gg)
 get derivative of base functions at Gauss pts More...
 
const MatrixAdaptor getDiffN (const int gg)
 get derivative of base functions at Gauss pts More...
 
const VectorAdaptor getN (const FieldApproximationBase base, const int gg, const int nb_base_functions)
 get base functions at Gauss pts More...
 
const VectorAdaptor getN (const int gg, const int nb_base_functions)
 get base functions at Gauss pts More...
 
const MatrixAdaptor getDiffN (const FieldApproximationBase base, const int gg, const int nb_base_functions)
 get derivatives of base functions at Gauss pts More...
 
const MatrixAdaptor getDiffN (const int gg, const int nb_base_functions)
 get derivatives of base functions at Gauss pts More...
 

Get base functions for vectorial approximation basese, i.e.

Hdiv/Hcurl

template<int DIM>
const MatrixAdaptor getVectorN (const FieldApproximationBase base, const int gg)
 get Hdiv of base functions at Gauss pts More...
 
template<int DIM>
const MatrixAdaptor getVectorN (const int gg)
 get Hdiv of base functions at Gauss pts More...
 
template<int DIM0, int DIM1>
const MatrixAdaptor getVectorDiffN (FieldApproximationBase base, const int gg)
 get DiffHdiv of base functions at Gauss pts More...
 
template<int DIM0, int DIM1>
const MatrixAdaptor getVectorDiffN (const int gg)
 get DiffHdiv of base functions at Gauss pts More...
 
template<int DIM0, int DIM1>
const MatrixAdaptor getVectorDiffN (const FieldApproximationBase base, const int dof, const int gg)
 get DiffHdiv of base functions at Gauss pts More...
 
template<int DIM0, int DIM1>
const MatrixAdaptor getVectorDiffN (const int dof, const int gg)
 get DiffHdiv of base functions at Gauss pts More...
 

Get base functions with FTensor

FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N (const FieldApproximationBase base)
 Get base function as Tensor0. More...
 
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N ()
 Get base function as Tensor0. More...
 
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N (const FieldApproximationBase base, const int gg, const int bb)
 Get base function as Tensor0 (Loop by integration points) More...
 
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N (const int gg, const int bb)
 Get base function as Tensor0 (Loop by integration points) More...
 
template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN (const FieldApproximationBase base)
 Get derivatives of base functions. More...
 
template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN ()
 Get derivatives of base functions. More...
 
template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN (const FieldApproximationBase base, const int gg, const int bb)
 Get derivatives of base functions (Loop by integration points) More...
 
template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN (const int gg, const int bb)
 Get derivatives of base functions (Loop by integration points) More...
 
template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N (FieldApproximationBase base)
 Get base functions for Hdiv/Hcurl spaces. More...
 
template<int Tensor_Dim>
auto getFTensor1N ()
 Get base functions for Hdiv space. More...
 
template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN (FieldApproximationBase base)
 Get derivatives of base functions for Hdiv space. More...
 
template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN (FieldApproximationBase base, const int gg, const int bb)
 Get derivatives of base functions for Hdiv space at integration pts. More...
 
template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN ()
 Get derivatives of base functions for Hdiv space. More...
 
template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN (const int gg, const int bb)
 Get derivatives of base functions for Hdiv space at integration pts. More...
 
template<int Tensor_Dim0, int Tensor_Dim1, int Tensor_Dim2>
FTensor::Tensor3< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 *Tensor_Dim2 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3Diff2N (FieldApproximationBase base)
 Get second derivatives of base functions for Hvec space. More...
 
template<int Tensor_Dim0, int Tensor_Dim1, int Tensor_Dim2>
FTensor::Tensor3< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 *Tensor_Dim2 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3Diff2N ()
 Get second derivatives of base functions for Hvec space. More...
 
template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N (FieldApproximationBase base, const int gg, const int bb)
 Get Hdiv base functions at integration point. More...
 
template<int Tensor_Dim>
auto getFTensor1N (const int gg, const int bb)
 Get Hdiv base functions at integration point. More...
 
template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2N (FieldApproximationBase base)
 Get base functions for Hdiv/Hcurl spaces. More...
 
template<int Tensor_Dim0, int Tensor_Dim1>
auto getFTensor2N ()
 Get base functions for Hdiv space. More...
 
template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2N (FieldApproximationBase base, const int gg, const int bb)
 Get base functions for tensor Hdiv/Hcurl spaces. More...
 
template<int Tensor_Dim0, int Tensor_Dim1>
auto getFTensor2N (const int gg, const int bb)
 Get base functions for Hdiv space. More...
 

Auxiliary functions

std::ostream & operator<< (std::ostream &os, const EntitiesFieldData::EntData &e)
 
MoFEMErrorCode resetFieldDependentData ()
 

Bernstein-Bezier base only functions

static constexpr size_t MaxBernsteinBezierOrder = BITFEID_SIZE
 
VectorIntgetBBNodeOrder ()
 Get orders at the nodes. More...
 
MatrixIntgetBBAlphaIndices ()
 Get file BB indices. More...
 
virtual boost::shared_ptr< MatrixInt > & getBBAlphaIndicesSharedPtr (const std::string &field_name)
 
virtual boost::shared_ptr< MatrixDouble > & getBBNSharedPtr (const std::string &field_name)
 
virtual const boost::shared_ptr< MatrixDouble > & getBBNSharedPtr (const std::string &field_name) const
 
virtual boost::shared_ptr< MatrixDouble > & getBBDiffNSharedPtr (const std::string &field_name)
 
virtual const boost::shared_ptr< MatrixDouble > & getBBDiffNSharedPtr (const std::string &field_name) const
 
virtual std::map< std::string, boost::shared_ptr< MatrixInt > > & getBBAlphaIndicesMap ()
 
virtual std::map< std::string, boost::shared_ptr< MatrixDouble > > & getBBNMap ()
 get hash map of base function for BB base, key is a field name More...
 
virtual std::map< std::string, boost::shared_ptr< MatrixDouble > > & getBBDiffNMap ()
 get hash map of direvarives base function for BB base, key is a field name More...
 
virtual boost::shared_ptr< MatrixInt > & getBBAlphaIndicesByOrderSharedPtr (const size_t o)
 get ALpha indices for BB base by order More...
 
virtual boost::shared_ptr< MatrixDouble > & getBBNByOrderSharedPtr (const size_t o)
 get BB base by order More...
 
virtual boost::shared_ptr< MatrixDouble > & getBBDiffNByOrderSharedPtr (const size_t o)
 get BB base derivative by order More...
 
virtual std::array< boost::shared_ptr< MatrixInt >, MaxBernsteinBezierOrder > & getBBAlphaIndicesByOrderArray ()
 
virtual std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > & getBBNByOrderArray ()
 
virtual std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > & getBBDiffNByOrderArray ()
 
virtual MoFEMErrorCode baseSwap (const std::string &field_name, const FieldApproximationBase base)
 Swap bases functions. More...
 

Specializations for H1/L2

template<>
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getFTensor1DiffN (const FieldApproximationBase base, const int gg, const int bb)
 Get spatial derivative of base function tensor for dimension 3d. More...
 
template<>
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getFTensor1DiffN (const int gg, const int bb)
 Get spatial derivative of base function tensor for dimension 3d. More...
 

Specializations for HDiv/HCrul

template<>
FTensor::Tensor2< FTensor::PackPtr< double *, 9 >, 3, 3 > getFTensor2N (FieldApproximationBase base, const int gg, const int bb)
 

Detailed Description

Data on single entity (This is passed as argument to DataOperator::doWork)

Todo:
Hdiv and Hcurl functions should be accessed through common interface.
Examples
adolc_plasticity.cpp, approx_sphere.cpp, boundary_marker.cpp, build_large_problem.cpp, child_and_parent.cpp, ContactOps.hpp, continuity_check_on_contact_prism_side_ele.cpp, continuity_check_on_skeleton_3d.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, dg_projection.cpp, dynamic_first_order_con_law.cpp, eigen_elastic.cpp, EshelbianOperators.cpp, EshelbianPlasticity.cpp, field_evaluator.cpp, forces_and_sources_testing_edge_element.cpp, forces_and_sources_testing_flat_prism_element.cpp, forces_and_sources_testing_users_base.cpp, free_surface.cpp, hanging_node_approx.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hcurl_divergence_operator_2d.cpp, hdiv_divergence_operator.cpp, heat_equation.cpp, heat_method.cpp, hello_world.cpp, higher_derivatives.cpp, HookeElement.cpp, HookeElement.hpp, HookeInternalStressElement.hpp, level_set.cpp, lorentz_force.cpp, mesh_smoothing.cpp, mixed_poisson.cpp, NavierStokesElement.cpp, NavierStokesElement.hpp, nonlinear_elastic.cpp, operators_tests.cpp, phase.cpp, photon_diffusion.cpp, plastic.cpp, PlasticOpsGeneric.hpp, plate.cpp, plot_base.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.hpp, PoissonDiscontinousGalerkin.hpp, prism_elements_from_surface.cpp, prism_polynomial_approximation.cpp, quad_polynomial_approximation.cpp, reaction_diffusion.cpp, scalar_check_approximation.cpp, seepage.cpp, shallow_wave.cpp, simple_contact_thermal.cpp, simple_elasticity.cpp, simple_interface.cpp, simple_l2_only.cpp, tensor_divergence_operator.cpp, test_cache_on_entities.cpp, thermo_elastic.cpp, ThermoElasticOps.hpp, and wave_equation.cpp.

Definition at line 127 of file EntitiesFieldData.hpp.

Member Enumeration Documentation

◆ BaseDerivatives

Enumerator
ZeroDerivative 
FirstDerivative 
SecondDerivative 
ThirdDerivative 
ForthDerivative 
LastDerivative 

Definition at line 129 of file EntitiesFieldData.hpp.

129  {
130  ZeroDerivative = 0,
136  };

Constructor & Destructor Documentation

◆ EntData()

MoFEM::EntitiesFieldData::EntData::EntData ( const bool  allocate_base_matrices = true)

Definition at line 10 of file EntitiesFieldData.cpp.

14  if (allocate_base_matrices) {
15 
16  for (auto d = 0; d != LastDerivative; ++d) {
17  for (int b = 0; b != LASTBASE; ++b) {
19  }
20  }
21  }
22 }

◆ ~EntData()

virtual MoFEM::EntitiesFieldData::EntData::~EntData ( )
virtualdefault

Member Function Documentation

◆ baseSwap()

MoFEMErrorCode MoFEM::EntitiesFieldData::EntData::baseSwap ( const std::string &  field_name,
const FieldApproximationBase  base 
)
virtual

Swap bases functions.

Some base are not hierarchical and depene on approximation order. Such case demand special handling, that appropiate base order is set depending on field, such that is accessible in operator.

Note
Base is not swap on meshsets
Parameters
field_name
base
Returns
MoFEMErrorCode

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 136 of file EntitiesFieldData.cpp.

137  {
139  auto make_swap = [](boost::shared_ptr<MatrixDouble> &ptr,
140  boost::shared_ptr<MatrixDouble> &ptrBB,
141  boost::shared_ptr<MatrixDouble> &swap_ptr) {
142  if (swap_ptr) {
143  ptr = swap_ptr;
144  swap_ptr.reset();
145  } else {
146  swap_ptr = ptr;
147  ptr = ptrBB;
148  }
149  };
154 }

◆ getBase()

FieldApproximationBase & MoFEM::EntitiesFieldData::EntData::getBase ( )
inline

Get approximation base.

Returns
Approximation base
Examples
hanging_node_approx.cpp.

Definition at line 1300 of file EntitiesFieldData.hpp.

1300 { return bAse; }

◆ getBBAlphaIndices()

MatrixInt & MoFEM::EntitiesFieldData::EntData::getBBAlphaIndices ( )
inline

Get file BB indices.

Returns
MatrixInt&

Definition at line 1536 of file EntitiesFieldData.hpp.

1536  {
1538 }

◆ getBBAlphaIndicesByOrderArray()

std::array< boost::shared_ptr< MatrixInt >, EntitiesFieldData::EntData::MaxBernsteinBezierOrder > & MoFEM::EntitiesFieldData::EntData::getBBAlphaIndicesByOrderArray ( )
virtual

Definition at line 890 of file EntitiesFieldData.cpp.

890  {
891  return bbAlphaIndicesByOrder;
892 }

◆ getBBAlphaIndicesByOrderSharedPtr()

boost::shared_ptr< MatrixInt > & MoFEM::EntitiesFieldData::EntData::getBBAlphaIndicesByOrderSharedPtr ( const size_t  o)
virtual

get ALpha indices for BB base by order

Parameters
oapproximation order
Returns
boost::shared_ptr<MatrixInt>&

Definition at line 874 of file EntitiesFieldData.cpp.

874  {
875  return bbAlphaIndicesByOrder[o];
876 }

◆ getBBAlphaIndicesMap()

std::map< std::string, boost::shared_ptr< MatrixInt > > & MoFEM::EntitiesFieldData::EntData::getBBAlphaIndicesMap ( )
virtual

Definition at line 859 of file EntitiesFieldData.cpp.

859  {
860  return bbAlphaIndices;
861 }

◆ getBBAlphaIndicesSharedPtr()

boost::shared_ptr< MatrixInt > & MoFEM::EntitiesFieldData::EntData::getBBAlphaIndicesSharedPtr ( const std::string &  field_name)
virtual

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 822 of file EntitiesFieldData.cpp.

823  {
824  return bbAlphaIndices[field_name];
825 }

◆ getBBDiffNByOrderArray()

std::array< boost::shared_ptr< MatrixDouble >, EntitiesFieldData::EntData::MaxBernsteinBezierOrder > & MoFEM::EntitiesFieldData::EntData::getBBDiffNByOrderArray ( )
virtual

Definition at line 902 of file EntitiesFieldData.cpp.

902  {
903  return bbDiffNByOrder;
904 }

◆ getBBDiffNByOrderSharedPtr()

boost::shared_ptr< MatrixDouble > & MoFEM::EntitiesFieldData::EntData::getBBDiffNByOrderSharedPtr ( const size_t  o)
virtual

get BB base derivative by order

Parameters
o
Returns
boost::shared_ptr<MatrixDouble>&

Definition at line 884 of file EntitiesFieldData.cpp.

884  {
885  return bbDiffNByOrder[o];
886 }

◆ getBBDiffNMap()

std::map< std::string, boost::shared_ptr< MatrixDouble > > & MoFEM::EntitiesFieldData::EntData::getBBDiffNMap ( )
virtual

get hash map of direvarives base function for BB base, key is a field name

Returns
std::map<std::string, boost::shared_ptr<MatrixDouble>>&

Definition at line 869 of file EntitiesFieldData.cpp.

869  {
870  return bbDiffN;
871 }

◆ getBBDiffNSharedPtr() [1/2]

boost::shared_ptr< MatrixDouble > & MoFEM::EntitiesFieldData::EntData::getBBDiffNSharedPtr ( const std::string &  field_name)
virtual

Get shared pointer to BB derivatives of base base functions

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 845 of file EntitiesFieldData.cpp.

845  {
846  return bbDiffN[field_name];
847 }

◆ getBBDiffNSharedPtr() [2/2]

const boost::shared_ptr< MatrixDouble > & MoFEM::EntitiesFieldData::EntData::getBBDiffNSharedPtr ( const std::string &  field_name) const
virtual

Get shared pointer to derivatives of BB base base functions

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 853 of file EntitiesFieldData.cpp.

854  {
855  return bbDiffN.at(field_name);
856 }

◆ getBBNByOrderArray()

std::array< boost::shared_ptr< MatrixDouble >, EntitiesFieldData::EntData::MaxBernsteinBezierOrder > & MoFEM::EntitiesFieldData::EntData::getBBNByOrderArray ( )
virtual

Definition at line 896 of file EntitiesFieldData.cpp.

896  {
897  return bbNByOrder;
898 }

◆ getBBNByOrderSharedPtr()

boost::shared_ptr< MatrixDouble > & MoFEM::EntitiesFieldData::EntData::getBBNByOrderSharedPtr ( const size_t  o)
virtual

get BB base by order

Parameters
o
Returns
boost::shared_ptr<MatrixDouble>&

Definition at line 879 of file EntitiesFieldData.cpp.

879  {
880  return bbNByOrder[o];
881 }

◆ getBBNMap()

std::map< std::string, boost::shared_ptr< MatrixDouble > > & MoFEM::EntitiesFieldData::EntData::getBBNMap ( )
virtual

get hash map of base function for BB base, key is a field name

Returns
std::map<std::string, boost::shared_ptr<MatrixDouble>>&

Definition at line 864 of file EntitiesFieldData.cpp.

864  {
865  return bbN;
866 }

◆ getBBNodeOrder()

VectorInt & MoFEM::EntitiesFieldData::EntData::getBBNodeOrder ( )
inline

Get orders at the nodes.

Returns
VectorInt&

Definition at line 1534 of file EntitiesFieldData.hpp.

1534 { return bbNodeOrder; }

◆ getBBNSharedPtr() [1/2]

boost::shared_ptr< MatrixDouble > & MoFEM::EntitiesFieldData::EntData::getBBNSharedPtr ( const std::string &  field_name)
virtual

Get shared pointer to BB base base functions

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 828 of file EntitiesFieldData.cpp.

828  {
829  return bbN[field_name];
830 }

◆ getBBNSharedPtr() [2/2]

const boost::shared_ptr< MatrixDouble > & MoFEM::EntitiesFieldData::EntData::getBBNSharedPtr ( const std::string &  field_name) const
virtual

Get shared pointer to BB base base functions

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 836 of file EntitiesFieldData.cpp.

837  {
838  return bbN.at(field_name);
839 }

◆ getDiffN() [1/7]

MatrixDouble & MoFEM::EntitiesFieldData::EntData::getDiffN ( )
inline

get derivatives of base functions

Matrix at rows has nb. of Gauss pts, at columns it has derivative of base functions. Columns are structured as follows, [ dN1/dx, dN1/dy, dN1/dz, dN2/dx, dN2/dy, dN2/dz, ... ]

Scalar base functions: Note that base functions are calculated in file H1.c Above description not apply for derivatives of nodal functions, since derivative of nodal functions in case of simplexes, EDGES, TRIANGLES and TETS are constant. So that matrix rows represents nb. of base functions, columns are derivatives. Nb. of columns depend on element dimension, for EDGES is one, for TRIS is 2 and TETS is 3.

Note
Note that for node element this function make no sense.

Tonsorial base functions:

Note
Note: In rows ale integration pts, columns are formatted that that components of vectors and then derivatives, for example row for given integration points is formatted in array

\[ t_{0,0}, t_{1,0}, t_{1,0}, t_{0,1}, t_{1,1}, t_{1,1}, t_{0,2}, t_{1,2}, t_{1,2} \]

where comma express derivative, i.e. \(t_{2,1} = \frac{\partial t_2}{\partial \xi_1}\)

Definition at line 1339 of file EntitiesFieldData.hpp.

1339 { return getDiffN(bAse); }

◆ getDiffN() [2/7]

MatrixDouble & MoFEM::EntitiesFieldData::EntData::getDiffN ( const FieldApproximationBase  base)
inline

get derivatives of base functions

Matrix at rows has nb. of Gauss pts, at columns it has derivative of base functions. Columns are structured as follows, [ dN1/dx, dN1/dy, dN1/dz, dN2/dx, dN2/dy, dN2/dz, ... ]

Scalar base functions: Note that base functions are calculated in file H1.c Above description not apply for derivatives of nodal functions, since derivative of nodal functions in case of simplexes, EDGES, TRIANGLES and TETS are constant. So that matrix rows represents nb. of base functions, columns are derivatives. Nb. of columns depend on element dimension, for EDGES is one, for TRIS is 2 and TETS is 3.

Note
Note that for node element this function make no sense.

Tonsorial base functions:

Note
Note: In rows ale integration pts, columns are formatted that that components of vectors and then derivatives, for example row for given integration points is formatted in array

\[ t_{0,0}, t_{1,0}, t_{1,0}, t_{0,1}, t_{1,1}, t_{1,1}, t_{0,2}, t_{1,2}, t_{1,2} \]

where comma express derivative, i.e. \(t_{2,1} = \frac{\partial t_2}{\partial \xi_1}\)
Examples
hcurl_curl_operator.cpp, hdiv_divergence_operator.cpp, PoissonDiscontinousGalerkin.hpp, and scalar_check_approximation.cpp.

Definition at line 1316 of file EntitiesFieldData.hpp.

1316  {
1317  return *(getDiffNSharedPtr(base));
1318 }

◆ getDiffN() [3/7]

const MatrixAdaptor MoFEM::EntitiesFieldData::EntData::getDiffN ( const FieldApproximationBase  base,
const int  gg 
)
inline

get derivative of base functions at Gauss pts

returned matrix on rows has base functions, in column its derivatives.

Parameters
baseApproximation base
ggNb. of Gauss pts.

Definition at line 1359 of file EntitiesFieldData.hpp.

1360  {
1361  // FIXME: That is bug, it will not work if number of integration pts is
1362  // equal to number of nodes on entity. User who not implementing low
1363  // level DataOperator will not experience this.
1364  if (getN(base).size1() == getDiffN(base).size1()) {
1365  int size = getN(base).size2();
1366  int dim = getDiffN(base).size2() / size;
1367  double *data = &getDiffN(base)(gg, 0);
1368  return MatrixAdaptor(
1369  getN(base).size2(), dim,
1370  ublas::shallow_array_adaptor<double>(getDiffN(base).size2(), data));
1371  } else {
1372  // in some cases, f.e. for derivatives of nodal base functions at only
1373  // one gauss point is needed
1374  return MatrixAdaptor(
1375  getN(base).size1(), getN(base).size2(),
1376  ublas::shallow_array_adaptor<double>(getDiffN(base).data().size(),
1377  &getDiffN(base).data()[0]));
1378  }
1379 }

◆ getDiffN() [4/7]

const MatrixAdaptor MoFEM::EntitiesFieldData::EntData::getDiffN ( const FieldApproximationBase  base,
const int  gg,
const int  nb_base_functions 
)
inline

get derivatives of base functions at Gauss pts

Note that multi field element, two different field can have different approximation orders. Since we use hierarchical approximation basis, base functions are calculated once for element, using maximal approximation order on given entity.

Parameters
baseApproximation base
ggnb. of Gauss point
nb_base_functionsnumber of of base functions

Definition at line 1401 of file EntitiesFieldData.hpp.

1403  {
1404  // FIXME: That is bug, it will not work if number of integration pts is
1405  // equal to number of nodes on entity. User who not implementing low
1406  // level DataOperator will not experience this.
1407  if (getN(base).size1() == getDiffN(base).size1()) {
1408  (void)getN(base)(gg,
1409  nb_base_functions -
1410  1); // throw error if nb_base_functions is to big
1411  int dim = getDiffN(base).size2() / getN(base).size2();
1412  double *data = &getDiffN(base)(gg, 0);
1413  return MatrixAdaptor(
1414  nb_base_functions, dim,
1415  ublas::shallow_array_adaptor<double>(dim * nb_base_functions, data));
1416  } else {
1417  // in some cases, f.e. for derivatives of nodal base functions only one
1418  // gauss point is needed
1419  return MatrixAdaptor(
1420  getN(base).size1(), getN(base).size2(),
1421  ublas::shallow_array_adaptor<double>(getDiffN(base).data().size(),
1422  &getDiffN(base).data()[0]));
1423  }
1424 }

◆ getDiffN() [5/7]

const MatrixAdaptor MoFEM::EntitiesFieldData::EntData::getDiffN ( const int  gg)
inline

get derivative of base functions at Gauss pts

returned matrix on rows has base functions, in column its derivatives.

Parameters
ggnb. of Gauss pts.

Definition at line 1381 of file EntitiesFieldData.hpp.

1381  {
1382  return getDiffN(bAse, gg);
1383 }

◆ getDiffN() [6/7]

const MatrixAdaptor MoFEM::EntitiesFieldData::EntData::getDiffN ( const int  gg,
const int  nb_base_functions 
)
inline

get derivatives of base functions at Gauss pts

Note that multi field element, two different field can have different approximation orders. Since we use hierarchical approximation basis, base functions are calculated once for element, using maximal approximation order on given entity.

Parameters
ggnb. of Gauss point
nb_base_functionsnumber of of base functions

Definition at line 1427 of file EntitiesFieldData.hpp.

1428  {
1429  return getDiffN(bAse, gg, nb_base_functions);
1430 }

◆ getDiffN() [7/7]

MatrixDouble & MoFEM::EntitiesFieldData::EntData::getDiffN ( const std::string &  field_name)
inline

get derivatives of base functions

Matrix at rows has nb. of Gauss pts, at columns it has derivative of base functions. Columns are structured as follows, [ dN1/dx, dN1/dy, dN1/dz, dN2/dx, dN2/dy, dN2/dz, ... ]

Scalar base functions: Note that base functions are calculated in file H1.c Above description not apply for derivatives of nodal functions, since derivative of nodal functions in case of simplexes, EDGES, TRIANGLES and TETS are constant. So that matrix rows represents nb. of base functions, columns are derivatives. Nb. of columns depend on element dimension, for EDGES is one, for TRIS is 2 and TETS is 3.

Note
Note that for node element this function make no sense.

Tonsorial base functions:

Note
Note: In rows ale integration pts, columns are formatted that that components of vectors and then derivatives, for example row for given integration points is formatted in array

\[ t_{0,0}, t_{1,0}, t_{1,0}, t_{0,1}, t_{1,1}, t_{1,1}, t_{0,2}, t_{1,2}, t_{1,2} \]

where comma express derivative, i.e. \(t_{2,1} = \frac{\partial t_2}{\partial \xi_1}\)

Definition at line 1335 of file EntitiesFieldData.hpp.

1335  {
1336  return *(getBBDiffNSharedPtr(field_name));
1337 }

◆ getDiffNSharedPtr()

boost::shared_ptr< MatrixDouble > & MoFEM::EntitiesFieldData::EntData::getDiffNSharedPtr ( const FieldApproximationBase  base)
virtual

Get shared pointer to derivatives of base base functions

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 37 of file EntitiesFieldData.cpp.

38  {
39  return diffN[base];
40 }

◆ getEntDataBitRefLevel()

std::vector< BitRefLevel > & MoFEM::EntitiesFieldData::EntData::getEntDataBitRefLevel ( )
virtual

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 944 of file EntitiesFieldData.cpp.

944  {
945  return entDataBitRefLevel;
946 }

◆ getFieldData() [1/2]

VectorDouble& MoFEM::EntitiesFieldData::EntData::getFieldData ( )
inline

get dofs data stature FEDofEntity

◆ getFieldData() [2/2]

VectorDouble & MoFEM::EntitiesFieldData::EntData::getFieldData ( ) const
inline

◆ getFieldDataUpToOrder()

const VectorAdaptor MoFEM::EntitiesFieldData::EntData::getFieldDataUpToOrder ( int  order)
inline

get dofs values up to given order

Definition at line 1246 of file EntitiesFieldData.hpp.

1246  {
1247  unsigned int size = 0;
1248  if (auto dof = dOfs[0]) {
1249  size = dof->getOrderNbDofs(order) * dof->getNbOfCoeffs();
1250  size = size < fieldData.size() ? size : fieldData.size();
1251  }
1252  double *data = &*fieldData.data().begin();
1253  return getVectorAdaptor(data, size);
1254 }

◆ getFieldDofs() [1/2]

VectorDofs& MoFEM::EntitiesFieldData::EntData::getFieldDofs ( )
inline

◆ getFieldDofs() [2/2]

VectorDofs & MoFEM::EntitiesFieldData::EntData::getFieldDofs ( ) const
inline

get dofs data stature FEDofEntity

Examples
EshelbianOperators.cpp, HookeElement.cpp, and prism_elements_from_surface.cpp.

Definition at line 1256 of file EntitiesFieldData.hpp.

1256  {
1257  return dOfs;
1258 }

◆ getFieldEntities() [1/2]

VectorFieldEntities& MoFEM::EntitiesFieldData::EntData::getFieldEntities ( )
inline

get field entities

◆ getFieldEntities() [2/2]

const VectorFieldEntities & MoFEM::EntitiesFieldData::EntData::getFieldEntities ( ) const
inline

get field entities

Examples
hanging_node_approx.cpp, and test_cache_on_entities.cpp.

Definition at line 1264 of file EntitiesFieldData.hpp.

1264  {
1265  return fieldEntities;
1266 }

◆ getFTensor0FieldData()

FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > MoFEM::EntitiesFieldData::EntData::getFTensor0FieldData ( )

Resturn scalar files as a FTensor of rank 0.

Returns
FTensor::Tensor0<FTensor::PackPtr<double *,1> >
Examples
continuity_check_on_skeleton_3d.cpp, and scalar_check_approximation.cpp.

Definition at line 506 of file EntitiesFieldData.cpp.

506  {
507 #ifndef NDEBUG
508  for (auto &d : dOfs) {
509  if (d) {
510  if (d->getNbOfCoeffs() != 1) {
511  std::stringstream s;
512  s << "Wrong number of coefficients is " << d->getNbOfCoeffs();
513  s << " but expected scalar field, tensor of rank 0";
515  }
516  break;
517  }
518  }
519 #endif
521  &*fieldData.data().begin());
522 }

◆ getFTensor0N() [1/4]

FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > MoFEM::EntitiesFieldData::EntData::getFTensor0N ( )
inline

Get base function as Tensor0.

Return base functions for field base

Returns
Tensor0

Definition at line 1495 of file EntitiesFieldData.hpp.

1495  {
1496  return getFTensor0N(bAse);
1497 };

◆ getFTensor0N() [2/4]

FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > MoFEM::EntitiesFieldData::EntData::getFTensor0N ( const FieldApproximationBase  base)
inline

◆ getFTensor0N() [3/4]

FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > MoFEM::EntitiesFieldData::EntData::getFTensor0N ( const FieldApproximationBase  base,
const int  gg,
const int  bb 
)
inline

Get base function as Tensor0 (Loop by integration points)

Parameters
base
ggintegration points
bbbase function
Returns
Tensor0

Note that:

t0 = data.getFTensor0N(base,bb);
++t0

Increment in above code will move pointer to base function in next integration point.

Definition at line 1500 of file EntitiesFieldData.hpp.

1501  {
1502  double *ptr = &getN(base)(gg, bb);
1504 };

◆ getFTensor0N() [4/4]

FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > MoFEM::EntitiesFieldData::EntData::getFTensor0N ( const int  gg,
const int  bb 
)
inline

Get base function as Tensor0 (Loop by integration points)

Return base functions for field base

Parameters
bbbase function
Returns
Tensor0

Definition at line 1507 of file EntitiesFieldData.hpp.

1507  {
1508  return getFTensor0N(bAse, gg, bb);
1509 };

◆ getFTensor1DiffN() [1/8]

template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN ( )

Get derivatives of base functions.

For volume element like tetrahedral or prism,

Tensor1<double*,3> diff_base = data.getFTensor1DiffN<3>();

For face element like triangle or quad

Tensor1<double*,2> diff_base = data.getFTensor1DiffN<2>();
Returns
Tensor rank 1 (vector)

Definition at line 536 of file EntitiesFieldData.cpp.

536  {
537  return getFTensor1DiffN<Tensor_Dim>(bAse);
538 }

◆ getFTensor1DiffN() [2/8]

FTensor::Tensor1< FTensor::PackPtr< double *, 2 >, 2 > MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN< 2 >

Get spatial derivative of base function tensor for dimension 3d.

Get spatial derivative of base function tensor for dimension 2d.

Definition at line 557 of file EntitiesFieldData.cpp.

557  {
558  return getFTensor1DiffN<3>(bAse);
559 }

◆ getFTensor1DiffN() [3/8]

template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN ( const FieldApproximationBase  base)

Get derivatives of base functions.

For volume element like tetrahedral or prism,

Tensor1<double*,3> diff_base = data.getFTensor1DiffN<3>();

For face element like triangle or quad

Tensor1<double*,2> diff_base = data.getFTensor1DiffN<2>();
Parameters
basefunctions
Returns
Tensor rank 1 (vector)
Examples
heat_method.cpp, HookeElement.cpp, HookeElement.hpp, HookeInternalStressElement.hpp, level_set.cpp, NavierStokesElement.cpp, reaction_diffusion.cpp, scalar_check_approximation.cpp, shallow_wave.cpp, and simple_elasticity.cpp.

Definition at line 526 of file EntitiesFieldData.cpp.

527  {
528  std::stringstream s;
529  s << "Template for tensor dimension " << Tensor_Dim << " not implemented";
530  THROW_MESSAGE(s.str());
532 }

◆ getFTensor1DiffN() [4/8]

FTensor::Tensor1< FTensor::PackPtr< double *, 2 >, 2 > MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN< 2 > ( const FieldApproximationBase  base)

Get spatial derivative of base function tensor for dimension 3d.

Get spatial derivative of base function tensor for dimension 2d.

Definition at line 545 of file EntitiesFieldData.cpp.

546  {
547  double *ptr = &*getDiffN(base).data().begin();
548  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(ptr, &ptr[1],
549  &ptr[2]);
550 }

◆ getFTensor1DiffN() [5/8]

FTensor::Tensor1< FTensor::PackPtr< double *, 2 >, 2 > MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN< 2 > ( const FieldApproximationBase  base,
const int  gg,
const int  bb 
)

Get spatial derivative of base function tensor for dimension 3d.

Get spatial derivative of base function tensor for dimension 2d.

Definition at line 596 of file EntitiesFieldData.cpp.

597  {
598  double *ptr = &getDiffN(base)(gg, 3 * bb);
599  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(ptr, &ptr[1],
600  &ptr[2]);
601 }

◆ getFTensor1DiffN() [6/8]

template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN ( const FieldApproximationBase  base,
const int  gg,
const int  bb 
)

Get derivatives of base functions (Loop by integration points)

For volume element like tetrahedral or prism,

Tensor1<double*,3> diff_base = data.getFTensor1DiffN<3>(base,gg,bb);

where bb is base function and gg is integration pt. Operator ++diff_base will move tensor pointer to next integration point.

For face element like triangle or quad

Tensor1<double*,2> diff_base = data.getFTensor1DiffN<2>(base,gg,bb);
Returns
Tensor rank 1 (vector)

Definition at line 583 of file EntitiesFieldData.cpp.

584  {
585  std::stringstream s;
586  s << "Template for tensor dimension " << Tensor_Dim << " not implemented";
587  THROW_MESSAGE(s.str());
589 }

◆ getFTensor1DiffN() [7/8]

FTensor::Tensor1< FTensor::PackPtr< double *, 2 >, 2 > MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN< 2 > ( const int  gg,
const int  bb 
)

Get spatial derivative of base function tensor for dimension 3d.

Get spatial derivative of base function tensor for dimension 2d.

Definition at line 608 of file EntitiesFieldData.cpp.

608  {
609  return getFTensor1DiffN<3>(bAse, gg, bb);
610 }

◆ getFTensor1DiffN() [8/8]

template<int Tensor_Dim>
FTensor::Tensor1<FTensor::PackPtr<double *, Tensor_Dim>, Tensor_Dim> MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN ( const int  gg,
const int  bb 
)

Get derivatives of base functions (Loop by integration points)

For volume element like tetrahedral or prism,

Tensor1<double*,3> diff_base = data.getFTensor1DiffN<3>(gg,bb);

where bb is base function and gg is integration pt. Operator ++diff_base will move tensor pointer to next base function.

For face element like triangle or quad

Tensor1<double*,2> diff_base = data.getFTensor1DiffN<2>(gg,bb);
Returns
Tensor rank 1 (vector)

◆ getFTensor1FieldData() [1/2]

template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > MoFEM::EntitiesFieldData::EntData::getFTensor1FieldData ( )

Return FTensor of rank 1, i.e. vector from filed data coeffinects.

auto t_vec = data.getFTensor1FieldData<3>();
Template Parameters
Tensor_Dimsize of vector
Returns
FTensor::Tensor1<FTensor::PackPtr<double *, Tensor_Dim>, Tensor_Dim>

Definition at line 1275 of file EntitiesFieldData.hpp.

1275  {
1276  std::stringstream s;
1277  s << "Not implemented for this dimension dim = " << Tensor_Dim;
1278  THROW_MESSAGE(s.str());
1279 }

◆ getFTensor1FieldData() [2/2]

FTensor::Tensor1< FTensor::PackPtr< double *, 1 >, 1 > MoFEM::EntitiesFieldData::EntData::getFTensor1FieldData< 1 >

Definition at line 287 of file EntitiesFieldData.cpp.

287  {
288 #ifndef NDEBUG
289  for (auto &d : dOfs) {
290  if (d) {
291  if (d->getNbOfCoeffs() != 3) {
292  std::stringstream s;
293  s << "Wrong number of coefficients is " << d->getNbOfCoeffs();
294  s << " but you ask for tensor rank 1 dimension 3";
296  }
297  break;
298  }
299  }
300 #endif
301  double *ptr = &*fieldData.data().begin();
302  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(ptr, &ptr[1],
303  &ptr[2]);
304 }

◆ getFTensor1N() [1/6]

template<int Tensor_Dim>
auto MoFEM::EntitiesFieldData::EntData::getFTensor1N ( )

Get base functions for Hdiv space.

Example:

int nb_dofs = data.getFieldData().size();
auto t_n_hdiv = data.getFTensor1N<3>();
for(int gg = 0;gg!=nb_gauss_pts;gg++) {
int ll = 0;
for(;ll!=nb_dofs;ll++) {
double dot_product = t_n_hdiv(i)*t_n_hdiv(i);
++t_n_hdiv;
}
for(;ll!=data.getVectorN().size2()/3;ll++) {
++t_n_hdiv;
}
}

Definition at line 1511 of file EntitiesFieldData.hpp.

1511  {
1512  return getFTensor1N<Tensor_Dim>(bAse);
1513 }

◆ getFTensor1N() [2/6]

template<int Tensor_Dim>
auto MoFEM::EntitiesFieldData::EntData::getFTensor1N ( const int  gg,
const int  bb 
)
inline

Get Hdiv base functions at integration point.

for(int gg = 0;gg!=nb_gauss_pts;gg++) {
auto t_base = data.getFTensor1N(gg,0);
for(int bb = 0;bb!=nb_base_functions;bb++) {
double dot = t_base(i)*t_base(i);
}
}

Definition at line 1516 of file EntitiesFieldData.hpp.

1516  {
1517  return getFTensor1N<Tensor_Dim>(bAse, gg, bb);
1518 }

◆ getFTensor1N() [3/6]

template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > MoFEM::EntitiesFieldData::EntData::getFTensor1N ( FieldApproximationBase  base)

Get base functions for Hdiv/Hcurl spaces.

Note
You probably like to use getFTensor1N(), in typical use base is set automatically based on base set to field.
Parameters
baseApproximation base

Example:

int nb_dofs = data.getFieldData().size();
auto t_n_hdiv = data.getFTensor1N<3>();
for(int gg = 0;gg!=nb_gauss_pts;gg++) {
int ll = 0;
for(;ll!=nb_dofs;ll++) {
double dot_product = t_n_hdiv(i)*t_n_hdiv(i);
++t_n_hdiv;
}
for(;ll!=data.getVectorN().size2()/3;ll++) {
++t_n_hdiv;
}
}
Examples
continuity_check_on_skeleton_3d.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, EshelbianOperators.cpp, hcurl_curl_operator.cpp, and hcurl_divergence_operator_2d.cpp.

Definition at line 640 of file EntitiesFieldData.cpp.

640  {
641  std::stringstream s;
642  s << "Template for tensor dimension " << Tensor_Dim << " not implemented";
643  THROW_MESSAGE(s.str());
645 }

◆ getFTensor1N() [4/6]

FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > MoFEM::EntitiesFieldData::EntData::getFTensor1N< 3 > ( FieldApproximationBase  base)

Definition at line 682 of file EntitiesFieldData.cpp.

682  {
683  double *t_n_ptr = &*getN(base).data().begin();
684  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(t_n_ptr, // HVEC0
685  &t_n_ptr[HVEC1],
686  &t_n_ptr[HVEC2]);
687 }

◆ getFTensor1N() [5/6]

template<int Tensor_Dim>
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > MoFEM::EntitiesFieldData::EntData::getFTensor1N ( FieldApproximationBase  base,
const int  gg,
const int  bb 
)

Get Hdiv base functions at integration point.

for(int gg = 0;gg!=nb_gauss_pts;gg++) {
auto t_base = data.getFTensor1N(base,gg,bb);
for(int bb = 0;bb!=nb_base_functions;bb++) {
auto dot = t_base(i)*t_base(i);
}
}

Definition at line 649 of file EntitiesFieldData.cpp.

650  {
651  std::stringstream s;
652  s << "Template for tensor dimension " << Tensor_Dim << " not implemented";
653  THROW_MESSAGE(s.str());
655 }

◆ getFTensor1N() [6/6]

FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > MoFEM::EntitiesFieldData::EntData::getFTensor1N< 3 > ( FieldApproximationBase  base,
const int  gg,
const int  bb 
)

Definition at line 691 of file EntitiesFieldData.cpp.

692  {
693  double *t_n_ptr = &getN(base)(gg, 3 * bb);
694  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(t_n_ptr, // HVEC0
695  &t_n_ptr[HVEC1],
696  &t_n_ptr[HVEC2]);
697 }

◆ getFTensor2DiffN() [1/6]

template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2<FTensor::PackPtr<double *, Tensor_Dim0 * Tensor_Dim1>, Tensor_Dim0, Tensor_Dim1> MoFEM::EntitiesFieldData::EntData::getFTensor2DiffN ( )
inline

Get derivatives of base functions for Hdiv space.

Definition at line 751 of file EntitiesFieldData.hpp.

751  {
752  return getFTensor2DiffN<Tensor_Dim0, Tensor_Dim1>(bAse);
753  }

◆ getFTensor2DiffN() [2/6]

template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2<FTensor::PackPtr<double *, Tensor_Dim0 * Tensor_Dim1>, Tensor_Dim0, Tensor_Dim1> MoFEM::EntitiesFieldData::EntData::getFTensor2DiffN ( const int  gg,
const int  bb 
)
inline

Get derivatives of base functions for Hdiv space at integration pts.

Definition at line 761 of file EntitiesFieldData.hpp.

761  {
762  return getFTensor2DiffN<Tensor_Dim0, Tensor_Dim1>(bAse, gg, bb);
763  }

◆ getFTensor2DiffN() [3/6]

template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > MoFEM::EntitiesFieldData::EntData::getFTensor2DiffN ( FieldApproximationBase  base)

Get derivatives of base functions for Hdiv space.

Examples
EshelbianOperators.cpp, hcurl_curl_operator.cpp, hcurl_divergence_operator_2d.cpp, and hdiv_divergence_operator.cpp.

Definition at line 660 of file EntitiesFieldData.cpp.

660  {
661  std::stringstream s;
662  s << "Template for tensor dimension " << Tensor_Dim0 << "x" << Tensor_Dim1
663  << " not implemented";
664  THROW_MESSAGE(s.str());
666 }

◆ getFTensor2DiffN() [4/6]

FTensor::Tensor2< FTensor::PackPtr< double *, 6 >, 3, 2 > MoFEM::EntitiesFieldData::EntData::getFTensor2DiffN< 3, 2 > ( FieldApproximationBase  base)

Definition at line 701 of file EntitiesFieldData.cpp.

702  {
703  double *t_diff_n_ptr = &*getDiffN(base).data().begin();
705  t_diff_n_ptr, &t_diff_n_ptr[HVEC0_1], &t_diff_n_ptr[HVEC0_2],
706  &t_diff_n_ptr[HVEC1_0], &t_diff_n_ptr[HVEC1_1], &t_diff_n_ptr[HVEC1_2],
707  &t_diff_n_ptr[HVEC2_0], &t_diff_n_ptr[HVEC2_1], &t_diff_n_ptr[HVEC2_2]);
708 }

◆ getFTensor2DiffN() [5/6]

template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > MoFEM::EntitiesFieldData::EntData::getFTensor2DiffN ( FieldApproximationBase  base,
const int  gg,
const int  bb 
)

Get derivatives of base functions for Hdiv space at integration pts.

Definition at line 671 of file EntitiesFieldData.cpp.

672  {
673  std::stringstream s;
674  s << "Template for tensor dimension " << Tensor_Dim0 << "x" << Tensor_Dim1
675  << " not implemented";
676  THROW_MESSAGE(s.str());
678 }

◆ getFTensor2DiffN() [6/6]

FTensor::Tensor2< FTensor::PackPtr< double *, 6 >, 3, 2 > MoFEM::EntitiesFieldData::EntData::getFTensor2DiffN< 3, 2 > ( FieldApproximationBase  base,
const int  gg,
const int  bb 
)

Definition at line 712 of file EntitiesFieldData.cpp.

713  {
714  double *t_diff_n_ptr = &getDiffN(base)(gg, 9 * bb);
716  t_diff_n_ptr, &t_diff_n_ptr[HVEC0_1], &t_diff_n_ptr[HVEC0_2],
717  &t_diff_n_ptr[HVEC1_0], &t_diff_n_ptr[HVEC1_1], &t_diff_n_ptr[HVEC1_2],
718  &t_diff_n_ptr[HVEC2_0], &t_diff_n_ptr[HVEC2_1], &t_diff_n_ptr[HVEC2_2]);
719 }

◆ getFTensor2FieldData() [1/2]

template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > MoFEM::EntitiesFieldData::EntData::getFTensor2FieldData ( )

Return FTensor rank 2, i.e. matrix from filed data coeffinects.

auto t_mat = data.getFTensor2FieldData<3,3>();
Template Parameters
Tensor_Dim0
Tensor_Dim1
Returns
FTensor::Tensor2<FTensor::PackPtr<double *, Tensor_Dim0 * Tensor_Dim1>, Tensor_Dim0, Tensor_Dim1>

Definition at line 1284 of file EntitiesFieldData.hpp.

1284  {
1285  std::stringstream s;
1286  s << "Not implemented for this dimension dim0 = " << Tensor_Dim0;
1287  s << " and dim1 " << Tensor_Dim1;
1288  THROW_MESSAGE(s.str());
1289 }

◆ getFTensor2FieldData() [2/2]

FTensor::Tensor2< FTensor::PackPtr< double *, 9 >, 3, 3 > MoFEM::EntitiesFieldData::EntData::getFTensor2FieldData< 3, 3 >

Definition at line 348 of file EntitiesFieldData.cpp.

348  {
349 #ifndef NDEBUG
350  for (auto &d : dOfs) {
351  if (d) {
352  if (d->getNbOfCoeffs() != 1) {
353  std::stringstream s;
354  s << "Wrong number of coefficients is " << d->getNbOfCoeffs();
355  s << " but you ask for tensor rank 2 dimensions 1 by 1 so 1 "
356  "coefficients "
357  "is expected";
359  }
360  break;
361  }
362  }
363 #endif
364  double *ptr = &*fieldData.data().begin();
366 }

◆ getFTensor2N() [1/6]

template<int Tensor_Dim0, int Tensor_Dim1>
auto MoFEM::EntitiesFieldData::EntData::getFTensor2N ( )

Get base functions for Hdiv space.

Example:

int nb_dofs = data.getFieldData().size();
auto t_n_hdiv = data.getFTensor2N<3,3>();
for(int gg = 0;gg!=nb_gauss_pts;gg++) {
int ll = 0;
for(;ll!=nb_dofs;ll++) {
double dot_product = t_n_hdiv(i,j)*t_n_hdiv(i,j);
++t_n_hdiv;
}
for(;ll!=data.getVectorN().size2()/3;ll++) {
++t_n_hdiv;
}
}

Definition at line 1521 of file EntitiesFieldData.hpp.

1521  {
1522  return getFTensor2N<Tensor_Dim0, Tensor_Dim1>(bAse);
1523 }

◆ getFTensor2N() [2/6]

template<int Tensor_Dim0, int Tensor_Dim1>
auto MoFEM::EntitiesFieldData::EntData::getFTensor2N ( const int  gg,
const int  bb 
)

Get base functions for Hdiv space.

Example:

int nb_dofs = data.getFieldData().size();
for(int gg = 0;gg!=nb_gauss_pts;++gg) {
int ll = 0;
auto t_n_hdiv = data.getFTensor2N<3,3>(gg,0);
for(;ll!=nb_dofs;ll++) {
double dot_product = t_n_hdiv(i)*t_n_hdiv(i);
++t_n_hdiv;
}
for(;ll!=data.getVectorN().size2()/3;ll++) {
++t_n_hdiv;
}
}

Definition at line 1526 of file EntitiesFieldData.hpp.

1526  {
1527  return getFTensor2N<Tensor_Dim0, Tensor_Dim1>(bAse, gg, bb);
1528 }

◆ getFTensor2N() [3/6]

template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > MoFEM::EntitiesFieldData::EntData::getFTensor2N ( FieldApproximationBase  base)

Get base functions for Hdiv/Hcurl spaces.

Note
You probably like to use getFTensor1N(), in typical use base is set automatically based on base set to field.
Parameters
baseApproximation base

Example:

int nb_dofs = data.getFieldData().size();
auto t_n_hdiv = data.getFTensor2N<3,3>();
for(int gg = 0;gg!=nb_gauss_pts;gg++) {
int ll = 0;
for(;ll!=nb_dofs;ll++) {
double dot_product = t_n_hdiv(i,j)*t_n_hdiv(i,j);
++t_n_hdiv;
}
for(;ll!=data.getVectorN().size2()/3;ll++) {
++t_n_hdiv;
}
}
Examples
EshelbianOperators.cpp.

Definition at line 762 of file EntitiesFieldData.cpp.

762  {
763  std::stringstream s;
764  s << "Template for tensor dimension " << Tensor_Dim0 << ", " << Tensor_Dim1
765  << " not implemented";
766  THROW_MESSAGE(s.str());
768  Tensor_Dim0, Tensor_Dim1>();
769 }

◆ getFTensor2N() [4/6]

FTensor::Tensor2< FTensor::PackPtr< double *, 9 >, 3, 3 > MoFEM::EntitiesFieldData::EntData::getFTensor2N< 3, 3 > ( FieldApproximationBase  base)

Definition at line 773 of file EntitiesFieldData.cpp.

773  {
774  double *t_n_ptr = &*(getN(base).data().begin());
776 
777  &t_n_ptr[HVEC0], &t_n_ptr[HVEC1], &t_n_ptr[HVEC2],
778 
779  &t_n_ptr[3 + HVEC0], &t_n_ptr[3 + HVEC1], &t_n_ptr[3 + HVEC2],
780 
781  &t_n_ptr[6 + HVEC0], &t_n_ptr[6 + HVEC1], &t_n_ptr[6 + HVEC2]
782 
783  );
784 }

◆ getFTensor2N() [5/6]

template<>
FTensor::Tensor2<FTensor::PackPtr<double *, 9>, 3, 3> MoFEM::EntitiesFieldData::EntData::getFTensor2N ( FieldApproximationBase  base,
const int  gg,
const int  bb 
)

Definition at line 801 of file EntitiesFieldData.cpp.

802  {
803  double *t_n_ptr = &getN(base)(gg, 9 * bb);
805 
806  &t_n_ptr[HVEC0], &t_n_ptr[HVEC1], &t_n_ptr[HVEC2],
807 
808  &t_n_ptr[3 + HVEC0], &t_n_ptr[3 + HVEC1], &t_n_ptr[3 + HVEC2],
809 
810  &t_n_ptr[6 + HVEC0], &t_n_ptr[6 + HVEC1], &t_n_ptr[6 + HVEC2]
811 
812  );
813 }

◆ getFTensor2N() [6/6]

template<int Tensor_Dim0, int Tensor_Dim1>
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > MoFEM::EntitiesFieldData::EntData::getFTensor2N ( FieldApproximationBase  base,
const int  gg,
const int  bb 
)

Get base functions for tensor Hdiv/Hcurl spaces.

Note
You probably like to use getFTensor2N(), in typical use base is set automatically based on base set to field.
Parameters
baseApproximation base

Example:

int nb_dofs = data.getFieldData().size();
for(int gg = 0;gg!=nb_gauss_pts;gg++) {
auto t_n_hdiv = data.getFTensor2N<3>(base,gg,bb);
int ll = 0;
for(;ll!=nb_dofs;ll++) {
double dot_product = t_n_hdiv(i,j)*t_n_hdiv(i,j);
++t_n_hdiv;
}
for(;ll!=data.getVectorN().size2()/3;ll++) {
++t_n_hdiv;
}
}

Definition at line 789 of file EntitiesFieldData.cpp.

790  {
791  std::stringstream s;
792  s << "Template for tensor dimension " << Tensor_Dim0 << ", " << Tensor_Dim1
793  << " not implemented";
794  THROW_MESSAGE(s.str());
796  Tensor_Dim0, Tensor_Dim1>();
797 }

◆ getFTensor2SymmetricFieldData() [1/2]

template<int Tensor_Dim>
FTensor::Tensor2_symmetric< FTensor::PackPtr< double *,(Tensor_Dim *(Tensor_Dim+1))/2 >, Tensor_Dim > MoFEM::EntitiesFieldData::EntData::getFTensor2SymmetricFieldData ( )

Return symmetric FTensor rank 2, i.e. matrix from filed data coeffinects.

auto t_mat = data.getFTensor2SymmetricFieldData<3>();
Template Parameters
Tensor_Dimdimension of the tensor
Returns
FTensor::Tensor2_symmetric<FTensor::PackPtr<double *, (Tensor_Dim
  • (Tensor_Dim + 1)) / 2>, Tensor_Dim>

Definition at line 1294 of file EntitiesFieldData.hpp.

1294  {
1295  std::stringstream s;
1296  s << "Not implemented for this dimension dim = " << Tensor_Dim;
1297  THROW_MESSAGE(s.str());
1298 }

◆ getFTensor2SymmetricFieldData() [2/2]

FTensor::Tensor2_symmetric< FTensor::PackPtr< double *, 3 >, 2 > MoFEM::EntitiesFieldData::EntData::getFTensor2SymmetricFieldData< 2 >

Definition at line 461 of file EntitiesFieldData.cpp.

461  {
462 #ifndef NDEBUG
463  for (auto &d : dOfs) {
464  if (d) {
465  if (d->getNbOfCoeffs() != 6) {
466  std::stringstream s;
467  s << "Wrong number of coefficients is " << d->getNbOfCoeffs();
468  s << " but you ask for symmetric tensor rank 2 dimensions 3 by 3 so 6 "
469  "coefficients "
470  "is expected";
472  }
473  break;
474  }
475  }
476 #endif
477  double *ptr = &*fieldData.data().begin();
479  ptr, &ptr[1], &ptr[2], &ptr[3], &ptr[4], &ptr[5]);
480 }

◆ getFTensor3Diff2N() [1/3]

template<int Tensor_Dim0, int Tensor_Dim1, int Tensor_Dim2>
FTensor::Tensor3< FTensor::PackPtr<double *, Tensor_Dim0 * Tensor_Dim1 * Tensor_Dim2>, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2> MoFEM::EntitiesFieldData::EntData::getFTensor3Diff2N ( )
inline

Get second derivatives of base functions for Hvec space.

Definition at line 779 of file EntitiesFieldData.hpp.

779  {
780  return getFTensor3Diff2N<Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>(bAse);
781  }

◆ getFTensor3Diff2N() [2/3]

template<int Tensor_Dim0, int Tensor_Dim1, int Tensor_Dim2>
FTensor::Tensor3< FTensor::PackPtr<double *, Tensor_Dim0 * Tensor_Dim1 * Tensor_Dim2>, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2> MoFEM::EntitiesFieldData::EntData::getFTensor3Diff2N ( FieldApproximationBase  base)

Get second derivatives of base functions for Hvec space.

◆ getFTensor3Diff2N() [3/3]

FTensor::Tensor3< FTensor::PackPtr< double *, 12 >, 3, 2, 2 > MoFEM::EntitiesFieldData::EntData::getFTensor3Diff2N ( FieldApproximationBase  base)

Definition at line 743 of file EntitiesFieldData.cpp.

743  {
744  double *ptr = &(getN(base, BaseDerivatives::SecondDerivative))(0, 0);
746 
747  &ptr[2 * HVEC0_0 + 0], &ptr[2 * HVEC0_0 + 1], &ptr[2 * HVEC0_1 + 0],
748  &ptr[2 * HVEC0_1 + 1],
749 
750  &ptr[2 * HVEC1_0 + 0], &ptr[2 * HVEC1_0 + 1], &ptr[2 * HVEC1_1 + 0],
751  &ptr[2 * HVEC1_1 + 1],
752 
753  &ptr[2 * HVEC2_0 + 0], &ptr[2 * HVEC2_0 + 1], &ptr[2 * HVEC2_1 + 0],
754  &ptr[2 * HVEC2_1 + 1]
755 
756  };
757 }

◆ getIndices() [1/2]

VectorInt& MoFEM::EntitiesFieldData::EntData::getIndices ( )
inline

◆ getIndices() [2/2]

VectorInt & MoFEM::EntitiesFieldData::EntData::getIndices ( ) const
inline

◆ getIndicesUpToOrder()

const VectorIntAdaptor MoFEM::EntitiesFieldData::EntData::getIndicesUpToOrder ( int  order)
inline

get global indices of dofs on entity up to given order

Definition at line 1206 of file EntitiesFieldData.hpp.

1206  {
1207  unsigned int size = 0;
1208  if (auto dof = dOfs[0]) {
1209  size = dof->getOrderNbDofs(order) * dof->getNbOfCoeffs();
1210  size = size < iNdices.size() ? size : iNdices.size();
1211  }
1212  int *data = &*iNdices.data().begin();
1213  return VectorIntAdaptor(size, ublas::shallow_array_adaptor<int>(size, data));
1214 }

◆ getLocalIndices() [1/2]

VectorInt& MoFEM::EntitiesFieldData::EntData::getLocalIndices ( )
inline

◆ getLocalIndices() [2/2]

VectorInt & MoFEM::EntitiesFieldData::EntData::getLocalIndices ( ) const
inline

get local indices of dofs on entity

Examples
boundary_marker.cpp, and prism_elements_from_surface.cpp.

Definition at line 1216 of file EntitiesFieldData.hpp.

1216  {
1217  return localIndices;
1218 }

◆ getLocalIndicesUpToOrder()

const VectorIntAdaptor MoFEM::EntitiesFieldData::EntData::getLocalIndicesUpToOrder ( int  order)
inline

get local indices of dofs on entity up to given order

Definition at line 1221 of file EntitiesFieldData.hpp.

1221  {
1222  unsigned int size = 0;
1223  if (auto dof = dOfs[0]) {
1224  size = dof->getOrderNbDofs(order) * dof->getNbOfCoeffs();
1225  size = size < localIndices.size() ? size : localIndices.size();
1226  }
1227  int *data = &*localIndices.data().begin();
1228  return VectorIntAdaptor(size, ublas::shallow_array_adaptor<int>(size, data));
1229 }

◆ getN() [1/9]

MatrixDouble & MoFEM::EntitiesFieldData::EntData::getN ( )
inline

get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb. of columns is equal to number of base functions on this entity.

Note
Note that for vectorial base, like Hdiv or Hcurl, in columns are vectorial base functions. For tonsorial would be tonsorial base functions. Interpretation depends on type of base, scalar, vectorial or tonsorial and dimension fo problem.

Definition at line 1313 of file EntitiesFieldData.hpp.

1313 { return getN(bAse); }

◆ getN() [2/9]

MatrixDouble & MoFEM::EntitiesFieldData::EntData::getN ( const BaseDerivatives  derivative)
inline

Get base function derivative.

Parameters
derivative
Returns
MatrixDouble&

Definition at line 1342 of file EntitiesFieldData.hpp.

1342  {
1343  return getN(bAse, derivative);
1344 }

◆ getN() [3/9]

MatrixDouble & MoFEM::EntitiesFieldData::EntData::getN ( const FieldApproximationBase  base)
inline

◆ getN() [4/9]

MatrixDouble & MoFEM::EntitiesFieldData::EntData::getN ( const FieldApproximationBase  base,
const BaseDerivatives  derivative 
)
inline

Get base function derivative.

Parameters
basebase
derivativederivative
Returns
MatrixDouble&

Definition at line 1321 of file EntitiesFieldData.hpp.

1322  {
1323 #ifndef NDEBUG
1324  if (!getNSharedPtr(base, derivative)) {
1325  MOFEM_LOG_C("SELF", Sev::error,
1326  "Ptr to base %s functions derivative %d is null",
1327  ApproximationBaseNames[base], derivative);
1328  THROW_MESSAGE("Null pointer");
1329  }
1330 #endif
1331  return *(getNSharedPtr(base, derivative));
1332 }

◆ getN() [5/9]

const VectorAdaptor MoFEM::EntitiesFieldData::EntData::getN ( const FieldApproximationBase  base,
const int  gg 
)
inline

get base functions at Gauss pts

Definition at line 1347 of file EntitiesFieldData.hpp.

1348  {
1349  int size = getN(base).size2();
1350  double *data = &getN(base)(gg, 0);
1351  return VectorAdaptor(size, ublas::shallow_array_adaptor<double>(size, data));
1352 }

◆ getN() [6/9]

const VectorAdaptor MoFEM::EntitiesFieldData::EntData::getN ( const FieldApproximationBase  base,
const int  gg,
const int  nb_base_functions 
)
inline

get base functions at Gauss pts

Note that multi field element, two different field can have different approximation orders. Since we use hierarchical approximation basis, base functions are calculated once for element, using maximal approximation order on given entity.

Parameters
baseApproximation base
ggnumber of Gauss point
nb_base_functionsnumber of of base functions returned

Definition at line 1386 of file EntitiesFieldData.hpp.

1387  {
1388  (void)getN()(gg, nb_base_functions -
1389  1); // throw error if nb_base_functions is to big
1390  double *data = &getN(base)(gg, 0);
1391  return VectorAdaptor(nb_base_functions, ublas::shallow_array_adaptor<double>(
1392  nb_base_functions, data));
1393 }

◆ getN() [7/9]

const VectorAdaptor MoFEM::EntitiesFieldData::EntData::getN ( const int  gg)
inline

get base functions at Gauss pts

Definition at line 1354 of file EntitiesFieldData.hpp.

1354  {
1355  return getN(bAse, gg);
1356 }

◆ getN() [8/9]

const VectorAdaptor MoFEM::EntitiesFieldData::EntData::getN ( const int  gg,
const int  nb_base_functions 
)
inline

get base functions at Gauss pts

Note that multi field element, two different field can have different approximation orders. Since we use hierarchical approximation basis, base functions are calculated once for element, using maximal approximation order on given entity.

Parameters
ggnumber of Gauss point
nb_base_functionsnumber of of base functions returned

Definition at line 1396 of file EntitiesFieldData.hpp.

1396  {
1397  return getN(bAse, gg, nb_base_functions);
1398 }

◆ getN() [9/9]

MatrixDouble & MoFEM::EntitiesFieldData::EntData::getN ( const std::string &  field_name)
inline

get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb. of columns is equal to number of base functions on this entity.

Note
Note that for vectorial base, like Hdiv or Hcurl, in columns are vectorial base functions. For tonsorial would be tonsorial base functions. Interpretation depends on type of base, scalar, vectorial or tonsorial and dimension fo problem.

Definition at line 1309 of file EntitiesFieldData.hpp.

1309  {
1310  return *(getBBNSharedPtr(field_name));
1311 }

◆ getNSharedPtr() [1/2]

boost::shared_ptr< MatrixDouble > & MoFEM::EntitiesFieldData::EntData::getNSharedPtr ( const FieldApproximationBase  base)
virtual

Get shared pointer to base base functions

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 33 of file EntitiesFieldData.cpp.

33  {
34  return N[base];
35 }

◆ getNSharedPtr() [2/2]

boost::shared_ptr< MatrixDouble > & MoFEM::EntitiesFieldData::EntData::getNSharedPtr ( const FieldApproximationBase  base,
const BaseDerivatives  direvatie 
)
virtual

Get shared pointer to base base functions

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 27 of file EntitiesFieldData.cpp.

28  {
29  return baseFunctionsAndBaseDerivatives[direvatie][base];
30 }

◆ getOrder() [1/2]

ApproximationOrder& MoFEM::EntitiesFieldData::EntData::getOrder ( )
inline

◆ getOrder() [2/2]

ApproximationOrder & MoFEM::EntitiesFieldData::EntData::getOrder ( ) const
inline

get approximation order

Definition at line 1197 of file EntitiesFieldData.hpp.

1197  {
1198  return oRder;
1199 }

◆ getSense() [1/2]

int& MoFEM::EntitiesFieldData::EntData::getSense ( )
inline

◆ getSense() [2/2]

int MoFEM::EntitiesFieldData::EntData::getSense ( ) const
virtual

get entity sense, need to calculate base functions with conforming approximation fields

Reimplemented in MoFEM::DerivedEntitiesFieldData::DerivedEntData.

Definition at line 1231 of file EntitiesFieldData.hpp.

1231 { return sEnse; }

◆ getSpace()

FieldSpace & MoFEM::EntitiesFieldData::EntData::getSpace ( )
inline

Get field space.

Returns
Field space
Examples
continuity_check_on_contact_prism_side_ele.cpp, continuity_check_on_skeleton_3d.cpp, and hanging_node_approx.cpp.

Definition at line 1302 of file EntitiesFieldData.hpp.

1302 { return sPace; }

◆ getVectorDiffN() [1/4]

template<int DIM0, int DIM1>
const MatrixAdaptor MoFEM::EntitiesFieldData::EntData::getVectorDiffN ( const FieldApproximationBase  base,
const int  dof,
const int  gg 
)
inline

get DiffHdiv of base functions at Gauss pts

Parameters
baseApproximation base
ggnb. of Gauss point
numberof of base functions

Definition at line 1474 of file EntitiesFieldData.hpp.

1475  {
1476  double *data =
1477  &EntitiesFieldData::EntData::getDiffN(base)(gg, DIM0 * DIM1 * dof);
1478  return MatrixAdaptor(DIM0, DIM1,
1479  ublas::shallow_array_adaptor<double>(DIM0 * DIM1, data));
1480 }

◆ getVectorDiffN() [2/4]

template<int DIM0, int DIM1>
const MatrixAdaptor MoFEM::EntitiesFieldData::EntData::getVectorDiffN ( const int  dof,
const int  gg 
)
inline

get DiffHdiv of base functions at Gauss pts

Parameters
ggnb. of Gauss point
numberof of base functions

Definition at line 1483 of file EntitiesFieldData.hpp.

1484  {
1485  return getVectorDiffN<DIM0, DIM1>(bAse, dof, gg);
1486 }

◆ getVectorDiffN() [3/4]

template<int DIM0, int DIM1>
const MatrixAdaptor MoFEM::EntitiesFieldData::EntData::getVectorDiffN ( const int  gg)
inline

get DiffHdiv of base functions at Gauss pts

Parameters
ggnb. of Gauss point
numberof of base functions

Definition at line 1468 of file EntitiesFieldData.hpp.

1468  {
1469  return getVectorDiffN<DIM0, DIM1>(bAse, gg);
1470 }

◆ getVectorDiffN() [4/4]

template<int DIM0, int DIM1>
const MatrixAdaptor MoFEM::EntitiesFieldData::EntData::getVectorDiffN ( FieldApproximationBase  base,
const int  gg 
)
inline

get DiffHdiv of base functions at Gauss pts

Parameters
baseApproximation base
ggnb. of Gauss point
numberof of base functions

Definition at line 1454 of file EntitiesFieldData.hpp.

1455  {
1456  if (PetscUnlikely(getDiffN(base).size2() % (DIM0 * DIM1))) {
1457  THROW_MESSAGE("Wrong dimension");
1458  }
1459 
1460  const int nb_base_functions = getN(base).size2() / (DIM0 * DIM1);
1461  double *data = &getN(base)(gg, 0);
1462  return MatrixAdaptor(nb_base_functions, DIM0 * DIM1,
1463  ublas::shallow_array_adaptor<double>(
1464  DIM0 * DIM1 * nb_base_functions, data));
1465 }

◆ getVectorN() [1/2]

template<int DIM>
const MatrixAdaptor MoFEM::EntitiesFieldData::EntData::getVectorN ( const FieldApproximationBase  base,
const int  gg 
)
inline

get Hdiv of base functions at Gauss pts

Parameters
baseApproximation base
ggnb. of Gauss point
Examples
hdiv_divergence_operator.cpp.

Definition at line 1434 of file EntitiesFieldData.hpp.

1435  {
1436  if (PetscUnlikely(getN(base).size2() % DIM)) {
1437  THROW_MESSAGE("Wrong dimension");
1438  }
1439 
1440  const int nb_base_functions = getN(base).size2() / DIM;
1441  double *data = &getN(base)(gg, 0);
1442  return MatrixAdaptor(
1443  nb_base_functions, DIM,
1444  ublas::shallow_array_adaptor<double>(DIM * nb_base_functions, data));
1445 }

◆ getVectorN() [2/2]

template<int DIM>
const MatrixAdaptor MoFEM::EntitiesFieldData::EntData::getVectorN ( const int  gg)
inline

get Hdiv of base functions at Gauss pts

Parameters
ggnb. of Gauss point
numberof of base functions

Definition at line 1448 of file EntitiesFieldData.hpp.

1448  {
1449  return getVectorN<DIM>(bAse, gg);
1450 }

◆ resetFieldDependentData()

MoFEMErrorCode MoFEM::EntitiesFieldData::EntData::resetFieldDependentData ( )
inline

Reset data associated with particular field name

Returns
error code

Definition at line 115 of file EntitiesFieldData.cpp.

115  {
117  sPace = NOSPACE;
118  bAse = NOBASE;
119  fieldEntities.resize(0, false);
120  iNdices.resize(0, false);
121  localIndices.resize(0, false);
122  dOfs.resize(0, false);
123  fieldData.resize(0, false);
125 }

Friends And Related Function Documentation

◆ OpAddParentEntData

friend struct OpAddParentEntData
friend

Definition at line 1114 of file EntitiesFieldData.hpp.

◆ operator<<

std::ostream& operator<< ( std::ostream &  os,
const EntitiesFieldData::EntData e 
)
friend

Definition at line 240 of file EntitiesFieldData.cpp.

241  {
242  os << "sEnse: " << e.getSense() << std::endl
243  << "oRder: " << e.getOrder() << std::endl
244  << "global indices: " << e.getIndices() << std::endl
245  << "local indices: " << e.getLocalIndices() << std::endl;
246  // FIXME: precision should not be set here
247  os << "fieldData: " << std::fixed << std::setprecision(2) << e.getFieldData()
248  << std::endl;
249  MatrixDouble base = const_cast<EntitiesFieldData::EntData &>(e).getN();
250  MatrixDouble diff_base =
251  const_cast<EntitiesFieldData::EntData &>(e).getDiffN();
252  const double eps = 1e-6;
253  for (unsigned int ii = 0; ii != base.size1(); ii++) {
254  for (unsigned int jj = 0; jj != base.size2(); jj++) {
255  if (fabs(base(ii, jj)) < eps)
256  base(ii, jj) = 0;
257  }
258  }
259  for (unsigned int ii = 0; ii != diff_base.size1(); ii++) {
260  for (unsigned int jj = 0; jj != diff_base.size2(); jj++) {
261  if (fabs(diff_base(ii, jj)) < eps)
262  diff_base(ii, jj) = 0;
263  }
264  }
265  os << "N: " << std::fixed << base << std::endl
266  << "diffN: " << std::fixed << diff_base;
267  return os;
268 }

Member Data Documentation

◆ bAse

FieldApproximationBase MoFEM::EntitiesFieldData::EntData::bAse
protected

Field approximation base.

Definition at line 1069 of file EntitiesFieldData.hpp.

◆ baseFunctionsAndBaseDerivatives

std::array<std::array<boost::shared_ptr<MatrixDouble>, LASTBASE>, LastDerivative> MoFEM::EntitiesFieldData::EntData::baseFunctionsAndBaseDerivatives
protected

Definition at line 1079 of file EntitiesFieldData.hpp.

◆ bbAlphaIndices

std::map<std::string, boost::shared_ptr<MatrixInt> > MoFEM::EntitiesFieldData::EntData::bbAlphaIndices
protected

Indices for Bernstein-Bezier (BB) base.

Definition at line 1090 of file EntitiesFieldData.hpp.

◆ bbAlphaIndicesByOrder

std::array<boost::shared_ptr<MatrixInt>, MaxBernsteinBezierOrder> MoFEM::EntitiesFieldData::EntData::bbAlphaIndicesByOrder
protected

BB alpha indices by order.

Definition at line 1097 of file EntitiesFieldData.hpp.

◆ bbDiffN

std::map<std::string, boost::shared_ptr<MatrixDouble> > MoFEM::EntitiesFieldData::EntData::bbDiffN
protected

Definition at line 1088 of file EntitiesFieldData.hpp.

◆ bbDiffNByOrder

std::array<boost::shared_ptr<MatrixDouble>, MaxBernsteinBezierOrder> MoFEM::EntitiesFieldData::EntData::bbDiffNByOrder
protected

BB base functions direvatives by order.

Definition at line 1095 of file EntitiesFieldData.hpp.

◆ bbFieldName

std::string MoFEM::EntitiesFieldData::EntData::bbFieldName
protected

field name

Definition at line 1085 of file EntitiesFieldData.hpp.

◆ bbN

std::map<std::string, boost::shared_ptr<MatrixDouble> > MoFEM::EntitiesFieldData::EntData::bbN
protected

Definition at line 1087 of file EntitiesFieldData.hpp.

◆ bbNByOrder

std::array<boost::shared_ptr<MatrixDouble>, MaxBernsteinBezierOrder> MoFEM::EntitiesFieldData::EntData::bbNByOrder
protected

BB base functions by order.

Definition at line 1093 of file EntitiesFieldData.hpp.

◆ bbNodeOrder

VectorInt MoFEM::EntitiesFieldData::EntData::bbNodeOrder
protected

order of nodes

Definition at line 1086 of file EntitiesFieldData.hpp.

◆ diffN

std::array<boost::shared_ptr<MatrixDouble>, LASTBASE>& MoFEM::EntitiesFieldData::EntData::diffN
protected

Derivatives of base functions.

Definition at line 1083 of file EntitiesFieldData.hpp.

◆ dOfs

VectorDofs MoFEM::EntitiesFieldData::EntData::dOfs
protected

DoFs on entity.

Definition at line 1072 of file EntitiesFieldData.hpp.

◆ entDataBitRefLevel

std::vector<BitRefLevel> MoFEM::EntitiesFieldData::EntData::entDataBitRefLevel
protected

Bit ref level in entity.

Definition at line 1075 of file EntitiesFieldData.hpp.

◆ fieldData

VectorDouble MoFEM::EntitiesFieldData::EntData::fieldData
protected

Field data on entity.

Definition at line 1074 of file EntitiesFieldData.hpp.

◆ fieldEntities

VectorFieldEntities MoFEM::EntitiesFieldData::EntData::fieldEntities
protected

Field entities.

Definition at line 1073 of file EntitiesFieldData.hpp.

◆ iNdices

VectorInt MoFEM::EntitiesFieldData::EntData::iNdices
protected

Global indices on entity.

Definition at line 1070 of file EntitiesFieldData.hpp.

◆ localIndices

VectorInt MoFEM::EntitiesFieldData::EntData::localIndices
protected

Local indices on entity.

Definition at line 1071 of file EntitiesFieldData.hpp.

◆ MaxBernsteinBezierOrder

constexpr size_t MoFEM::EntitiesFieldData::EntData::MaxBernsteinBezierOrder = BITFEID_SIZE
staticconstexpr

Definition at line 1036 of file EntitiesFieldData.hpp.

◆ N

std::array<boost::shared_ptr<MatrixDouble>, LASTBASE>& MoFEM::EntitiesFieldData::EntData::N
protected

Base functions.

Definition at line 1081 of file EntitiesFieldData.hpp.

◆ oRder

ApproximationOrder MoFEM::EntitiesFieldData::EntData::oRder
protected

Entity order.

Definition at line 1067 of file EntitiesFieldData.hpp.

◆ sEnse

int MoFEM::EntitiesFieldData::EntData::sEnse
protected

Entity sense (orientation)

Definition at line 1066 of file EntitiesFieldData.hpp.

◆ sPace

FieldSpace MoFEM::EntitiesFieldData::EntData::sPace
protected

Entity space.

Definition at line 1068 of file EntitiesFieldData.hpp.

◆ swapBaseDiffNPtr

boost::shared_ptr<MatrixDouble> MoFEM::EntitiesFieldData::EntData::swapBaseDiffNPtr
protected

Used by Bernstein base to keep temporally pointer.

Swap approximation base. Bernstein-Bezier (BB) base is not hierarchical, and is calculated for particular field, since it all shape functions change with the order. BB base is precalculated for every field, and when user push operator with paricular field using BB base, pointers to shape funtions and BaseDerivatives of shape functions are set to particular location, once operator is executed, pointers are switch back to its oroginal position.

getNSharedPtr(base) <=== getBBNSharedPtr(field_name); // DO OPERATOR WORK getNSharedPtr(base) ==> getBBNSharedPtr(field_name);

Parameters
field_name
base
Returns
MoFEMErrorCode

Definition at line 1112 of file EntitiesFieldData.hpp.

◆ swapBaseNPtr

boost::shared_ptr<MatrixDouble> MoFEM::EntitiesFieldData::EntData::swapBaseNPtr
protected

Used by Bernstein base to keep temporally pointer.

Swap approximation base. Bernstein-Bezier (BB) base is not hierarchical, and is calculated for particular field, since it all shape functions change with the order. BB base is precalculated for every field, and when user push operator with paricular field using BB base, pointers to shape funtions and BaseDerivatives of shape functions are set to particular location, once operator is executed, pointers are switch back to its oroginal position.

getNSharedPtr(base) <=== getBBNSharedPtr(field_name); // DO OPERATOR WORK getNSharedPtr(base) ==> getBBNSharedPtr(field_name);

Parameters
field_name
base
Returns
MoFEMErrorCode

Definition at line 1105 of file EntitiesFieldData.hpp.


The documentation for this struct was generated from the following files:
NOSPACE
@ NOSPACE
Definition: definitions.h:83
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::EntitiesFieldData::EntData::getBBAlphaIndicesSharedPtr
virtual boost::shared_ptr< MatrixInt > & getBBAlphaIndicesSharedPtr(const std::string &field_name)
Definition: EntitiesFieldData.cpp:822
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N()
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1495
MoFEM::EntitiesFieldData::EntData::N
std::array< boost::shared_ptr< MatrixDouble >, LASTBASE > & N
Base functions.
Definition: EntitiesFieldData.hpp:1081
sdf_hertz.d
float d
Definition: sdf_hertz.py:5
LASTBASE
@ LASTBASE
Definition: definitions.h:69
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
HVEC0_2
@ HVEC0_2
Definition: definitions.h:198
MoFEM::EntitiesFieldData::EntData::sPace
FieldSpace sPace
Entity space.
Definition: EntitiesFieldData.hpp:1068
MoFEM::EntitiesFieldData::EntData::ForthDerivative
@ ForthDerivative
Definition: EntitiesFieldData.hpp:134
NOBASE
@ NOBASE
Definition: definitions.h:59
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:596
ApproximationBaseNames
const static char *const ApproximationBaseNames[]
Definition: definitions.h:72
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
HVEC0_1
@ HVEC0_1
Definition: definitions.h:195
MoFEM::EntitiesFieldData::EntData::bbFieldName
std::string bbFieldName
field name
Definition: EntitiesFieldData.hpp:1085
MoFEM::Types::VectorIntAdaptor
VectorShallowArrayAdaptor< int > VectorIntAdaptor
Definition: Types.hpp:116
FTensor::Tensor2_symmetric
Definition: Tensor2_symmetric_value.hpp:13
HVEC0_0
@ HVEC0_0
Definition: definitions.h:192
HVEC1_1
@ HVEC1_1
Definition: definitions.h:196
THROW_MESSAGE
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:561
MoFEM::EntitiesFieldData::EntData::sEnse
int sEnse
Entity sense (orientation)
Definition: EntitiesFieldData.hpp:1066
MoFEM::EntitiesFieldData::EntData::swapBaseNPtr
boost::shared_ptr< MatrixDouble > swapBaseNPtr
Used by Bernstein base to keep temporally pointer.
Definition: EntitiesFieldData.hpp:1105
MoFEM::EntitiesFieldData::EntData::getDiffN
MatrixDouble & getDiffN()
get derivatives of base functions
Definition: EntitiesFieldData.hpp:1339
HVEC1
@ HVEC1
Definition: definitions.h:186
FTensor::Tensor2
Definition: Tensor2_value.hpp:16
order
constexpr int order
Definition: dg_projection.cpp:18
HVEC2_1
@ HVEC2_1
Definition: definitions.h:197
MoFEM::EntitiesFieldData::EntData::getBBNSharedPtr
virtual boost::shared_ptr< MatrixDouble > & getBBNSharedPtr(const std::string &field_name)
Definition: EntitiesFieldData.cpp:828
FTensor::Tensor3
Definition: Tensor3_value.hpp:12
HVEC2_2
@ HVEC2_2
Definition: definitions.h:200
HVEC1_2
@ HVEC1_2
Definition: definitions.h:199
MoFEM::EntitiesFieldData::EntData::bbAlphaIndices
std::map< std::string, boost::shared_ptr< MatrixInt > > bbAlphaIndices
Indices for Bernstein-Bezier (BB) base.
Definition: EntitiesFieldData.hpp:1090
MoFEM::EntitiesFieldData::EntData::bAse
FieldApproximationBase bAse
Field approximation base.
Definition: EntitiesFieldData.hpp:1069
MoFEM::getVectorAdaptor
auto getVectorAdaptor(T1 ptr, const size_t n)
Get Vector adaptor.
Definition: Templates.hpp:31
MoFEM::EntitiesFieldData::EntData::ThirdDerivative
@ ThirdDerivative
Definition: EntitiesFieldData.hpp:133
MoFEM::EntitiesFieldData::EntData::iNdices
VectorInt iNdices
Global indices on entity.
Definition: EntitiesFieldData.hpp:1070
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
MoFEM::EntitiesFieldData::EntData::FirstDerivative
@ FirstDerivative
Definition: EntitiesFieldData.hpp:131
MoFEM::EntitiesFieldData::EntData::bbDiffN
std::map< std::string, boost::shared_ptr< MatrixDouble > > bbDiffN
Definition: EntitiesFieldData.hpp:1088
DIM1
constexpr int DIM1
Definition: level_set.cpp:21
MoFEM::EntitiesFieldData::EntData::fieldEntities
VectorFieldEntities fieldEntities
Field entities.
Definition: EntitiesFieldData.hpp:1073
MoFEM::EntitiesFieldData::EntData::bbNodeOrder
VectorInt bbNodeOrder
order of nodes
Definition: EntitiesFieldData.hpp:1086
MoFEM::EntitiesFieldData::EntData::getBBDiffNSharedPtr
virtual boost::shared_ptr< MatrixDouble > & getBBDiffNSharedPtr(const std::string &field_name)
Definition: EntitiesFieldData.cpp:845
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN()
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1313
MoFEM::EntitiesFieldData::EntData::getDiffNSharedPtr
virtual boost::shared_ptr< MatrixDouble > & getDiffNSharedPtr(const FieldApproximationBase base)
Definition: EntitiesFieldData.cpp:37
MoFEM::EntitiesFieldData::EntData::bbAlphaIndicesByOrder
std::array< boost::shared_ptr< MatrixInt >, MaxBernsteinBezierOrder > bbAlphaIndicesByOrder
BB alpha indices by order.
Definition: EntitiesFieldData.hpp:1097
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::EntitiesFieldData::EntData::bbNByOrder
std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > bbNByOrder
BB base functions by order.
Definition: EntitiesFieldData.hpp:1093
EntData
EntitiesFieldData::EntData EntData
Definition: child_and_parent.cpp:37
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', 3 >
MoFEM::EntitiesFieldData::EntData::bbDiffNByOrder
std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > bbDiffNByOrder
BB base functions direvatives by order.
Definition: EntitiesFieldData.hpp:1095
MoFEM::EntitiesFieldData::EntData::entDataBitRefLevel
std::vector< BitRefLevel > entDataBitRefLevel
Bit ref level in entity.
Definition: EntitiesFieldData.hpp:1075
FTensor::Tensor0
Definition: Tensor0.hpp:16
MoFEM::EntitiesFieldData::EntData::localIndices
VectorInt localIndices
Local indices on entity.
Definition: EntitiesFieldData.hpp:1071
MoFEM::EntitiesFieldData::EntData::ZeroDerivative
@ ZeroDerivative
Definition: EntitiesFieldData.hpp:130
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::EntitiesFieldData::EntData::LastDerivative
@ LastDerivative
Definition: EntitiesFieldData.hpp:135
Tensor1
Definition: single.cpp:40
HVEC2_0
@ HVEC2_0
Definition: definitions.h:194
MoFEM::EntitiesFieldData::EntData::baseFunctionsAndBaseDerivatives
std::array< std::array< boost::shared_ptr< MatrixDouble >, LASTBASE >, LastDerivative > baseFunctionsAndBaseDerivatives
Definition: EntitiesFieldData.hpp:1079
MoFEM::EntitiesFieldData::EntData::bbN
std::map< std::string, boost::shared_ptr< MatrixDouble > > bbN
Definition: EntitiesFieldData.hpp:1087
MoFEM::EntitiesFieldData::EntData::swapBaseDiffNPtr
boost::shared_ptr< MatrixDouble > swapBaseDiffNPtr
Used by Bernstein base to keep temporally pointer.
Definition: EntitiesFieldData.hpp:1112
MoFEM::EntitiesFieldData::EntData::oRder
ApproximationOrder oRder
Entity order.
Definition: EntitiesFieldData.hpp:1067
MoFEM::Types::MatrixAdaptor
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
Definition: Types.hpp:132
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::EntitiesFieldData::EntData::fieldData
VectorDouble fieldData
Field data on entity.
Definition: EntitiesFieldData.hpp:1074
HVEC2
@ HVEC2
Definition: definitions.h:186
MoFEM::EntitiesFieldData::EntData::dOfs
VectorDofs dOfs
DoFs on entity.
Definition: EntitiesFieldData.hpp:1072
MoFEM::EntitiesFieldData::EntData::SecondDerivative
@ SecondDerivative
Definition: EntitiesFieldData.hpp:132
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
HVEC1_0
@ HVEC1_0
Definition: definitions.h:193
HVEC0
@ HVEC0
Definition: definitions.h:186
MoFEM::EntitiesFieldData::EntData::diffN
std::array< boost::shared_ptr< MatrixDouble >, LASTBASE > & diffN
Derivatives of base functions.
Definition: EntitiesFieldData.hpp:1083
MoFEM::EntitiesFieldData::EntData::getNSharedPtr
virtual boost::shared_ptr< MatrixDouble > & getNSharedPtr(const FieldApproximationBase base, const BaseDerivatives direvatie)
Definition: EntitiesFieldData.cpp:27