11#ifndef __POISSON3DHOMOGENEOUS_HPP__
12#define __POISSON3DHOMOGENEOUS_HPP__
41 const int nb_row_dofs = row_data.
getIndices().size();
42 const int nb_col_dofs = col_data.
getIndices().size();
44 if (nb_row_dofs && nb_col_dofs) {
46 locLhs.resize(nb_row_dofs, nb_col_dofs,
false);
53 const int nb_integration_points =
getGaussPts().size2();
61 for (
int gg = 0; gg != nb_integration_points; gg++) {
62 const double a = t_w * volume;
64 for (
int rr = 0; rr != nb_row_dofs; ++rr) {
68 for (
int cc = 0; cc != nb_col_dofs; cc++) {
69 locLhs(rr, cc) += t_row_diff_base(
i) * t_col_diff_base(
i) *
a;
86 if (row_side != col_side || row_type != col_type) {
87 transLocLhs.resize(nb_col_dofs, nb_row_dofs,
false);
112 locRhs.resize(nb_dofs,
false);
119 const int nb_integration_points =
getGaussPts().size2();
127 for (
int gg = 0; gg != nb_integration_points; gg++) {
128 const double a = t_w * volume;
130 for (
int rr = 0; rr != nb_dofs; rr++) {
144 CHKERR VecSetOption(
getKSPf(), VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
#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.
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.
FTensor::Index< 'i', 3 > i
constexpr auto field_name
bool sYmm
If true assume that matrix is symmetric structure.
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.
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
Volume finite element base.
OpDomainLhsMatrixK(std::string row_field_name, std::string col_field_name)
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.
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpDomainRhsVectorF(std::string field_name)