v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
MinimalSurfaceEquation::MinimalSurfaceElement::OpGetDataAtGaussPts Struct Reference

Evaluate function values and gradients at Gauss Pts. More...

#include <users_modules/minimal_surface_equation/src/MinimalSurfaceElement.hpp>

Inheritance diagram for MinimalSurfaceEquation::MinimalSurfaceElement::OpGetDataAtGaussPts:
[legend]
Collaboration diagram for MinimalSurfaceEquation::MinimalSurfaceElement::OpGetDataAtGaussPts:
[legend]

Public Member Functions

 OpGetDataAtGaussPts (const string field_name, CommonData &common_data)
 
PetscErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 

Public Attributes

CommonDatacommonData
 

Detailed Description

Evaluate function values and gradients at Gauss Pts.

Examples
minimal_surface_area.cpp.

Definition at line 190 of file MinimalSurfaceElement.hpp.

Constructor & Destructor Documentation

◆ OpGetDataAtGaussPts()

MinimalSurfaceEquation::MinimalSurfaceElement::OpGetDataAtGaussPts::OpGetDataAtGaussPts ( const string  field_name,
CommonData common_data 
)
inline

Definition at line 193 of file MinimalSurfaceElement.hpp.

194 : FaceElementForcesAndSourcesCore::UserDataOperator(
195 field_name, UserDataOperator::OPCOL),
196 commonData(common_data) {}
constexpr auto field_name

Member Function Documentation

◆ doWork()

PetscErrorCode MinimalSurfaceEquation::MinimalSurfaceElement::OpGetDataAtGaussPts::doWork ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData &  data 
)
inline

Definition at line 198 of file MinimalSurfaceElement.hpp.

199 {
201
202 int nb_dofs = data.getFieldData().size();
203 if (nb_dofs == 0) {
205 }
206 int nb_gauss_pts = data.getN().size1();
207
208 // Element loops over entities, it start from vertices
209 if (type == MBVERTEX) {
210 // Setting size of common data vectors
211 commonData.gradU.resize(nb_gauss_pts, 2);
212 commonData.gradU.clear();
213 }
214
215 for (int gg = 0; gg != nb_gauss_pts; gg++) {
216 MatrixAdaptor grad_shape_fun = data.getDiffN(gg);
217 ublas::matrix_row<MatrixDouble> grad_at_gauss_pt(commonData.gradU, gg);
218 noalias(grad_at_gauss_pt) +=
219 prod(trans(grad_shape_fun), data.getFieldData());
220 }
221
223 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
Definition: Types.hpp:132
MatrixDouble gradU
nb_gauss_pts x 2 gradients at integration pts

Member Data Documentation

◆ commonData

CommonData& MinimalSurfaceEquation::MinimalSurfaceElement::OpGetDataAtGaussPts::commonData

Definition at line 192 of file MinimalSurfaceElement.hpp.


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