#include <users_modules/tutorials/scl-1/src/poisson_2d_homogeneous.hpp>
|
| OpDomainRhsVectorF (std::string field_name) |
|
MoFEMErrorCode | iNtegrate (EntitiesFieldData::EntData &data) |
|
| OpBaseImpl (const std::string row_field_name, const std::string col_field_name, const OpType type, boost::shared_ptr< Range > ents_ptr=nullptr) |
|
MoFEMErrorCode | doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data) |
| Do calculations for the left hand side. More...
|
|
MoFEMErrorCode | doWork (int row_side, EntityType row_type, EntData &row_data) |
| Do calculations for the right hand side. More...
|
|
|
using | OpType = typename EleOp::OpType |
|
using | EntData = EntitiesFieldData::EntData |
|
using | MatSetValuesHook = boost::function< MoFEMErrorCode(ForcesAndSourcesCore::UserDataOperator *op_ptr, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, MatrixDouble &m)> |
|
TimeFun | timeScalingFun |
| assumes that time variable is set More...
|
|
FEFun | feScalingFun |
| assumes that time variable is set More...
|
|
boost::shared_ptr< Range > | entsPtr |
| Entities on which element is run. More...
|
|
static MatSetValuesHook | matSetValuesHook |
|
template<int DIM> |
FTensor::Tensor1< FTensor::PackPtr< double *, DIM >, DIM > | getNf () |
|
template<int DIM> |
FTensor::Tensor2< FTensor::PackPtr< double *, DIM >, DIM, DIM > | getLocMat (const int rr) |
|
virtual MoFEMErrorCode | iNtegrate (EntData &row_data, EntData &col_data) |
| Integrate grad-grad operator. More...
|
|
virtual MoFEMErrorCode | aSsemble (EntData &row_data, EntData &col_data, const bool trans) |
|
virtual MoFEMErrorCode | iNtegrate (EntData &data) |
| Class dedicated to integrate operator. More...
|
|
virtual MoFEMErrorCode | aSsemble (EntData &data) |
|
virtual size_t | getNbOfBaseFunctions (EntitiesFieldData::EntData &data) |
| Get number of base functions. More...
|
|
int | nbRows |
| number of dofs on rows More...
|
|
int | nbCols |
| number if dof on column More...
|
|
int | nbIntegrationPts |
| number of integration points More...
|
|
int | nbRowBaseFunctions |
| number or row base functions More...
|
|
int | rowSide |
| row side number More...
|
|
int | colSide |
| column side number More...
|
|
EntityType | rowType |
| row type More...
|
|
EntityType | colType |
| column type More...
|
|
bool | assembleTranspose |
|
bool | onlyTranspose |
|
MatrixDouble | locMat |
| local entity block matrix More...
|
|
MatrixDouble | locMatTranspose |
| local entity block matrix More...
|
|
VectorDouble | locF |
| local entity vector More...
|
|
◆ OpDomainRhsVectorF()
Poisson2DHomogeneousOperators::OpDomainRhsVectorF::OpDomainRhsVectorF |
( |
std::string |
field_name | ) |
|
|
inline |
Definition at line 94 of file poisson_2d_homogeneous.hpp.
constexpr auto field_name
@ OPROW
operator doWork function is executed on FE rows
FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase AssemblyDomainEleOp
◆ iNtegrate()
MoFEMErrorCode Poisson2DHomogeneousOperators::OpDomainRhsVectorF::iNtegrate |
( |
EntitiesFieldData::EntData & |
data | ) |
|
|
inline |
- Examples
- poisson_2d_homogeneous.hpp.
Definition at line 97 of file poisson_2d_homogeneous.hpp.
97 {
99
100 const int nb_dofs = data.getIndices().size();
101
102 this->
locF.resize(nb_dofs,
false);
104
105
106 const double area = getMeasure();
107
108
109 const int nb_integration_points = getGaussPts().size2();
110
111 auto t_w = getFTensor0IntegrationWeight();
112
113
114 auto t_base = data.getFTensor0N();
115
116
117 for (int gg = 0; gg != nb_integration_points; gg++) {
118 const double a = t_w * area;
119
120 for (int rr = 0; rr != nb_dofs; rr++) {
122
123
124 ++t_base;
125 }
126
127
128 ++t_w;
129 }
130
132 }
#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()
VectorDouble locF
local entity vector
The documentation for this struct was generated from the following file: