v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
AtomTest::OpError< FIELD_DIM > Struct Template Reference

Operator to evaluate errors. More...

Inheritance diagram for AtomTest::OpError< FIELD_DIM >:
[legend]
Collaboration diagram for AtomTest::OpError< FIELD_DIM >:
[legend]

Public Member Functions

 OpError (boost::shared_ptr< CommonData > &common_data_ptr)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object 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)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is 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 call this function to loop over parent elements. This function calls finite element with is 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 call this function to loop over parent elements. This function calls finite element with is 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 doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
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...
 

Public Attributes

boost::shared_ptr< CommonDatacommonDataPtr
 
- 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...
 

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...
 
- 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)>
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

template<int FIELD_DIM>
struct AtomTest::OpError< FIELD_DIM >

Operator to evaluate errors.

Examples
child_and_parent.cpp, hanging_node_approx.cpp, and higher_derivatives.cpp.

Definition at line 115 of file higher_derivatives.cpp.

Constructor & Destructor Documentation

◆ OpError()

template<int FIELD_DIM>
AtomTest::OpError< FIELD_DIM >::OpError ( boost::shared_ptr< CommonData > &  common_data_ptr)
inline

Definition at line 117 of file higher_derivatives.cpp.

118 : DomainEleOp(FIELD_NAME, OPROW), commonDataPtr(common_data_ptr) {
119 std::fill(doEntities.begin(), doEntities.end(), false);
120 doEntities[MBVERTEX] = true;
121 }
constexpr char FIELD_NAME[]
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
boost::shared_ptr< CommonData > commonDataPtr
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
@ OPROW
operator doWork function is executed on FE rows

Member Function Documentation

◆ doWork()

template<int FIELD_DIM>
MoFEMErrorCode AtomTest::OpError< FIELD_DIM >::doWork ( int  side,
EntityType  type,
EntData data 
)
inlinevirtual

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

Reimplemented from MoFEM::DataOperator.

Examples
child_and_parent.cpp, hanging_node_approx.cpp, and higher_derivatives.cpp.

Definition at line 122 of file higher_derivatives.cpp.

122 {
124
125 const int nb_integration_pts = getGaussPts().size2();
126 auto t_w = getFTensor0IntegrationWeight(); // ger integration weights
127 auto t_val = getFTensor0FromVec(*(
128 commonDataPtr->approxVals)); // get function approximation at gauss pts
129 auto t_grad_val = getFTensor1FromMat<SPACE_DIM>(
131 ->approxGradVals)); // get gradient of approximation at gauss pts
132 auto t_hessian_val = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(
133 *(commonDataPtr)->approxHessianVals); // get hessian of approximation
134 // at integration pts
135
136 auto t_inv_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(
137 *(commonDataPtr->invJacPtr)); // get inverse of element jacobian
138 auto t_coords = getFTensor1CoordsAtGaussPts(); // get coordinates of
139 // integration points
140
141 // Indices used for tensor operations
146
147 const double volume = getMeasure(); // get finite element area
148
149
150 std::array<double, 3> error = {0, 0,
151 0}; // array for storing operator errors
152
153 for (int gg = 0; gg != nb_integration_pts; ++gg) {
154
155 const double alpha = t_w * volume;
156
157 // Calculate errors
158
159 double diff = t_val - fun(t_coords(0), t_coords(1), t_coords(2));
160 error[0] += alpha * pow(diff, 2);
161 auto t_diff_grad = diff_fun(t_coords(0), t_coords(1), t_coords(2));
162 t_diff_grad(i) -= t_grad_val(i);
163
164 error[1] += alpha * t_diff_grad(i) *
165 t_diff_grad(i); // note push forward derivatives
166
168 MOFEM_LOG("SELF", Sev::noisy) << "t_hessian_val " << t_hessian_push;
169
170 // hessian expected to have symmetry
171 if (std::abs(t_hessian_val(0, 1) - t_hessian_val(1, 0)) >
172 std::numeric_limits<float>::epsilon()) {
173 MOFEM_LOG("SELF", Sev::error) << "t_hessian_val " << t_hessian_val;
174 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
175 "Hessian should be symmetric");
176 }
177
178 auto t_diff_hessian = diff2_fun(t_coords(0), t_coords(1), t_coords(2));
179 t_diff_hessian(i, j) -= t_hessian_val(i, j);
180 error[2] = t_diff_hessian(i, j) * t_diff_hessian(i, j);
181
182 // move data to next integration point
183 ++t_w;
184 ++t_val;
185 ++t_grad_val;
186 ++t_hessian_val;
187 ++t_inv_jac;
188 ++t_coords;
189 }
190
191 // assemble error vector
192 std::array<int, 3> index = {0, 1, 2};
193 CHKERR VecSetValues(commonDataPtr->L2Vec, 3, index.data(), error.data(),
194 ADD_VALUES);
195
197 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
FTensor::Index< 'i', SPACE_DIM > i
auto diff2_fun
Function second derivative.
auto fun
Function to approximate.
auto diff_fun
Function derivative.
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

Member Data Documentation

◆ commonDataPtr

template<int FIELD_DIM>
boost::shared_ptr<CommonData> AtomTest::OpError< FIELD_DIM >::commonDataPtr

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