v0.14.0
Public Member Functions | Public Attributes | List of all members
MinimalSurfaceEquation::MinimalSurfaceElement::OpCalculateCoefficientsAtGaussPts 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::OpCalculateCoefficientsAtGaussPts:
[legend]
Collaboration diagram for MinimalSurfaceEquation::MinimalSurfaceElement::OpCalculateCoefficientsAtGaussPts:
[legend]

Public Member Functions

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

Public Attributes

CommonDatacommonData
 
bool rhsAndNotLhs
 

Detailed Description

Evaluate function values and gradients at Gauss Pts.

Examples
minimal_surface_area.cpp.

Definition at line 228 of file MinimalSurfaceElement.hpp.

Constructor & Destructor Documentation

◆ OpCalculateCoefficientsAtGaussPts()

MinimalSurfaceEquation::MinimalSurfaceElement::OpCalculateCoefficientsAtGaussPts::OpCalculateCoefficientsAtGaussPts ( const string  field_name,
CommonData common_data,
bool  rhs_and_not_lhs 
)
inline

Definition at line 232 of file MinimalSurfaceElement.hpp.

Member Function Documentation

◆ doWork()

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

Definition at line 239 of file MinimalSurfaceElement.hpp.

240  {
242 
243  // Element loops over entities, it start from vertices
244  if (type == MBVERTEX) {
245  int nb_gauss_pts = data.getN().size1();
246  commonData.normGradU2.resize(nb_gauss_pts, false);
247  for (int gg = 0; gg != nb_gauss_pts; gg++) {
248  ublas::matrix_row<MatrixDouble> grad_at_gauss_pt(commonData.gradU,
249  gg);
250  commonData.normGradU2[gg] =
251  inner_prod(grad_at_gauss_pt, grad_at_gauss_pt);
252  }
253  commonData.aN.resize(nb_gauss_pts, false);
254  for (int gg = 0; gg != nb_gauss_pts; gg++) {
255  commonData.aN[gg] = 1. / sqrt(1 + commonData.normGradU2[gg]);
256  }
257  if (rhsAndNotLhs) {
258  // needed at right hand side when residual is calculated
259  commonData.aNbyGradU.resize(nb_gauss_pts, 2, false);
260  for (int gg = 0; gg != nb_gauss_pts; gg++) {
261  ublas::matrix_row<MatrixDouble> grad_at_gauss_pt(commonData.gradU,
262  gg);
263  for (int rr = 0; rr != 2; rr++) {
264  commonData.aNbyGradU(gg, rr) =
265  commonData.aN[gg] * grad_at_gauss_pt[rr];
266  }
267  }
268  } else {
269  // needed at left hand side when matrix is calculated
270  commonData.aNpow3byGradU.resize(nb_gauss_pts, 2, false);
271  for (int gg = 0; gg != nb_gauss_pts; gg++) {
272  ublas::matrix_row<MatrixDouble> grad_at_gauss_pt(commonData.gradU,
273  gg);
274  for (int rr = 0; rr != 2; rr++) {
275  commonData.aNpow3byGradU(gg, rr) =
276  pow(commonData.aN[gg], 3) * grad_at_gauss_pt[rr];
277  }
278  }
279  }
280  }
281 
282  // doVerticesRow = true;
283  // doEdgesRow = false;
284  // doQuadsRow = false;
285  // doTrisRow = false;
286  // doTetsRow = false;
287  // doPrismsRow = false;
288 
290  }

Member Data Documentation

◆ commonData

CommonData& MinimalSurfaceEquation::MinimalSurfaceElement::OpCalculateCoefficientsAtGaussPts::commonData

Definition at line 230 of file MinimalSurfaceElement.hpp.

◆ rhsAndNotLhs

bool MinimalSurfaceEquation::MinimalSurfaceElement::OpCalculateCoefficientsAtGaussPts::rhsAndNotLhs

Definition at line 231 of file MinimalSurfaceElement.hpp.


The documentation for this struct was generated from the following file:
MinimalSurfaceEquation::MinimalSurfaceElement::OpCalculateCoefficientsAtGaussPts::commonData
CommonData & commonData
Definition: MinimalSurfaceElement.hpp:230
MinimalSurfaceEquation::MinimalSurfaceElement::CommonData::gradU
MatrixDouble gradU
nb_gauss_pts x 2 gradients at integration pts
Definition: MinimalSurfaceElement.hpp:76
MinimalSurfaceEquation::MinimalSurfaceElement::CommonData::aNpow3byGradU
MatrixDouble aNpow3byGradU
nb_gauss_pts x 2,
Definition: MinimalSurfaceElement.hpp:81
MinimalSurfaceEquation::MinimalSurfaceElement::CommonData::normGradU2
VectorDouble normGradU2
size of nb_gauss_pts, norm of gradient
Definition: MinimalSurfaceElement.hpp:77
MinimalSurfaceEquation::MinimalSurfaceElement::OpCalculateCoefficientsAtGaussPts::rhsAndNotLhs
bool rhsAndNotLhs
Definition: MinimalSurfaceElement.hpp:231
convert.type
type
Definition: convert.py:64
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPCOL
@ OPCOL
operator doWork function is executed on FE columns
Definition: ForcesAndSourcesCore.hpp:568
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MinimalSurfaceEquation::MinimalSurfaceElement::CommonData::aNbyGradU
MatrixDouble aNbyGradU
nb_gauss_pts x 2,
Definition: MinimalSurfaceElement.hpp:80
MinimalSurfaceEquation::MinimalSurfaceElement::CommonData::aN
VectorDouble aN
size of nb_gauss_pts,
Definition: MinimalSurfaceElement.hpp:79
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346