v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
Poisson2DiscontGalerkinOperators::OpL2LhsPenalty Struct Reference

Operator the left hand side matrix. More...

#include <users_modules/tutorials/scl-11/src/PoissonDiscontinousGalerkin.hpp>

Inheritance diagram for Poisson2DiscontGalerkinOperators::OpL2LhsPenalty:
[legend]
Collaboration diagram for Poisson2DiscontGalerkinOperators::OpL2LhsPenalty:
[legend]

Public Member Functions

 OpL2LhsPenalty (boost::shared_ptr< FaceSideEle > side_fe_ptr)
 Construct a new OpL2LhsPenalty. More...
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 
- Public Member Functions inherited from MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
double getArea ()
 get area of face More...
 
double getMeasure ()
 get measure of element More...
 
VectorDoublegetNormal ()
 get triangle normal More...
 
VectorDoublegetTangent1 ()
 get triangle tangent 1 More...
 
VectorDoublegetTangent2 ()
 get triangle tangent 2 More...
 
auto getFTensor1Normal ()
 get normal as tensor More...
 
auto getFTensor1Tangent1 ()
 get tangentOne as tensor More...
 
auto getFTensor1Tangent2 ()
 get tangentTwo as tensor More...
 
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
VectorDoublegetCoords ()
 get triangle coordinates More...
 
auto getFTensor1Coords ()
 get get coords at gauss points More...
 
MatrixDoublegetNormalsAtGaussPts ()
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
MatrixDoublegetTangent1AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
auto getFTensor1NormalsAtGaussPts ()
 get normal at integration points More...
 
auto getFTensor1Tangent1AtGaussPts ()
 get tangent 1 at integration points More...
 
auto getFTensor1Tangent2AtGaussPts ()
 get tangent 2 at integration points More...
 
FaceElementForcesAndSourcesCoregetFaceFE ()
 return pointer to Generic Triangle Finite Element object More...
 
MoFEMErrorCode loopSideVolumes (const string fe_name, VolumeElementForcesAndSourcesCoreOnSide &fe_method)
 
- 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)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is 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 call this function to loop over parent elements. This function calls finite element with is 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 call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
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...
 

Private Attributes

boost::shared_ptr< FaceSideElesideFEPtr
 pointer to element to get data on edge/face sides More...
 
MatrixDouble locMat
 local operator matrix More...
 

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...
 
- 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 []
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Operator the left hand side matrix.

Examples
poisson_2d_dis_galerkin.cpp.

Definition at line 115 of file PoissonDiscontinousGalerkin.hpp.

Constructor & Destructor Documentation

◆ OpL2LhsPenalty()

Poisson2DiscontGalerkinOperators::OpL2LhsPenalty::OpL2LhsPenalty ( boost::shared_ptr< FaceSideEle side_fe_ptr)
inline

Construct a new OpL2LhsPenalty.

Parameters
side_fe_ptrpointer to FE to evaluate side elements

Definition at line 123 of file PoissonDiscontinousGalerkin.hpp.

BoundaryEle::UserDataOperator BoundaryEleOp
@ NOSPACE
Definition: definitions.h:83
@ OPSPACE
operator do Work is execute on space data
boost::shared_ptr< FaceSideEle > sideFEPtr
pointer to element to get data on edge/face sides

Member Function Documentation

◆ doWork()

MoFEMErrorCode Poisson2DiscontGalerkinOperators::OpL2LhsPenalty::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData &  data 
)
inline
Examples
PoissonDiscontinousGalerkin.hpp.

Definition at line 126 of file PoissonDiscontinousGalerkin.hpp.

127 {
129
130 // Collect data from side domian elements
131 CHKERR loopSideFaces("dFE", *sideFEPtr);
132 const auto in_the_loop =
133 sideFEPtr->nInTheLoop; // return number of elements on the side
134
135#ifndef NDEBUG
136 const std::array<std::string, 2> ele_type_name = {"BOUNDARY", "SKELETON"};
137 MOFEM_LOG("SELF", Sev::noisy)
138 << "OpL2LhsPenalty inTheLoop " << ele_type_name[in_the_loop];
139 MOFEM_LOG("SELF", Sev::noisy)
140 << "OpL2LhsPenalty sense " << senseMap[0] << " " << senseMap[1];
141#endif
142
143 // calculate penalty
144 const double s = getMeasure() / (areaMap[0] + areaMap[1]);
145 const double p = penalty * s;
146
147 // get normal of the face or edge
148 auto t_normal = getFTensor1Normal();
149 t_normal.normalize();
150
151 // get number of integration points on face
152 const size_t nb_integration_pts = getGaussPts().size2();
153
154 // beta paramter is zero, when penalty method is used, also, takes value 1,
155 // when boundary edge/face is evaluated, and 2 when is skeleton edge/face.
156 const double beta = static_cast<double>(nitsche) / (in_the_loop + 1);
157
158 // iterate the sides rows
159 for (auto s0 : {LEFT_SIDE, RIGHT_SIDE}) {
160
161 // gent number of DOFs on the right side.
162 const auto nb_rows = indicesRowSideMap[s0].size();
163
164 if (nb_rows) {
165
166 // get orientation of the local element edge
167 const auto sense_row = senseMap[s0];
168
169 // iterate the side cols
170 const auto nb_row_base_functions = rowBaseSideMap[s0].size2();
171 for (auto s1 : {LEFT_SIDE, RIGHT_SIDE}) {
172
173 // get orientation of the edge
174 const auto sense_col = senseMap[s1];
175
176 // number of dofs, for homogenous approximation this value is
177 // constant.
178 const auto nb_cols = indicesColSideMap[s1].size();
179
180 // resize local element matrix
181 locMat.resize(nb_rows, nb_cols, false);
182 locMat.clear();
183
184 // get base functions, and integration weights
185 auto t_row_base = get_ntensor(rowBaseSideMap[s0]);
186 auto t_diff_row_base = get_diff_ntensor(rowDiffBaseSideMap[s0]);
187 auto t_w = getFTensor0IntegrationWeight();
188
189 // iterate integration points on face/edge
190 for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
191
192 // t_w is integration weight, and measure is element area, or
193 // volume, depending if problem is in 2d/3d.
194 const double alpha = getMeasure() * t_w;
195 auto t_mat = locMat.data().begin();
196
197 // iterate rows
198 size_t rr = 0;
199 for (; rr != nb_rows; ++rr) {
200
201 // calculate tetting function
203 t_vn_plus(i) = beta * (phi * t_diff_row_base(i) / p);
205 t_vn(i) = t_row_base * t_normal(i) * sense_row - t_vn_plus(i);
206
207 // get base functions on columns
208 auto t_col_base = get_ntensor(colBaseSideMap[s1], gg, 0);
209 auto t_diff_col_base =
211
212 // iterate columns
213 for (size_t cc = 0; cc != nb_cols; ++cc) {
214
215 // calculate variance of tested function
217 t_un(i) = -p * (t_col_base * t_normal(i) * sense_col -
218 beta * t_diff_col_base(i) / p);
219
220 // assemble matrix
221 *t_mat -= alpha * (t_vn(i) * t_un(i));
222 *t_mat -= alpha * (t_vn_plus(i) * (beta * t_diff_col_base(i)));
223
224 // move to next column base and element of matrix
225 ++t_col_base;
226 ++t_diff_col_base;
227 ++t_mat;
228 }
229
230 // move to next row base
231 ++t_row_base;
232 ++t_diff_row_base;
233 }
234
235 // this is obsolete for this particular example, we keep it for
236 // generality. in case of multi-physics problems different fields can
237 // chare hierarchical base but use different approx. order, so is
238 // possible to have more base functions than DOFs on element.
239 for (; rr < nb_row_base_functions; ++rr) {
240 ++t_row_base;
241 ++t_diff_row_base;
242 }
243
244 ++t_w;
245 }
246
247 // assemble system
248 CHKERR ::MatSetValues(getKSPB(), indicesRowSideMap[s0].size(),
249 &*indicesRowSideMap[s0].begin(),
250 indicesColSideMap[s1].size(),
251 &*indicesColSideMap[s1].begin(),
252 &*locMat.data().begin(), ADD_VALUES);
253
254 // stop of boundary element
255 if (!in_the_loop)
257 }
258 }
259 }
260
262 }
static Index< 'p', 3 > p
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
static double penalty
static double phi
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
std::array< MatrixDouble, 2 > rowBaseSideMap
std::array< VectorInt, 2 > indicesColSideMap
indices on columns for left hand-side
FTensor::Index< 'i', SPACE_DIM > i
std::array< MatrixDouble, 2 > colDiffBaseSideMap
std::array< MatrixDouble, 2 > colBaseSideMap
std::array< VectorInt, 2 > indicesRowSideMap
indices on rows for left hand-side
std::array< MatrixDouble, 2 > rowDiffBaseSideMap
static double nitsche
auto getFTensor0IntegrationWeight()
Get integration weights.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

Member Data Documentation

◆ locMat

MatrixDouble Poisson2DiscontGalerkinOperators::OpL2LhsPenalty::locMat
private

local operator matrix

Examples
PoissonDiscontinousGalerkin.hpp.

Definition at line 267 of file PoissonDiscontinousGalerkin.hpp.

◆ sideFEPtr

boost::shared_ptr<FaceSideEle> Poisson2DiscontGalerkinOperators::OpL2LhsPenalty::sideFEPtr
private

pointer to element to get data on edge/face sides

Examples
PoissonDiscontinousGalerkin.hpp.

Definition at line 266 of file PoissonDiscontinousGalerkin.hpp.


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