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

Calculate the grad-grad operator and assemble matrix. More...

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

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

Public Member Functions

 OpK (bool symm=true)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Do calculations for give operator. More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Protected Member Functions

virtual MoFEMErrorCode iNtegrate (EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Integrate grad-grad operator. More...
 
virtual MoFEMErrorCode aSsemble (EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Assemble local entity block matrix. More...
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 

Protected Attributes

int nbRows
 < error code More...
 
int nbCols
 number if dof on column More...
 
int nbIntegrationPts
 number of integration points More...
 
bool isDiag
 true if this block is on diagonal More...
 
FTensor::Index< 'i', 3 > i
 summit Index More...
 
MatrixDouble locMat
 local entity block matrix More...
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0, OPCOL = 1 << 1, OPROWCOL = 1 << 2, OPSPACE = 1 << 3,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > >>
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 

Detailed Description

Calculate the grad-grad operator and assemble matrix.

Calculate

\[ \mathbf{K}=\int_\Omega \nabla \boldsymbol\phi \cdot \nabla \boldsymbol\phi \textrm{d}\Omega \]

and assemble to global matrix.

This operator is executed on element for each unique combination of entities.

Definition at line 28 of file PoissonOperators.hpp.

Constructor & Destructor Documentation

◆ OpK()

PoissonExample::OpK::OpK ( bool  symm = true)
inline

Definition at line 30 of file PoissonOperators.hpp.

Member Function Documentation

◆ aSsemble()

virtual MoFEMErrorCode PoissonExample::OpK::aSsemble ( EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)
inlineprotectedvirtual

Assemble local entity block matrix.

Parameters
row_datarow data (consist base functions on row entity)
col_datacolumn data (consist base functions on column entity)
Returns
error code

Definition at line 137 of file PoissonOperators.hpp.

138  {
140  // get pointer to first global index on row
141  const int *row_indices = &*row_data.getIndices().data().begin();
142  // get pointer to first global index on column
143  const int *col_indices = &*col_data.getIndices().data().begin();
144  Mat B = getFEMethod()->ksp_B != PETSC_NULL ? getFEMethod()->ksp_B
145  : getFEMethod()->snes_B;
146  // assemble local matrix
147  CHKERR MatSetValues(B, nbRows, row_indices, nbCols, col_indices,
148  &*locMat.data().begin(), ADD_VALUES);
149 
150  if (!isDiag && sYmm) {
151  // if not diagonal term and since global matrix is symmetric assemble
152  // transpose term.
153  locMat = trans(locMat);
154  CHKERR MatSetValues(B, nbCols, col_indices, nbRows, row_indices,
155  &*locMat.data().begin(), ADD_VALUES);
156  }
158  }

◆ doWork()

MoFEMErrorCode PoissonExample::OpK::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)
inlinevirtual

Do calculations for give operator.

Parameters
row_siderow side number (local number) of entity on element
col_sidecolumn side number (local number) of entity on element
row_typetype of row entity MBVERTEX, MBEDGE, MBTRI or MBTET
col_typetype of column entity MBVERTEX, MBEDGE, MBTRI or MBTET
row_datadata for row
col_datadata for column
Returns
error code

Reimplemented from MoFEM::DataOperator.

Definition at line 44 of file PoissonOperators.hpp.

47  {
49  // get number of dofs on row
50  nbRows = row_data.getIndices().size();
51  // if no dofs on row, exit that work, nothing to do here
52  if (!nbRows)
54  // get number of dofs on column
55  nbCols = col_data.getIndices().size();
56  // if no dofs on Columbia, exit nothing to do here
57  if (!nbCols)
59  // get number of integration points
60  nbIntegrationPts = getGaussPts().size2();
61  // check if entity block is on matrix diagonal
62  if (row_side == col_side && row_type == col_type) {
63  isDiag = true; // yes, it is
64  } else {
65  isDiag = false;
66  }
67  // integrate local matrix for entity block
68  CHKERR iNtegrate(row_data, col_data);
69  // assemble local matrix
70  CHKERR aSsemble(row_data, col_data);
72  }

◆ iNtegrate()

virtual MoFEMErrorCode PoissonExample::OpK::iNtegrate ( EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)
inlineprotectedvirtual

Integrate grad-grad operator.

Parameters
row_datarow data (consist base functions on row entity)
col_datacolumn data (consist base functions on column entity)
Returns
error code

Reimplemented in PoissonExample::OpKt.

Definition at line 91 of file PoissonOperators.hpp.

92  {
94  // set size of local entity bock
95  locMat.resize(nbRows, nbCols, false);
96  // clear matrix
97  locMat.clear();
98  // get element volume
99  double vol = getVolume();
100  // get integration weights
101  auto t_w = getFTensor0IntegrationWeight();
102  // get base function gradient on rows
103  auto t_row_grad = row_data.getFTensor1DiffN<3>();
104  // loop over integration points
105  for (int gg = 0; gg != nbIntegrationPts; gg++) {
106  // take into account Jacobean
107  const double alpha = t_w * vol;
108  // noalias(locMat) +=
109  // alpha*prod(row_data.getDiffN(gg),trans(col_data.getDiffN(gg))); take
110  // fist element to local matrix
111  FTensor::Tensor0<double *> a(&*locMat.data().begin());
112  // loop over rows base functions
113  for (int rr = 0; rr != nbRows; rr++) {
114  // get column base functions gradient at gauss point gg
115  auto t_col_grad = col_data.getFTensor1DiffN<3>(gg, 0);
116  // loop over columns
117  for (int cc = 0; cc != nbCols; cc++) {
118  // calculate element of local matrix
119  a += alpha * (t_row_grad(i) * t_col_grad(i));
120  ++t_col_grad; // move to another gradient of base function on column
121  ++a; // move to another element of local matrix in column
122  }
123  ++t_row_grad; // move to another element of gradient of base function on
124  // row
125  }
126  ++t_w; // move to another integration weight
127  }
129  }

Member Data Documentation

◆ i

FTensor::Index<'i', 3> PoissonExample::OpK::i
protected

summit Index

Definition at line 82 of file PoissonOperators.hpp.

◆ isDiag

bool PoissonExample::OpK::isDiag
protected

true if this block is on diagonal

Definition at line 80 of file PoissonOperators.hpp.

◆ locMat

MatrixDouble PoissonExample::OpK::locMat
protected

local entity block matrix

Definition at line 83 of file PoissonOperators.hpp.

◆ nbCols

int PoissonExample::OpK::nbCols
protected

number if dof on column

Definition at line 78 of file PoissonOperators.hpp.

◆ nbIntegrationPts

int PoissonExample::OpK::nbIntegrationPts
protected

number of integration points

Definition at line 79 of file PoissonOperators.hpp.

◆ nbRows

int PoissonExample::OpK::nbRows
protected

< error code

number of dofs on rows

Definition at line 77 of file PoissonOperators.hpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::MatSetValues
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
Definition: EntitiesFieldData.hpp:1644
MoFEM::KspMethod::ksp_B
Mat & ksp_B
Definition: LoopMethods.hpp:91
MoFEM::SnesMethod::snes_B
Mat & snes_B
preconditioner of jacobian matrix
Definition: LoopMethods.hpp:124
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getVolume
double getVolume() const
element volume (linear geometry)
Definition: VolumeElementForcesAndSourcesCore.hpp:161
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeight
auto getFTensor0IntegrationWeight()
Get integration weights.
Definition: ForcesAndSourcesCore.hpp:1240
MoFEM::ForcesAndSourcesCore::UserDataOperator::getGaussPts
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Definition: ForcesAndSourcesCore.hpp:1236
PoissonExample::OpK::nbIntegrationPts
int nbIntegrationPts
number of integration points
Definition: PoissonOperators.hpp:79
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
a
constexpr double a
Definition: approx_sphere.cpp:30
PoissonExample::OpK::locMat
MatrixDouble locMat
local entity block matrix
Definition: PoissonOperators.hpp:83
PoissonExample::OpK::nbCols
int nbCols
number if dof on column
Definition: PoissonOperators.hpp:78
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1214
PoissonExample::OpK::i
FTensor::Index< 'i', 3 > i
summit Index
Definition: PoissonOperators.hpp:82
PoissonExample::OpK::nbRows
int nbRows
< error code
Definition: PoissonOperators.hpp:77
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEMethod
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Definition: ForcesAndSourcesCore.hpp:1042
PoissonExample::OpK::aSsemble
virtual MoFEMErrorCode aSsemble(EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Assemble local entity block matrix.
Definition: PoissonOperators.hpp:137
MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
Definition: EntitiesFieldData.cpp:526
FTensor::Tensor0
Definition: Tensor0.hpp:16
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
PoissonExample::OpK::isDiag
bool isDiag
true if this block is on diagonal
Definition: PoissonOperators.hpp:80
PoissonExample::OpK::iNtegrate
virtual MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Integrate grad-grad operator.
Definition: PoissonOperators.hpp:91
MoFEM::DataOperator::sYmm
bool sYmm
If true assume that matrix is symmetric structure.
Definition: DataOperators.hpp:82
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359