1#ifndef __POISSON2DNONHOMOGENEOUS_HPP__
2#define __POISSON2DNONHOMOGENEOUS_HPP__
18typedef boost::function<
double(
const double,
const double,
const double)>
23 OpDomainLhs(std::string row_field_name, std::string col_field_name)
29 EntityType col_type,
EntData &row_data,
33 const int nb_row_dofs = row_data.
getIndices().size();
34 const int nb_col_dofs = col_data.
getIndices().size();
36 if (nb_row_dofs && nb_col_dofs) {
38 locLhs.resize(nb_row_dofs, nb_col_dofs,
false);
45 const int nb_integration_points =
getGaussPts().size2();
53 for (
int gg = 0; gg != nb_integration_points; gg++) {
55 const double a = t_w * area;
57 for (
int rr = 0; rr != nb_row_dofs; ++rr) {
61 for (
int cc = 0; cc != nb_col_dofs; cc++) {
63 locLhs(rr, cc) += t_row_diff_base(
i) * t_col_diff_base(
i) *
a;
84 if (row_side != col_side || row_type != col_type) {
85 transLocLhs.resize(nb_col_dofs, nb_row_dofs,
false);
115 locRhs.resize(nb_dofs,
false);
122 const int nb_integration_points =
getGaussPts().size2();
132 for (
int gg = 0; gg != nb_integration_points; gg++) {
134 const double a = t_w * area;
138 for (
int rr = 0; rr != nb_dofs; rr++) {
140 locRhs[rr] += t_base * body_source *
a;
175 EntityType col_type,
EntData &row_data,
179 const int nb_row_dofs = row_data.
getIndices().size();
180 const int nb_col_dofs = col_data.
getIndices().size();
182 if (nb_row_dofs && nb_col_dofs) {
191 const int nb_integration_points =
getGaussPts().size2();
199 for (
int gg = 0; gg != nb_integration_points; gg++) {
200 const double a = t_w * edge;
202 for (
int rr = 0; rr != nb_row_dofs; ++rr) {
206 for (
int cc = 0; cc != nb_col_dofs; cc++) {
224 if (row_side != col_side || row_type != col_type) {
262 const int nb_integration_points =
getGaussPts().size2();
272 for (
int gg = 0; gg != nb_integration_points; gg++) {
274 const double a = t_w * edge;
275 double boundary_term =
278 for (
int rr = 0; rr != nb_dofs; rr++) {
#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.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
boost::function< double(const double, const double, const double)> ScalarFunc
FTensor::Index< 'i', 2 > i
constexpr auto field_name
bool sYmm
If true assume that matrix is symmetric structure.
default operator for EDGE element
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
const VectorInt & getIndices() const
Get global indices of dofs on entity.
default operator for TRI element
friend class UserDataOperator
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
@ OPROW
operator doWork function is executed on FE rows
@ OPROWCOL
operator doWork is executed on FE rows &columns
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
friend class UserDataOperator
MatrixDouble transLocBoundaryLhs
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
MatrixDouble locBoundaryLhs
OpBoundaryLhs(std::string row_field_name, std::string col_field_name)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
VectorDouble locBoundaryRhs
OpBoundaryRhs(std::string field_name, ScalarFunc boundary_function)
OpDomainLhs(std::string row_field_name, std::string col_field_name)
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
ScalarFunc sourceTermFunc
OpDomainRhs(std::string field_name, ScalarFunc source_term_function)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.