v0.14.0
Public Member Functions | Protected Attributes | List of all members
MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl< Tensor_Dim, CTX > Struct Template Reference

Approximate field values for given petsc vector. More...

#include <src/finite_elements/UserDataOperators.hpp>

Inheritance diagram for MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl< Tensor_Dim, CTX >:
[legend]
Collaboration diagram for MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl< Tensor_Dim, CTX >:
[legend]

Public Member Functions

 OpCalculateVectorFieldValuesFromPetscVecImpl (const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_at_type=MBVERTEX)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
- 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. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
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 More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, 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. More...
 
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. More...
 
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. More...
 
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. More...
 
- 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. More...
 
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. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Protected Attributes

boost::shared_ptr< MatrixDoubledataPtr
 
const EntityHandle zeroAtType
 
VectorDouble dotVector
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

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 = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > >>
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- 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. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 
- 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)
 

Detailed Description

template<int Tensor_Dim, PetscData::DataContext CTX>
struct MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl< Tensor_Dim, CTX >

Approximate field values for given petsc vector.

Note
Look at PetscData to see what vectors could be extracted with that user data operator.
Examples
approx_sphere.cpp, EshelbianPlasticity.cpp, free_surface.cpp, operators_tests.cpp, plastic.cpp, and shallow_wave.cpp.

Definition at line 595 of file UserDataOperators.hpp.

Constructor & Destructor Documentation

◆ OpCalculateVectorFieldValuesFromPetscVecImpl()

template<int Tensor_Dim, PetscData::DataContext CTX>
MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl< Tensor_Dim, CTX >::OpCalculateVectorFieldValuesFromPetscVecImpl ( const std::string  field_name,
boost::shared_ptr< MatrixDouble data_ptr,
const EntityType  zero_at_type = MBVERTEX 
)
inline

Definition at line 598 of file UserDataOperators.hpp.

603  dataPtr(data_ptr), zeroAtType(zero_at_type) {
604  if (!dataPtr)
605  THROW_MESSAGE("Pointer is not set");
606  }

Member Function Documentation

◆ doWork()

template<int Tensor_Dim, PetscData::DataContext CTX>
MoFEMErrorCode MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl< Tensor_Dim, CTX >::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)
inlinevirtual

Operator for linear form, usually to calculate values on right hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 608 of file UserDataOperators.hpp.

609  {
611  auto &local_indices = data.getLocalIndices();
612  const size_t nb_dofs = local_indices.size();
613  const size_t nb_gauss_pts = getGaussPts().size2();
614 
615  MatrixDouble &mat = *dataPtr;
616  if (type == zeroAtType) {
617  mat.resize(Tensor_Dim, nb_gauss_pts, false);
618  mat.clear();
619  }
620  if (!nb_dofs)
622 
623  const double *array;
624 
625  auto get_array = [&](const auto ctx, auto vec) {
627 #ifndef NDEBUG
628  if ((getFEMethod()->data_ctx & ctx).none()) {
629  MOFEM_LOG_CHANNEL("SELF");
630  MOFEM_LOG("SELF", Sev::error)
631  << "In this case filed degrees of freedom are read from vector. "
632  "That usually happens when time solver is used, and acces to "
633  "first or second rates is needed. You probably not set ts_u, "
634  "ts_u_t, or ts_u_tt and associated data structure, i.e. "
635  "data_ctx to CTX_SET_X, CTX_SET_X_T, or CTX_SET_X_TT "
636  "respectively";
637  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Vector not set");
638  }
639 #endif
640  CHKERR VecGetArrayRead(vec, &array);
642  };
643 
644  auto restore_array = [&](auto vec) {
645  return VecRestoreArrayRead(vec, &array);
646  };
647 
648  switch (CTX) {
650  CHKERR get_array(PetscData::CtxSetX, getFEMethod()->ts_u);
651  break;
653  CHKERR get_array(PetscData::CtxSetX_T, getFEMethod()->ts_u_t);
654  break;
656  CHKERR get_array(PetscData::CtxSetX_TT, getFEMethod()->ts_u_tt);
657  break;
658  default:
659  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
660  "That case is not implemented");
661  }
662 
663  dotVector.resize(local_indices.size());
664  for (int i = 0; i != local_indices.size(); ++i)
665  if (local_indices[i] != -1)
666  dotVector[i] = array[local_indices[i]];
667  else
668  dotVector[i] = 0;
669 
670  switch (CTX) {
672  CHKERR restore_array(getFEMethod()->ts_u);
673  break;
675  CHKERR restore_array(getFEMethod()->ts_u_t);
676  break;
678  CHKERR restore_array(getFEMethod()->ts_u_tt);
679  break;
680  default:
681  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
682  "That case is not implemented");
683  }
684 
685  const size_t nb_base_functions = data.getN().size2();
686  auto base_function = data.getFTensor0N();
687  auto values_at_gauss_pts = getFTensor1FromMat<Tensor_Dim>(mat);
688 
690  const size_t size = nb_dofs / Tensor_Dim;
691  if (nb_dofs % Tensor_Dim) {
692  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Data inconsistency");
693  }
694  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
695  auto field_data = getFTensor1FromArray<Tensor_Dim, Tensor_Dim>(dotVector);
696  size_t bb = 0;
697  for (; bb != size; ++bb) {
698  values_at_gauss_pts(I) += field_data(I) * base_function;
699  ++field_data;
700  ++base_function;
701  }
702  // Number of dofs can be smaller than number of Tensor_Dim x base
703  // functions
704  for (; bb < nb_base_functions; ++bb)
705  ++base_function;
706  ++values_at_gauss_pts;
707  }
709  }

Member Data Documentation

◆ dataPtr

template<int Tensor_Dim, PetscData::DataContext CTX>
boost::shared_ptr<MatrixDouble> MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl< Tensor_Dim, CTX >::dataPtr
protected

Definition at line 712 of file UserDataOperators.hpp.

◆ dotVector

template<int Tensor_Dim, PetscData::DataContext CTX>
VectorDouble MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl< Tensor_Dim, CTX >::dotVector
protected

Definition at line 714 of file UserDataOperators.hpp.

◆ zeroAtType

template<int Tensor_Dim, PetscData::DataContext CTX>
const EntityHandle MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl< Tensor_Dim, CTX >::zeroAtType
protected

Definition at line 713 of file UserDataOperators.hpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
THROW_MESSAGE
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:574
I
constexpr IntegrationType I
Definition: operators_tests.cpp:31
MoFEM::ForcesAndSourcesCore::UserDataOperator::getGaussPts
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Definition: ForcesAndSourcesCore.hpp:1236
MoFEM::PetscData::CTX_SET_X_TT
@ CTX_SET_X_TT
Definition: LoopMethods.hpp:29
MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl::dataPtr
boost::shared_ptr< MatrixDouble > dataPtr
Definition: UserDataOperators.hpp:712
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl::dotVector
VectorDouble dotVector
Definition: UserDataOperators.hpp:714
convert.type
type
Definition: convert.py:64
MoFEM::PetscData::CtxSetX_T
static constexpr Switches CtxSetX_T
Definition: LoopMethods.hpp:40
MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl::zeroAtType
const EntityHandle zeroAtType
Definition: UserDataOperators.hpp:713
MoFEM::ForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: ForcesAndSourcesCore.hpp:482
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index
Definition: Index.hpp:23
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEMethod
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Definition: ForcesAndSourcesCore.hpp:1042
MoFEM::PetscData::CtxSetX_TT
static constexpr Switches CtxSetX_TT
Definition: LoopMethods.hpp:41
MoFEM::PetscData::CtxSetX
static constexpr Switches CtxSetX
Definition: LoopMethods.hpp:39
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::PetscData::CTX_SET_X_T
@ CTX_SET_X_T
Definition: LoopMethods.hpp:28
MoFEM::PetscData::CTX_SET_X
@ CTX_SET_X
Definition: LoopMethods.hpp:27
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567