v0.14.0
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
PoissonExample::Op_g Struct Reference

Assemble constrains vector. More...

#include <users_modules/tutorials/cor-2to5/src/PoissonOperators.hpp>

Inheritance diagram for PoissonExample::Op_g:
[legend]
Collaboration diagram for PoissonExample::Op_g:
[legend]

Public Types

typedef boost::function< double(const double, const double, const double)> FVal
 

Public Member Functions

 Op_g (FVal f_value, const string field_name="L", const double beta=1)
 
- Public Member Functions inherited from PoissonExample::OpBaseRhs< FaceElementForcesAndSourcesCore::UserDataOperator >
 OpBaseRhs (const std::string field_name)
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
 This function is called by finite element. More...
 

Protected Member Functions

MoFEMErrorCode iNtegrate (EntitiesFieldData::EntData &data)
 Integrate local constrains vector. More...
 
MoFEMErrorCode aSsemble (EntitiesFieldData::EntData &data)
 assemble constrains vectors More...
 

Protected Attributes

FTensor::Number< 0 > NX
 x-direction index More...
 
FTensor::Number< 1 > NY
 y-direction index More...
 
FTensor::Number< 2 > NZ
 z-direction index More...
 
FVal fValue
 Function pointer evaluating values of "U" at the boundary. More...
 
VectorDouble locVec
 
const double bEta
 
- Protected Attributes inherited from PoissonExample::OpBaseRhs< FaceElementForcesAndSourcesCore::UserDataOperator >
int nbRows
 < error code More...
 
int nbIntegrationPts
 number of integration points More...
 

Detailed Description

Assemble constrains vector.

\[ \mathbf{g} = \int_{\partial\Omega} \boldsymbol\psi \overline{u} \textrm{d}\partial\Omega \]

Examples
analytical_poisson_field_split.cpp.

Definition at line 416 of file PoissonOperators.hpp.

Member Typedef Documentation

◆ FVal

typedef boost::function<double(const double, const double, const double)> PoissonExample::Op_g::FVal
Examples
PoissonOperators.hpp.

Definition at line 420 of file PoissonOperators.hpp.

Constructor & Destructor Documentation

◆ Op_g()

PoissonExample::Op_g::Op_g ( FVal  f_value,
const string  field_name = "L",
const double  beta = 1 
)
inline
Examples
PoissonOperators.hpp.

Definition at line 422 of file PoissonOperators.hpp.

423  : OpBaseRhs<FaceElementForcesAndSourcesCore::UserDataOperator>(
424  field_name),
425  fValue(f_value), bEta(beta) {}

Member Function Documentation

◆ aSsemble()

MoFEMErrorCode PoissonExample::Op_g::aSsemble ( EntitiesFieldData::EntData data)
inlineprotectedvirtual

assemble constrains vectors

Implements PoissonExample::OpBaseRhs< FaceElementForcesAndSourcesCore::UserDataOperator >.

Examples
PoissonOperators.hpp.

Definition at line 477 of file PoissonOperators.hpp.

477  {
479  const int *indices = &*data.getIndices().data().begin();
480  const double *vals = &*locVec.data().begin();
481  Vec f = getFEMethod()->ksp_f != PETSC_NULL ? getFEMethod()->ksp_f
482  : getFEMethod()->snes_f;
483  CHKERR VecSetValues(f, nbRows, indices, &*vals, ADD_VALUES);
485  }

◆ iNtegrate()

MoFEMErrorCode PoissonExample::Op_g::iNtegrate ( EntitiesFieldData::EntData data)
inlineprotectedvirtual

Integrate local constrains vector.

Implements PoissonExample::OpBaseRhs< FaceElementForcesAndSourcesCore::UserDataOperator >.

Reimplemented in PoissonExample::OpResF_Boundary, and PoissonExample::OpRes_g.

Examples
PoissonOperators.hpp.

Definition at line 439 of file PoissonOperators.hpp.

439  {
441  // set size to local vector
442  locVec.resize(nbRows, false);
443  // clear local vector
444  locVec.clear();
445  // get face area
446  const double area = getArea() * bEta;
447  // get integration weight
448  auto t_w = getFTensor0IntegrationWeight();
449  // get base function
450  auto t_l = data.getFTensor0N();
451  // get coordinates at integration point
452  auto t_coords = getFTensor1CoordsAtGaussPts();
453  // make loop over integration points
454  for (int gg = 0; gg != nbIntegrationPts; gg++) {
455  // evaluate function on boundary and scale it by area and integration
456  // weight
457  double alpha =
458  area * t_w * fValue(t_coords(NX), t_coords(NY), t_coords(NZ));
459  // get element of vector
461  &*locVec.data().begin());
462  //
463  for (int rr = 0; rr != nbRows; rr++) {
464  t_a += alpha * t_l;
465  ++t_a;
466  ++t_l;
467  }
468  ++t_w;
469  ++t_coords;
470  }
472  }

Member Data Documentation

◆ bEta

const double PoissonExample::Op_g::bEta
protected
Examples
PoissonOperators.hpp.

Definition at line 434 of file PoissonOperators.hpp.

◆ fValue

FVal PoissonExample::Op_g::fValue
protected

Function pointer evaluating values of "U" at the boundary.

Examples
PoissonOperators.hpp.

Definition at line 431 of file PoissonOperators.hpp.

◆ locVec

VectorDouble PoissonExample::Op_g::locVec
protected
Examples
PoissonOperators.hpp.

Definition at line 433 of file PoissonOperators.hpp.

◆ NX

FTensor::Number<0> PoissonExample::Op_g::NX
protected

x-direction index

Examples
PoissonOperators.hpp.

Definition at line 428 of file PoissonOperators.hpp.

◆ NY

FTensor::Number<1> PoissonExample::Op_g::NY
protected

y-direction index

Examples
PoissonOperators.hpp.

Definition at line 429 of file PoissonOperators.hpp.

◆ NZ

FTensor::Number<2> PoissonExample::Op_g::NZ
protected

z-direction index

Examples
PoissonOperators.hpp.

Definition at line 430 of file PoissonOperators.hpp.


The documentation for this struct was generated from the following file:
PoissonExample::Op_g::NY
FTensor::Number< 1 > NY
y-direction index
Definition: PoissonOperators.hpp:429
PoissonExample::Op_g::NX
FTensor::Number< 0 > NX
x-direction index
Definition: PoissonOperators.hpp:428
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1576
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
PoissonExample::Op_g::bEta
const double bEta
Definition: PoissonOperators.hpp:434
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
PoissonExample::Op_g::locVec
VectorDouble locVec
Definition: PoissonOperators.hpp:433
PoissonExample::Op_g::fValue
FVal fValue
Function pointer evaluating values of "U" at the boundary.
Definition: PoissonOperators.hpp:431
FTensor::Tensor0
Definition: Tensor0.hpp:16
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
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
PoissonExample::OpBaseRhs< FaceElementForcesAndSourcesCore::UserDataOperator >::nbRows
int nbRows
< error code
Definition: PoissonOperators.hpp:204
PoissonExample::OpBaseRhs< FaceElementForcesAndSourcesCore::UserDataOperator >::nbIntegrationPts
int nbIntegrationPts
number of integration points
Definition: PoissonOperators.hpp:205
PoissonExample::Op_g::NZ
FTensor::Number< 2 > NZ
z-direction index
Definition: PoissonOperators.hpp:430