11#ifndef __POISSON2DISGALERKIN_HPP__
12#define __POISSON2DISGALERKIN_HPP__
22std::array<VectorInt, 2>
24std::array<VectorInt, 2>
43 std::fill(&doEntities[MBVERTEX], &doEntities[MBMAXTYPE],
false);
45 for (
auto t = moab::CN::TypeDimensionMap[
SPACE_DIM].first;
46 t <= moab::CN::TypeDimensionMap[
SPACE_DIM].second; ++
t)
61 if ((CN::Dimension(row_type) ==
SPACE_DIM) &&
63 const auto nb_in_loop = getFEMethod()->nInTheLoop;
70 areaMap[nb_in_loop] = getMeasure();
71 senseMap[nb_in_loop] = getSkeletonSense();
92 &*base_mat.data().begin());
95template <
typename T>
inline auto get_ntensor(
T &base_mat,
int gg,
int bb) {
96 double *ptr = &base_mat(gg, bb);
101 double *ptr = &*base_mat.data().begin();
102 return getFTensor1FromPtr<2>(ptr);
107 double *ptr = &base_mat(gg, 2 * bb);
108 return getFTensor1FromPtr<2>(ptr);
127 EntitiesFieldData::EntData &data) {
132 const auto in_the_loop =
136 const std::array<std::string, 2> ele_type_name = {
"BOUNDARY",
"SKELETON"};
138 <<
"OpL2LhsPenalty inTheLoop " << ele_type_name[in_the_loop];
149 t_normal.normalize();
152 const size_t nb_integration_pts =
getGaussPts().size2();
156 const double beta =
static_cast<double>(
nitsche) / (in_the_loop + 1);
167 const auto sense_row =
senseMap[s0];
174 const auto sense_col =
senseMap[s1];
181 locMat.resize(nb_rows, nb_cols,
false);
190 for (
size_t gg = 0; gg != nb_integration_pts; ++gg) {
195 auto t_mat =
locMat.data().begin();
199 for (; rr != nb_rows; ++rr) {
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);
209 auto t_diff_col_base =
213 for (
size_t cc = 0; cc != nb_cols; ++cc) {
217 t_un(
i) = -
p * (t_col_base * t_normal(
i) * sense_col -
218 beta * t_diff_col_base(
i) /
p);
221 *t_mat -= alpha * (t_vn(
i) * t_un(
i));
222 *t_mat -= alpha * (t_vn_plus(
i) * (beta * t_diff_col_base(
i)));
239 for (; rr < nb_row_base_functions; ++rr) {
252 &*
locMat.data().begin(), ADD_VALUES);
265 boost::shared_ptr<FaceSideEle>
285 const auto in_the_loop =
294 t_normal.normalize();
304 rhsF.resize(nb_rows,
false);
307 const size_t nb_integration_pts =
getGaussPts().size2();
316 const double beta =
static_cast<double>(
nitsche) / (in_the_loop + 1);
318 for (
size_t gg = 0; gg != nb_integration_pts; ++gg) {
321 const double source_val =
322 -
p *
sourceFun(t_coords(0), t_coords(1), t_coords(2));
324 auto t_f =
rhsF.data().begin();
327 for (; rr != nb_rows; ++rr) {
330 t_vn_plus(
i) = beta * (
phi * t_diff_row_base(
i) /
p);
332 t_vn(
i) = t_row_base * t_normal(
i) * sense_row - t_vn_plus(
i);
335 *t_f -= alpha * t_vn(
i) * (source_val * t_normal(
i));
342 for (; rr < nb_row_base_functions; ++rr) {
360 boost::shared_ptr<FaceSideEle>
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MOFEM_LOG(channel, severity)
Log.
auto fun
Function to approximate.
auto get_ntensor(T &base_mat)
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< double, 2 > areaMap
std::array< MatrixDouble, 2 > colBaseSideMap
std::array< VectorInt, 2 > indicesRowSideMap
indices on rows for left hand-side
std::array< int, 2 > senseMap
auto get_diff_ntensor(T &base_mat)
std::array< MatrixDouble, 2 > rowDiffBaseSideMap
constexpr double t
plate stiffness
constexpr auto field_name
Data on single entity (This is passed as argument to DataOperator::doWork)
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of dofs on entity.
default operator for TRI element
double getMeasure()
get measure of element
auto getFTensor1Normal()
get normal as tensor
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
@ OPSPACE
operator do Work is execute on space data
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Operator tp collect data from elements on the side of Edge/Face.
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
OpCalculateSideData(std::string field_name, std::string col_field_name)
Operator to evaluate Dirichlet boundary conditions using DG.
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
VectorDouble rhsF
vector to store local operator right hand side
OpL2BoundaryRhs(boost::shared_ptr< FaceSideEle > side_fe_ptr, ScalarFun fun)
ScalarFun sourceFun
pointer to function to evaluate value of function on boundary
boost::shared_ptr< FaceSideEle > sideFEPtr
pointer to element to get data on edge/face sides
Operator the left hand side matrix.
MatrixDouble locMat
local operator matrix
OpL2LhsPenalty(boost::shared_ptr< FaceSideEle > side_fe_ptr)
Construct a new OpL2LhsPenalty.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
boost::shared_ptr< FaceSideEle > sideFEPtr
pointer to element to get data on edge/face sides