v0.8.23
Public Member Functions | Public Attributes | List of all members
MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcRhs Struct Reference

Integrate vector on rhs,. More...

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

Inherits UserDataOperator.

Public Member Functions

 OpAssmebleBcRhs (const string field_name, Vec v_f)
 
virtual double evalFunction (double x, double y)
 Function on boundary Inherit this class and overload this function to change bc. More...
 
PetscErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 

Public Attributes

Vec vF
 
VectorDouble nF
 

Detailed Description

Integrate vector on rhs,.

\[ \mathbf{F} = \int_L \mathbf{N}^\textrm{T} g(x,y) \textrm{d}L \]

Examples
minimal_surface_area.cpp.

Definition at line 92 of file MinimalSurfaceElement.hpp.

Constructor & Destructor Documentation

◆ OpAssmebleBcRhs()

MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcRhs::OpAssmebleBcRhs ( const string  field_name,
Vec  v_f 
)

Definition at line 95 of file MinimalSurfaceElement.hpp.

97  field_name, UserDataOperator::OPROW),
98  vF(v_f) {}
ForcesAndSourcesCore::UserDataOperator UserDataOperator

Member Function Documentation

◆ doWork()

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

Definition at line 109 of file MinimalSurfaceElement.hpp.

110  {
112 
113  int nb_dofs = data.getFieldData().size();
114  if (nb_dofs == 0) {
116  }
117  nF.resize(nb_dofs, false);
118  nF.clear();
119  int nb_gauss_pts = data.getN().size1();
120  for (int gg = 0; gg != nb_gauss_pts; gg++) {
121  double val = getLength() * getGaussPts()(1, gg);
122  double x = getCoordsAtGaussPts()(gg, 0);
123  double y = getCoordsAtGaussPts()(gg, 1);
124  noalias(nF) += val * evalFunction(x, y) * data.getN(gg);
125  }
126  CHKERR VecSetValues(vF, data.getIndices().size(),
127  &*data.getIndices().data().begin(), &nF[0],
128  ADD_VALUES);
130  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
MoFEMErrorCode VecSetValues(Vec V, const DataForcesAndSourcesCore::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual double evalFunction(double x, double y)
Function on boundary Inherit this class and overload this function to change bc.

◆ evalFunction()

virtual double MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcRhs::evalFunction ( double  x,
double  y 
)
virtual

Function on boundary Inherit this class and overload this function to change bc.

Definition at line 103 of file MinimalSurfaceElement.hpp.

103  {
104  return sin(2 * M_PI * (x + y));
105  // return sin(M_PI*x);
106  }

Member Data Documentation

◆ nF

VectorDouble MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcRhs::nF

Definition at line 108 of file MinimalSurfaceElement.hpp.

◆ vF

Vec MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcRhs::vF

Definition at line 94 of file MinimalSurfaceElement.hpp.


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