v0.15.0
Loading...
Searching...
No Matches
PoissonExample::Op_g Struct Reference

Assemble constrains vector. More...

#include "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.
 
virtual MoFEMErrorCode iNtegrate (EntitiesFieldData::EntData &data)=0
 Class dedicated to integrate operator.
 
virtual MoFEMErrorCode aSsemble (EntitiesFieldData::EntData &data)=0
 Class dedicated to assemble operator to global system vector.
 

Protected Member Functions

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

Protected Attributes

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

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 419 of file PoissonOperators.hpp.

Member Typedef Documentation

◆ FVal

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

Definition at line 423 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 425 of file PoissonOperators.hpp.

426 : OpBaseRhs<FaceElementForcesAndSourcesCore::UserDataOperator>(
427 field_name),
428 fValue(f_value), bEta(beta) {}
constexpr auto field_name
FVal fValue
Function pointer evaluating values of "U" at the boundary.

Member Function Documentation

◆ aSsemble()

MoFEMErrorCode PoissonExample::Op_g::aSsemble ( EntitiesFieldData::EntData & data)
inlineprotected

assemble constrains vectors

Examples
PoissonOperators.hpp.

Definition at line 480 of file PoissonOperators.hpp.

480 {
482 const int *indices = &*data.getIndices().data().begin();
483 const double *vals = &*locVec.data().begin();
484 Vec f = getFEMethod()->ksp_f != PETSC_NULLPTR ? getFEMethod()->ksp_f
485 : getFEMethod()->snes_f;
486 CHKERR VecSetValues(f, nbRows, indices, &*vals, ADD_VALUES);
488 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
const VectorInt & getIndices() const
Get global indices of dofs on entity.

◆ iNtegrate()

MoFEMErrorCode PoissonExample::Op_g::iNtegrate ( EntitiesFieldData::EntData & data)
inlineprotected

Integrate local constrains vector.

Examples
PoissonOperators.hpp.

Definition at line 442 of file PoissonOperators.hpp.

442 {
444 // set size to local vector
445 locVec.resize(nbRows, false);
446 // clear local vector
447 locVec.clear();
448 // get face area
449 const double area = getArea() * bEta;
450 // get integration weight
451 auto t_w = getFTensor0IntegrationWeight();
452 // get base function
453 auto t_l = data.getFTensor0N();
454 // get coordinates at integration point
455 auto t_coords = getFTensor1CoordsAtGaussPts();
456 // make loop over integration points
457 for (int gg = 0; gg != nbIntegrationPts; gg++) {
458 // evaluate function on boundary and scale it by area and integration
459 // weight
460 double alpha =
461 area * t_w * fValue(t_coords(NX), t_coords(NY), t_coords(NZ));
462 // get element of vector
464 &*locVec.data().begin());
465 //
466 for (int rr = 0; rr != nbRows; rr++) {
467 t_a += alpha * t_l;
468 ++t_a;
469 ++t_l;
470 }
471 ++t_w;
472 ++t_coords;
473 }
475 }
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
FTensor::Number< 2 > NZ
z-direction index
FTensor::Number< 0 > NX
x-direction index
FTensor::Number< 1 > NY
y-direction index

Member Data Documentation

◆ bEta

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

Definition at line 437 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 434 of file PoissonOperators.hpp.

◆ locVec

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

Definition at line 436 of file PoissonOperators.hpp.

◆ NX

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

x-direction index

Examples
PoissonOperators.hpp.

Definition at line 431 of file PoissonOperators.hpp.

◆ NY

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

y-direction index

Examples
PoissonOperators.hpp.

Definition at line 432 of file PoissonOperators.hpp.

◆ NZ

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

z-direction index

Examples
PoissonOperators.hpp.

Definition at line 433 of file PoissonOperators.hpp.


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