10#ifndef __ENTITIES_FIELD_DATA_HPP__
11#define __ENTITIES_FIELD_DATA_HPP__
13using namespace boost::numeric;
23using VectorDofs = ublas::vector<FEDofEntity *, DofsAllocator>;
44 std::bitset<LASTBASE>
bAse;
48 std::array<std::bitset<LASTSPACE>, MBMAXTYPE>
50 std::array<std::bitset<LASTBASE>, MBMAXTYPE>
52 std::array<std::bitset<LASTBASE>,
LASTSPACE>
54 std::array<std::bitset<LASTBASE>,
LASTSPACE>
56 std::array<boost::ptr_vector<EntData>, MBMAXTYPE>
132 const boost::shared_ptr<EntitiesFieldData> &data_ptr);
143 const boost::shared_ptr<EntitiesFieldData>
dataPtr;
172 EntData(
const bool allocate_base_matrices =
true);
323 template <
int Tensor_Dim>
339 template <
int Tensor_Dim0,
int Tensor_Dim1>
341 Tensor_Dim0, Tensor_Dim1>
356 template <
int Tensor_Dim>
396 virtual boost::shared_ptr<MatrixDouble> &
406 virtual boost::shared_ptr<MatrixDouble> &
415 virtual boost::shared_ptr<MatrixDouble> &
563 const int gg,
const int nb_base_functions);
588 const int nb_base_functions);
598 const int nb_base_functions);
632 template <
int DIM0,
int DIM1>
642 template <
int DIM0,
int DIM1>
652 template <
int DIM0,
int DIM1>
654 const int dof,
const int gg);
662 template <
int DIM0,
int DIM1>
741 template <
int Tensor_Dim>
761 template <
int Tensor_Dim>
783 template <
int Tensor_Dim>
806 template <
int Tensor_Dim>
812 template <
int Tensor_Dim>
814 Tensor_Dim, Tensor_Dim>
820 template <
int Tensor_Dim>
822 Tensor_Dim, Tensor_Dim>
828 template <
int Tensor_Dim>
830 Tensor_Dim, Tensor_Dim>
836 template <
int Tensor_Dim>
838 Tensor_Dim, Tensor_Dim>
867 template <
int Tensor_Dim>
895 template <
int Tensor_Dim0,
int Tensor_Dim1>
897 Tensor_Dim0, Tensor_Dim1>
903 template <
int Tensor_Dim0,
int Tensor_Dim1>
905 Tensor_Dim0, Tensor_Dim1>
910 template <
int Tensor_Dim0,
int Tensor_Dim1>
912 Tensor_Dim0, Tensor_Dim1>
914 return getFTensor2DiffN<Tensor_Dim0, Tensor_Dim1>(
bAse);
920 template <
int Tensor_Dim0,
int Tensor_Dim1>
922 Tensor_Dim0, Tensor_Dim1>
924 return getFTensor2DiffN<Tensor_Dim0, Tensor_Dim1>(
bAse, gg, bb);
929 template <
int Tensor_Dim0,
int Tensor_Dim1,
int Tensor_Dim2>
932 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
937 template <
int Tensor_Dim0,
int Tensor_Dim1,
int Tensor_Dim2>
940 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
942 return getFTensor3Diff2N<Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>(
bAse);
947 template <
int Tensor_Dim0,
int Tensor_Dim1,
int Tensor_Dim2>
950 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
956 template <
int Tensor_Dim0,
int Tensor_Dim1,
int Tensor_Dim2>
959 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
964 template <
int Tensor_Dim0,
int Tensor_Dim1,
int Tensor_Dim2>
967 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
969 return getFTensor3DiffN<Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>(
bAse);
975 template <
int Tensor_Dim0,
int Tensor_Dim1,
int Tensor_Dim2>
978 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
980 return getFTensor3DiffN<Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>(
bAse, gg,
998 template <
int Tensor_Dim>
1016 template <
int Tensor_Dim>
1045 template <
int Tensor_Dim0,
int Tensor_Dim1>
1047 Tensor_Dim0, Tensor_Dim1>
1072 template <
int Tensor_Dim0,
int Tensor_Dim1>
auto getFTensor2N();
1100 template <
int Tensor_Dim0,
int Tensor_Dim1>
1102 Tensor_Dim0, Tensor_Dim1>
1126 template <
int Tensor_Dim0,
int Tensor_Dim1>
1135 friend std::ostream &
operator<<(std::ostream &os,
1164 virtual boost::shared_ptr<MatrixInt> &
1170 virtual boost::shared_ptr<MatrixDouble> &
1176 virtual const boost::shared_ptr<MatrixDouble> &
1182 virtual boost::shared_ptr<MatrixDouble> &
1188 virtual const boost::shared_ptr<MatrixDouble> &
1191 virtual std::map<std::string, boost::shared_ptr<MatrixInt>> &
1199 virtual std::map<std::string, boost::shared_ptr<MatrixDouble>> &
getBBNMap();
1207 virtual std::map<std::string, boost::shared_ptr<MatrixDouble>> &
1216 virtual boost::shared_ptr<MatrixInt> &
1225 virtual boost::shared_ptr<MatrixDouble> &
1234 virtual boost::shared_ptr<MatrixDouble> &
1281 std::vector<EntityType>
1284 std::array<std::array<boost::shared_ptr<MatrixDouble>,
LASTBASE>,
1289 std::array<boost::shared_ptr<MatrixDouble>,
LASTBASE>
1294 std::map<std::string, boost::shared_ptr<MatrixDouble>>
bbN;
1295 std::map<std::string, boost::shared_ptr<MatrixDouble>>
bbDiffN;
1296 std::map<std::string, boost::shared_ptr<MatrixInt>>
1340 const boost::shared_ptr<EntitiesFieldData::EntData> &ent_data_ptr);
1347 boost::shared_ptr<MatrixDouble> &
1351 boost::shared_ptr<MatrixDouble> &
1354 boost::shared_ptr<MatrixDouble> &
1357 const boost::shared_ptr<MatrixDouble> &
1360 const boost::shared_ptr<MatrixDouble> &
1363 inline boost::shared_ptr<MatrixDouble> &
1366 inline boost::shared_ptr<MatrixDouble> &
1369 boost::shared_ptr<MatrixInt> &
1375 boost::shared_ptr<MatrixDouble> &
1381 const boost::shared_ptr<MatrixDouble> &
1387 boost::shared_ptr<MatrixDouble> &
1393 const boost::shared_ptr<MatrixDouble> &
1420 unsigned int size = 0;
1421 if (
auto dof = dOfs[0]) {
1422 size = dof->getOrderNbDofs(
order) * dof->getNbOfCoeffs();
1423 size = size < iNdices.size() ? size : iNdices.size();
1425 int *data = &*iNdices.data().begin();
1426 return VectorIntAdaptor(size, ublas::shallow_array_adaptor<int>(size, data));
1430 return localIndices;
1435 unsigned int size = 0;
1436 if (
auto dof = dOfs[0]) {
1437 size = dof->getOrderNbDofs(
order) * dof->getNbOfCoeffs();
1438 size = size < localIndices.size() ? size : localIndices.size();
1440 int *data = &*localIndices.data().begin();
1441 return VectorIntAdaptor(size, ublas::shallow_array_adaptor<int>(size, data));
1451 return localIndices;
1460 unsigned int size = 0;
1461 if (
auto dof = dOfs[0]) {
1462 size = dof->getOrderNbDofs(
order) * dof->getNbOfCoeffs();
1463 size = size < fieldData.size() ? size : fieldData.size();
1465 double *data = &*fieldData.data().begin();
1478 return fieldEntities;
1483 return fieldEntities;
1486template <
int Tensor_Dim>
1489 std::stringstream s;
1490 s <<
"Not implemented for this dimension dim = " << Tensor_Dim;
1494template <
int Tensor_Dim0,
int Tensor_Dim1>
1496 Tensor_Dim0, Tensor_Dim1>
1498 std::stringstream s;
1499 s <<
"Not implemented for this dimension dim0 = " << Tensor_Dim0;
1500 s <<
" and dim1 " << Tensor_Dim1;
1504template <
int Tensor_Dim>
1506 FTensor::PackPtr<
double *, (Tensor_Dim * (Tensor_Dim + 1)) / 2>, Tensor_Dim>
1508 std::stringstream s;
1509 s <<
"Not implemented for this dimension dim = " << Tensor_Dim;
1520 if (!getNSharedPtr(base)) {
1521 MOFEM_LOG_C(
"SELF", Sev::error,
"Ptr to base %s functions is null",
1524 "Null pointer to base functions");
1527 return *(getNSharedPtr(base));
1533 MOFEM_LOG_C(
"SELF", Sev::error,
"Ptr to field %s base functions is null ",
1536 "Null pointer to base functions");
1546 return *(getDiffNSharedPtr(base));
1553 if (!getNSharedPtr(base, derivative)) {
1555 "Ptr to base %s functions derivative %d is null",
1560 return *(getNSharedPtr(base, derivative));
1568 "Null pointer to base functions derivative");
1578 return getN(
bAse, derivative);
1584 int size = getN(base).size2();
1585 double *data = &getN(base)(gg, 0);
1586 return VectorAdaptor(size, ublas::shallow_array_adaptor<double>(size, data));
1590 return getN(
bAse, gg);
1599 if (getN(base).size1() == getDiffN(base).size1()) {
1600 int size = getN(base).size2();
1601 int dim = getDiffN(base).size2() / size;
1602 double *data = &getDiffN(base)(gg, 0);
1604 getN(base).size2(), dim,
1605 ublas::shallow_array_adaptor<double>(getDiffN(base).size2(), data));
1610 getN(base).size1(), getN(base).size2(),
1611 ublas::shallow_array_adaptor<double>(getDiffN(base).data().size(),
1612 &getDiffN(base).data()[0]));
1617 return getDiffN(
bAse, gg);
1622 const int gg,
const int nb_base_functions) {
1623 (void)getN()(gg, nb_base_functions -
1625 double *data = &getN(base)(gg, 0);
1626 return VectorAdaptor(nb_base_functions, ublas::shallow_array_adaptor<double>(
1627 nb_base_functions, data));
1632 return getN(
bAse, gg, nb_base_functions);
1638 const int nb_base_functions) {
1642 if (getN(base).size1() == getDiffN(base).size1()) {
1643 (void)getN(base)(gg,
1646 int dim = getDiffN(base).size2() / getN(base).size2();
1647 double *data = &getDiffN(base)(gg, 0);
1649 nb_base_functions, dim,
1650 ublas::shallow_array_adaptor<double>(dim * nb_base_functions, data));
1655 getN(base).size1(), getN(base).size2(),
1656 ublas::shallow_array_adaptor<double>(getDiffN(base).data().size(),
1657 &getDiffN(base).data()[0]));
1663 const int nb_base_functions) {
1664 return getDiffN(
bAse, gg, nb_base_functions);
1671 if (PetscUnlikely(getN(base).size2() % DIM)) {
1675 const int nb_base_functions = getN(base).size2() / DIM;
1676 double *data = &getN(base)(gg, 0);
1678 nb_base_functions, DIM,
1679 ublas::shallow_array_adaptor<double>(DIM * nb_base_functions, data));
1684 return getVectorN<DIM>(
bAse, gg);
1687template <
int DIM0,
int DIM1>
1691 if (PetscUnlikely(getDiffN(base).size2() % (DIM0 *
DIM1))) {
1695 const int nb_base_functions = getN(base).size2() / (DIM0 *
DIM1);
1696 double *data = &getN(base)(gg, 0);
1698 ublas::shallow_array_adaptor<double>(
1699 DIM0 *
DIM1 * nb_base_functions, data));
1702template <
int DIM0,
int DIM1>
1704 return getVectorDiffN<DIM0, DIM1>(
bAse, gg);
1707template <
int DIM0,
int DIM1>
1710 const int dof,
const int gg) {
1714 ublas::shallow_array_adaptor<double>(DIM0 *
DIM1, data));
1717template <
int DIM0,
int DIM1>
1720 return getVectorDiffN<DIM0, DIM1>(
bAse, dof, gg);
1725 double *ptr = &*getN(base).data().begin();
1731 return getFTensor0N(
bAse);
1736 const int gg,
const int bb) {
1737 double *ptr = &getN(base)(gg, bb);
1743 return getFTensor0N(
bAse, gg, bb);
1747 return getFTensor1N<Tensor_Dim>(
bAse);
1750template <
int Tensor_Dim>
1752 return getFTensor1N<Tensor_Dim>(
bAse, gg, bb);
1755template <
int Tensor_Dim0,
int Tensor_Dim1>
1757 return getFTensor2N<Tensor_Dim0, Tensor_Dim1>(
bAse);
1760template <
int Tensor_Dim0,
int Tensor_Dim1>
1762 return getFTensor2N<Tensor_Dim0, Tensor_Dim1>(
bAse, gg, bb);
1772 return *getBBAlphaIndicesSharedPtr(bbFieldName);
1781boost::shared_ptr<MatrixDouble> &
1787boost::shared_ptr<MatrixDouble> &
1810template <
typename T = EntityStorage>
1813 const double *ptr, InsertMode iora) {
1814 static_assert(!std::is_same<T, T>::value,
1815 "VecSetValues value for this data storage is not implemented");
1822 const double *ptr, InsertMode iora) {
1842template <
typename T = EntityStorage>
1846 return VecSetValues<T>(V, data, &*vec.data().begin(), iora);
1865template <
typename T = EntityStorage>
1869 const double *ptr, InsertMode iora) {
1870 static_assert(!std::is_same<T, T>::value,
1871 "MatSetValues value for this data storage is not implemented");
1891template <
typename T = EntityStorage>
1896 return MatSetValues<T>(M, row_data, col_data, &*mat.data().begin(), iora);
1903 const double *ptr, InsertMode iora) {
1920 const int gg,
const int bb);
1934EntitiesFieldData::EntData::getFTensor1DiffN<3>(
1938EntitiesFieldData::EntData::getFTensor1DiffN<3>();
1942EntitiesFieldData::EntData::getFTensor1DiffN<2>(
1946EntitiesFieldData::EntData::getFTensor1DiffN<2>();
1954 const int gg,
const int bb);
1961 const int gg,
const int bb);
1974 const int gg,
const int bb);
1978EntitiesFieldData::EntData::getFTensor2DiffN2<2>(
1983EntitiesFieldData::EntData::getFTensor2DiffN2<2>(
1988EntitiesFieldData::EntData::getFTensor2DiffN2<2>();
1992EntitiesFieldData::EntData::getFTensor2DiffN2<2>(
const int gg,
const int bb);
1996EntitiesFieldData::EntData::getFTensor2DiffN2<3>(
2001EntitiesFieldData::EntData::getFTensor2DiffN2<3>(
2006EntitiesFieldData::EntData::getFTensor2DiffN2<3>();
2010EntitiesFieldData::EntData::getFTensor2DiffN2<3>(
const int gg,
const int bb);
2020EntitiesFieldData::EntData::getFTensor1FieldData<3>();
2024EntitiesFieldData::EntData::getFTensor1FieldData<2>();
2028EntitiesFieldData::EntData::getFTensor1FieldData<1>();
2032EntitiesFieldData::EntData::getFTensor2FieldData<1, 1>();
2036EntitiesFieldData::EntData::getFTensor2FieldData<1, 2>();
2040EntitiesFieldData::EntData::getFTensor2FieldData<1, 3>();
2044EntitiesFieldData::EntData::getFTensor2FieldData<2, 2>();
2048EntitiesFieldData::EntData::getFTensor2FieldData<3, 3>();
2052EntitiesFieldData::EntData::getFTensor2SymmetricFieldData<3>();
2056EntitiesFieldData::EntData::getFTensor2SymmetricFieldData<2>();
#define MOFEM_LOG_C(channel, severity, format,...)
FieldApproximationBase
approximation base
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
FieldSpace
approximation spaces
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
@ MOFEM_DATA_INCONSISTENCY
static const char *const ApproximationBaseNames[]
#define BITFEID_SIZE
max number of finite elements
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VectorShallowArrayAdaptor< int > VectorIntAdaptor
int ApproximationOrder
Approximation on the entity.
VectorShallowArrayAdaptor< double > VectorAdaptor
UBlasMatrix< double > MatrixDouble
UBlasVector< double > VectorDouble
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
UBlasMatrix< int > MatrixInt
UBlasVector< int > VectorInt
implementation of Data Operators for Forces and Sources
auto getVectorAdaptor(T1 ptr, const size_t n)
Get Vector adaptor.
DEPRECATED typedef DerivedEntitiesFieldData DerivedDataForcesAndSourcesCore
ublas::vector< FieldEntity *, FieldEntAllocator > VectorFieldEntities
DEPRECATED typedef EntitiesFieldData DataForcesAndSourcesCore
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
ublas::unbounded_array< FieldEntity *, std::allocator< FieldEntity * > > FieldEntAllocator
MoFEMErrorCode VecSetValues< EntityStorage >(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
MoFEMErrorCode MatSetValues< EntityStorage >(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
ublas::unbounded_array< FEDofEntity *, std::allocator< FEDofEntity * > > DofsAllocator
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
constexpr auto field_name
Derived data on single entity (This is passed as argument to DataOperator::doWork)
boost::shared_ptr< MatrixDouble > & getBBNSharedPtr(const std::string &field_name)
boost::shared_ptr< MatrixDouble > & getDerivedDiffNSharedPtr(const FieldApproximationBase base)
const boost::shared_ptr< EntitiesFieldData::EntData > entDataPtr
MoFEMErrorCode baseSwap(const std::string &field_name, const FieldApproximationBase base)
Used by Bernstein base to keep temporary pointer.
std::vector< BitRefLevel > & getEntDataBitRefLevel()
Get entity bit reference level.
int getSense() const
Get entity sense for conforming approximation fields.
boost::shared_ptr< MatrixInt > & getBBAlphaIndicesSharedPtr(const std::string &field_name)
boost::shared_ptr< MatrixDouble > & getBBDiffNSharedPtr(const std::string &field_name)
boost::shared_ptr< MatrixDouble > & getDiffNSharedPtr(const FieldApproximationBase base)
Get shared pointer to derivatives of base functions.
boost::shared_ptr< MatrixDouble > & getNSharedPtr(const FieldApproximationBase base, const BaseDerivatives derivative)
Get shared pointer to base functions with derivatives.
boost::shared_ptr< MatrixDouble > & getDerivedNSharedPtr(const FieldApproximationBase base)
this class derives data from other data structure
MoFEMErrorCode setElementType(const EntityType type)
Set element type for derived data.
const boost::shared_ptr< EntitiesFieldData > dataPtr
Data on single entity (This is passed as argument to DataOperator::doWork)
std::vector< int > dofBrokenSideVec
Map side to dofs number.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim *Tensor_Dim >, Tensor_Dim, Tensor_Dim > getFTensor2DiffN2(const FieldApproximationBase base)
Get second derivatives of scalar base functions.
std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > bbDiffNByOrder
BB base functions derivatives by order.
boost::shared_ptr< MatrixDouble > swapBaseDiffNPtr
Used by Bernstein base to keep temporary pointer.
virtual std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > & getBBDiffNByOrderArray()
FieldSpace sPace
Entity space.
const VectorIntAdaptor getIndicesUpToOrder(int order)
Get global indices of DOFs on entity up to given order.
virtual std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > & getBBNByOrderArray()
int sEnse
Entity sense (orientation)
auto getFTensor2N()
Get base functions for Hdiv space.
VectorDouble fieldData
Field data on entity.
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)
MatrixDouble & getN()
Get base functions using default base.
FieldApproximationBase & getBase()
Get approximation base.
boost::shared_ptr< MatrixDouble > swapBaseNPtr
Used by Bernstein base to keep temporary pointer.
ApproximationOrder getOrder() const
Get approximation order.
std::map< std::string, boost::shared_ptr< MatrixInt > > bbAlphaIndices
Indices for Bernstein-Bezier (BB) base.
FTensor::Tensor3< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 *Tensor_Dim2 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3DiffN(const int gg, const int bb)
Get derivatives of base functions for tonsorial Hdiv space at integration pts.
friend std::ostream & operator<<(std::ostream &os, const EntitiesFieldData::EntData &e)
VectorInt localIndices
Local indices on entity.
const MatrixAdaptor getVectorDiffN(FieldApproximationBase base, const int gg)
get DiffHdiv of base functions at Gauss pts
const VectorFieldEntities & getFieldEntities() const
Get field entities (const version)
virtual std::array< boost::shared_ptr< MatrixInt >, MaxBernsteinBezierOrder > & getBBAlphaIndicesByOrderArray()
virtual std::map< std::string, boost::shared_ptr< MatrixInt > > & getBBAlphaIndicesMap()
const VectorIntAdaptor getLocalIndicesUpToOrder(int order)
Get local indices of DOFs on entity up to given order.
const VectorAdaptor getFieldDataUpToOrder(int order)
Get DOF values on entity up to given order.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN()
Get derivatives of base functions.
std::vector< BitRefLevel > entDataBitRefLevel
Bit ref level in entity.
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.
static constexpr size_t MaxBernsteinBezierOrder
VectorInt iNdices
Global indices on entity.
virtual boost::shared_ptr< MatrixDouble > & getNSharedPtr(const FieldApproximationBase base, const BaseDerivatives derivative)
Get shared pointer to base functions with derivatives.
virtual boost::shared_ptr< MatrixInt > & getBBAlphaIndicesSharedPtr(const std::string &field_name)
const VectorDouble & getFieldData() const
Get DOF values on entity.
VectorInt bbNodeOrder
order of nodes
virtual boost::shared_ptr< MatrixDouble > & getBBDiffNByOrderSharedPtr(const size_t o)
get BB base derivative by order
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.
virtual std::map< std::string, boost::shared_ptr< MatrixDouble > > & getBBNMap()
get hash map of base function for BB base, key is a field name
const MatrixAdaptor getVectorN(const FieldApproximationBase base, const int gg)
get Hdiv of base functions at Gauss pts
virtual MoFEMErrorCode baseSwap(const std::string &field_name, const FieldApproximationBase base)
Swap bases functions.
const VectorInt & getLocalIndices() const
Get local indices of degrees of freedom on entity.
ApproximationOrder oRder
Entity order.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1FieldData()
Return FTensor of rank 1, i.e. vector from field data coefficients.
std::vector< EntityType > dofBrokenTypeVec
Map type of entity to dof number.
FTensor::Tensor3< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 *Tensor_Dim2 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3DiffN()
Get derivatives of base functions for tonsorial Hdiv space.
virtual std::map< std::string, boost::shared_ptr< MatrixDouble > > & getBBDiffNMap()
get hash map of derivatives base function for BB base, key is a field name
virtual std::vector< BitRefLevel > & getEntDataBitRefLevel()
Get entity bit reference level.
virtual int getSense() const
Get entity sense for conforming approximation fields.
VectorFieldEntities fieldEntities
Field entities.
virtual ~EntData()=default
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N()
Get base function as Tensor0.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData()
Return scalar files as a FTensor of rank 0.
MatrixInt & getBBAlphaIndices()
Get file BB indices.
std::array< boost::shared_ptr< MatrixInt >, MaxBernsteinBezierOrder > bbAlphaIndicesByOrder
BB alpha indices by order.
MatrixDouble & getDiffN()
Get derivatives of base functions using default base.
VectorInt & getBBNodeOrder()
Get orders at the nodes.
std::string bbFieldName
field name
virtual boost::shared_ptr< MatrixInt > & getBBAlphaIndicesByOrderSharedPtr(const size_t o)
get ALpha indices for BB base by order
std::array< boost::shared_ptr< MatrixDouble >, LASTBASE > & diffN
Derivatives of base functions.
std::map< std::string, boost::shared_ptr< MatrixDouble > > bbDiffN
MoFEMErrorCode resetFieldDependentData()
const VectorDofs & getFieldDofs() const
Get DOF data structures (const version)
FTensor::Tensor2_symmetric< FTensor::PackPtr< double *,(Tensor_Dim *(Tensor_Dim+1))/2 >, Tensor_Dim > getFTensor2SymmetricFieldData()
Return symmetric FTensor rank 2, i.e. matrix from field data coefficients.
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.
FieldSpace & getSpace()
Get field space.
auto getFTensor1N()
Get base functions for Hdiv space.
virtual boost::shared_ptr< MatrixDouble > & getBBNSharedPtr(const std::string &field_name)
std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > bbNByOrder
BB base functions by order.
std::array< std::array< boost::shared_ptr< MatrixDouble >, LASTBASE >, LastDerivative > baseFunctionsAndBaseDerivatives
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN()
Get derivatives of base functions for Hdiv space.
std::map< std::string, boost::shared_ptr< MatrixDouble > > bbN
VectorDofs dOfs
DoFs on entity.
virtual boost::shared_ptr< MatrixDouble > & getDiffNSharedPtr(const FieldApproximationBase base)
Get shared pointer to derivatives of base functions.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim *Tensor_Dim >, Tensor_Dim, Tensor_Dim > getFTensor2DiffN2()
Get second derivatives of scalar base functions.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2FieldData()
Return FTensor rank 2, i.e. matrix from field data coefficients.
virtual boost::shared_ptr< MatrixDouble > & getBBNByOrderSharedPtr(const size_t o)
get BB base by order
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim *Tensor_Dim >, Tensor_Dim, Tensor_Dim > getFTensor2DiffN2(const int gg, const int bb)
Get second derivatives of scalar base functions at integration pts.
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
FieldApproximationBase bAse
Field approximation base.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim *Tensor_Dim >, Tensor_Dim, Tensor_Dim > getFTensor2DiffN2(const FieldApproximationBase base, const int gg, const int bb)
Get second derivatives of scalar base functions at integration pts.
virtual boost::shared_ptr< MatrixDouble > & getBBDiffNSharedPtr(const std::string &field_name)
std::array< boost::shared_ptr< MatrixDouble >, LASTBASE > & N
Base functions.
data structure for finite element entity
friend std::ostream & operator<<(std::ostream &os, const EntitiesFieldData &e)
MoFEMErrorCode resetFieldDependentData()
Reset data associated with particular field name.
MatrixInt facesNodes
nodes on finite element faces
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
std::bitset< LASTBASE > bAse
bases on element
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
std::array< std::bitset< LASTBASE >, LASTSPACE > basesOnSpaces
base on spaces
virtual ~EntitiesFieldData()=default
virtual MoFEMErrorCode baseSwap(const std::string &field_name, const FieldApproximationBase base)
Swap approximation base.
virtual MoFEMErrorCode setElementType(const EntityType type)
Set element type and initialize data structures.
std::array< std::bitset< LASTBASE >, LASTSPACE > brokenBasesOnSpaces
base on spaces
MatrixInt facesNodesOrder
order of face nodes on element
std::array< std::bitset< LASTBASE >, MBMAXTYPE > basesOnEntities
bases on entity types
keeps information about indexed dofs for the finite element
Struct keeps handle to entity in the field.
Operator to project base functions from parent entity to child.