v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | List of all members
MoFEM::OpDGProjectionEvaluation Struct Reference

#include "src/finite_elements/DGProjection.hpp"

Inheritance diagram for MoFEM::OpDGProjectionEvaluation:
[legend]
Collaboration diagram for MoFEM::OpDGProjectionEvaluation:
[legend]

Public Member Functions

 OpDGProjectionEvaluation (boost::shared_ptr< MatrixDouble > data_ptr, boost::shared_ptr< MatrixDouble > coeffs_ptr, boost::shared_ptr< EntitiesFieldData > data_l2, const FieldApproximationBase b, const FieldSpace s, const LogManager::SeverityLevel sev=Sev::noisy)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side.
 
- Public Member Functions inherited from MoFEM::OpDGProjectionCoefficients
 OpDGProjectionCoefficients (boost::shared_ptr< MatrixDouble > data_ptr, boost::shared_ptr< MatrixDouble > coeffs_ptr, boost::shared_ptr< MatrixDouble > mass_ptr, boost::shared_ptr< EntitiesFieldData > data_l2, const FieldApproximationBase b, const FieldSpace s, const LogManager::SeverityLevel sev=Sev::noisy)
 
- Public Member Functions inherited from MoFEM::OpBaseDerivativesBase
 OpBaseDerivativesBase (boost::shared_ptr< MatrixDouble > base_mass_ptr, boost::shared_ptr< EntitiesFieldData > data_l2, const FieldApproximationBase b, const FieldSpace s, int verb=QUIET, Sev sev=Sev::verbose)
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 Constructor for operators working on finite element spaces.
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 Constructor for operators working on a single field.
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 Constructor for operators working on two fields (bilinear forms)
 
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
 
MoFEM::InterfacegetMField ()
 
moab::Interface & getMoab ()
 
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 = 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.
 
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 Types inherited from MoFEM::OpBaseDerivativesBase
using GetOrderFun = boost::function< int()>
 
- Protected Member Functions inherited from MoFEM::OpBaseDerivativesBase
MoFEMErrorCode calculateBase (GetOrderFun get_order)
 
MoFEMErrorCode calculateMass ()
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::OpDGProjectionCoefficients
boost::shared_ptr< MatrixDoubledataPtr
 
boost::shared_ptr< MatrixDoublecoeffsPtr
 
- Protected Attributes inherited from MoFEM::OpBaseDerivativesBase
FieldApproximationBase base
 
int verbosity
 
Sev severityLevel
 
boost::shared_ptr< MatrixDoublebaseMassPtr
 
boost::shared_ptr< EntitiesFieldDatadataL2
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Examples
mofem/atom_tests/dg_projection.cpp, mofem/tutorials/adv-6/between_meshes_dg_projection.cpp, and mofem/users_modules/adolc-plasticity/adolc_plasticity.cpp.

Definition at line 166 of file DGProjection.hpp.

Constructor & Destructor Documentation

◆ OpDGProjectionEvaluation()

MoFEM::OpDGProjectionEvaluation::OpDGProjectionEvaluation ( boost::shared_ptr< MatrixDouble data_ptr,
boost::shared_ptr< MatrixDouble coeffs_ptr,
boost::shared_ptr< EntitiesFieldData data_l2,
const FieldApproximationBase  b,
const FieldSpace  s,
const LogManager::SeverityLevel  sev = Sev::noisy 
)

Definition at line 335 of file DGProjection.cpp.

341 : OpDGProjectionCoefficients(data_ptr, coeffs_ptr,
342 boost::make_shared<MatrixDouble>(), data_l2, b,
343 s, sev) {}
OpDGProjectionCoefficients(boost::shared_ptr< MatrixDouble > data_ptr, boost::shared_ptr< MatrixDouble > coeffs_ptr, boost::shared_ptr< MatrixDouble > mass_ptr, boost::shared_ptr< EntitiesFieldData > data_l2, const FieldApproximationBase b, const FieldSpace s, const LogManager::SeverityLevel sev=Sev::noisy)

Member Function Documentation

◆ doWork()

MoFEMErrorCode MoFEM::OpDGProjectionEvaluation::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)
virtual

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

Reimplemented from MoFEM::OpDGProjectionCoefficients.

Definition at line 346 of file DGProjection.cpp.

347 {
348
350
351 if (sPace != L2) {
352 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
353 "Space should be set to L2");
354 }
355
356 auto fe_type = getFEType();
357 constexpr auto side_number = 0;
358 auto order = dataL2->dataOnEntities[fe_type][side_number].getOrder();
360
361 [order]() { return order; }
362
363 );
364
365 auto &ent_data = dataL2->dataOnEntities[fe_type][side_number];
366 const auto nb_base_functions = ent_data.getN(base).size2();
367 auto nb_data_coeffs = dataPtr->size1();
368 if (!nb_data_coeffs) {
369 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
370 "Number of data coefficients should be set");
371 }
372 auto &data = *dataPtr;
373 auto &coeffs = *coeffsPtr;
374
375 auto nb_integration_pts = getGaussPts().size2();
376 data.resize(nb_integration_pts, nb_data_coeffs);
377 data.clear();
378
380
381 auto t_base = ent_data.getFTensor0N(base, 0, 0);
382 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
383 Tensor0 t_coeffs(&*coeffs.data().begin());
384 for (auto bb = 0; bb != nb_base_functions; ++bb) {
385 double alpha = t_base;
386 Tensor0 t_data(&data(gg, 0));
387 for (auto cc = 0; cc != nb_data_coeffs; ++cc) {
388 t_data += alpha * t_coeffs;
389 ++t_data;
390 ++t_coeffs;
391 }
392 ++t_base;
393 }
394 }
395
396 *dataPtr = trans(*dataPtr);
397
399}
@ L2
field with C-1 continuity
Definition definitions.h:88
#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.
constexpr int order
Order displacement.
EntityType getFEType() const
Get dimension of finite element.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
boost::shared_ptr< EntitiesFieldData > dataL2
MoFEMErrorCode calculateBase(GetOrderFun get_order)
boost::shared_ptr< MatrixDouble > coeffsPtr
boost::shared_ptr< MatrixDouble > dataPtr

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