v0.15.0
Loading...
Searching...
No Matches
MoFEM::OpCalculateVectorFieldHessian< Tensor_Dim0, Tensor_Dim1 > Struct Template Reference

#include "src/finite_elements/UserDataOperators.hpp"

Inheritance diagram for MoFEM::OpCalculateVectorFieldHessian< Tensor_Dim0, Tensor_Dim1 >:
[legend]
Collaboration diagram for MoFEM::OpCalculateVectorFieldHessian< Tensor_Dim0, Tensor_Dim1 >:
[legend]

Public Member Functions

Gradients and hessian of tensor fields at integration points

}

MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 calculate values of vector field at integration points
 
- Public Member Functions inherited from MoFEM::OpCalculateTensor2FieldValues_General< Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator >
 OpCalculateTensor2FieldValues_General (const std::string field_name, boost::shared_ptr< ublas::matrix< double, ublas::row_major, DoubleAllocator > > data_ptr, const EntityType zero_type=MBVERTEX)
 
 OpCalculateTensor2FieldValues_General (const std::string field_name, boost::shared_ptr< ublas::matrix< double, ublas::row_major, DoubleAllocator > > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBVERTEX)
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement.
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle.
 
int getFEDim () const
 Get dimension of finite element.
 
EntityType getFEType () const
 Get dimension of finite element.
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer.
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity.
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element.
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices.
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices.
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object.
 
int getOpType () const
 Get operator types.
 
void setOpType (const OpType type)
 Set operator type.
 
void addOpType (const OpType type)
 Add operator type.
 
int getNinTheLoop () const
 get number of finite element in the loop
 
int getLoopSize () const
 get size of elements in the loop
 
std::string getFEName () const
 Get name of the element.
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element
 
auto getFTensor0IntegrationWeight ()
 Get integration weights.
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3)
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry.
 
double getMeasure () const
 get measure of element
 
doublegetMeasure ()
 get measure of element
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, boost::shared_ptr< Range > fe_range=nullptr, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopRange (const string &fe_name, ForcesAndSourcesCore *range_fe, boost::shared_ptr< Range > fe_range, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 Iterate over range of elements.
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side.
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not.
 
void setSymm ()
 set if operator is executed taking in account symmetry
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType
 
using DoWorkRhsHookFunType
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure.
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity.
 
booldoVertices
 \deprectaed If false skip vertices
 
booldoEdges
 \deprectaed If false skip edges
 
booldoQuads
 \deprectaed
 
booldoTris
 \deprectaed
 
booldoTets
 \deprectaed
 
booldoPrisms
 \deprectaed
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::OpCalculateTensor2FieldValues_General< Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator >
boost::shared_ptr< ublas::matrix< double, ublas::row_major, DoubleAllocator > > dataPtr
 
const EntityHandle zeroType
 
SmartPetscObj< Vec > dataVec
 
VectorDouble dotVector
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

template<int Tensor_Dim0, int Tensor_Dim1>
struct MoFEM::OpCalculateVectorFieldHessian< Tensor_Dim0, Tensor_Dim1 >

Definition at line 1784 of file UserDataOperators.hpp.

Member Function Documentation

◆ doWork()

template<int Tensor_Dim0, int Tensor_Dim1>
MoFEMErrorCode MoFEM::OpCalculateVectorFieldHessian< Tensor_Dim0, Tensor_Dim1 >::doWork ( int side,
EntityType type,
EntitiesFieldData::EntData & data )
virtual

calculate values of vector field at integration points

Member function specialization calculating vector field gradients for tenor field rank 2.

Parameters
sideside entity number
typeside entity type
dataentity data
Returns
error code

Reimplemented from MoFEM::OpCalculateTensor2FieldValues_General< Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator >.

Definition at line 1809 of file UserDataOperators.hpp.

1810 {
1812 if (!this->dataPtr)
1813 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1814 "Data pointer not allocated");
1815
1816 const size_t nb_gauss_pts = this->getGaussPts().size2();
1817 auto &mat = *this->dataPtr;
1818 if (type == this->zeroType) {
1819 mat.resize(Tensor_Dim0 * Tensor_Dim1 * Tensor_Dim1, nb_gauss_pts, false);
1820 mat.clear();
1821 }
1822
1823 if (nb_gauss_pts) {
1824 const size_t nb_dofs = data.getFieldData().size();
1825
1826 if (nb_dofs) {
1827
1828 if (this->dataVec.use_count()) {
1829 this->dotVector.resize(nb_dofs, false);
1830 const double *array;
1831 CHKERR VecGetArrayRead(this->dataVec, &array);
1832 const auto &local_indices = data.getLocalIndices();
1833 for (int i = 0; i != local_indices.size(); ++i)
1834 if (local_indices[i] != -1)
1835 this->dotVector[i] = array[local_indices[i]];
1836 else
1837 this->dotVector[i] = 0;
1838 CHKERR VecRestoreArrayRead(this->dataVec, &array);
1839 data.getFieldData().swap(this->dotVector);
1840 }
1841
1842 const int nb_base_functions = data.getN().size2();
1843
1844 auto &hessian_base = data.getN(BaseDerivatives::SecondDerivative);
1845 #ifndef NDEBUG
1846 if (hessian_base.size1() != nb_gauss_pts) {
1847 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1848 "Wrong number of integration pts (%d != %d)",
1849 hessian_base.size1(), nb_gauss_pts);
1850 }
1851 if (hessian_base.size2() !=
1852 nb_base_functions * Tensor_Dim1 * Tensor_Dim1) {
1853 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1854 "Wrong number of base functions (%d != %d)",
1855 hessian_base.size2() / (Tensor_Dim1 * Tensor_Dim1),
1856 nb_base_functions);
1857 }
1858 if (hessian_base.size2() < nb_dofs * Tensor_Dim1 * Tensor_Dim1) {
1859 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1860 "Wrong number of base functions (%d < %d)",
1861 hessian_base.size2(), nb_dofs * Tensor_Dim1 * Tensor_Dim1);
1862 }
1863 #endif
1864
1865 auto t_diff2_base_function = getFTensor2FromPtr<Tensor_Dim1, Tensor_Dim1>(
1866 &*hessian_base.data().begin());
1867
1868 auto t_hessian_at_gauss_pts =
1870
1871 FTensor::Index<'I', Tensor_Dim0> I;
1872 FTensor::Index<'J', Tensor_Dim1> J;
1873 FTensor::Index<'K', Tensor_Dim1> K;
1874
1875 int size = nb_dofs / Tensor_Dim0;
1876 #ifndef NDEBUG
1877 if (nb_dofs % Tensor_Dim0) {
1878 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1879 "Data inconsistency");
1880 }
1881 #endif
1882
1883 for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1884 auto field_data = data.getFTensor1FieldData<Tensor_Dim0>();
1885 int bb = 0;
1886 for (; bb < size; ++bb) {
1887 t_hessian_at_gauss_pts(I, J, K) +=
1888 field_data(I) * t_diff2_base_function(J, K);
1889 ++field_data;
1890 ++t_diff2_base_function;
1891 }
1892 // Number of dofs can be smaller than number of Tensor_Dim0 x base
1893 // functions
1894 for (; bb != nb_base_functions; ++bb)
1895 ++t_diff2_base_function;
1896 ++t_hessian_at_gauss_pts;
1897 }
1898
1899 if (this->dataVec.use_count()) {
1900 data.getFieldData().swap(this->dotVector);
1901 }
1902 }
1903 }
1905}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
static FTensor::Tensor3< FTensor::PackPtr< T *, 1 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3FromMat(ublas::matrix< T, L, A > &data)
Get tensor rank 3 (non symmetries) form data matrix.
auto getFTensor2FromPtr(double *ptr)
Make Tensor2 from pointer.
constexpr IntegrationType I
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

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