![Logo](MoFEMLogo.png) |
| v0.14.0
|
Go to the documentation of this file. 1 #ifndef __POISSON2DLAGRANGEMULTIPLIER_HPP__
2 #define __POISSON2DLAGRANGEMULTIPLIER_HPP__
20 typedef boost::function<
double(
const double,
const double,
const double)>
25 OpDomainLhsK(std::string row_field_name, std::string col_field_name,
26 boost::shared_ptr<std::vector<unsigned char>> boundary_marker =
nullptr)
33 EntityType col_type,
EntData &row_data,
37 const int nb_row_dofs = row_data.
getIndices().size();
38 const int nb_col_dofs = col_data.
getIndices().size();
40 if (nb_row_dofs && nb_col_dofs) {
42 locLhs.resize(nb_row_dofs, nb_col_dofs,
false);
49 const int nb_integration_points =
getGaussPts().size2();
57 for (
int gg = 0; gg != nb_integration_points; gg++) {
59 const double a = t_w * area;
61 for (
int rr = 0; rr != nb_row_dofs; ++rr) {
65 for (
int cc = 0; cc != nb_col_dofs; cc++) {
67 locLhs(rr, cc) += t_row_diff_base(
i) * t_col_diff_base(
i) *
a;
87 if (row_side != col_side || row_type != col_type) {
88 transLocLhs.resize(nb_col_dofs, nb_row_dofs,
false);
107 boost::shared_ptr<std::vector<unsigned char>> boundary_marker =
nullptr)
118 locRhs.resize(nb_dofs,
false);
125 const int nb_integration_points =
getGaussPts().size2();
135 for (
int gg = 0; gg != nb_integration_points; gg++) {
137 const double a = t_w * area;
141 for (
int rr = 0; rr != nb_dofs; rr++) {
177 EntityType col_type,
EntData &row_data,
181 const int nb_row_dofs = row_data.
getIndices().size();
182 const int nb_col_dofs = col_data.
getIndices().size();
184 if (nb_row_dofs && nb_col_dofs) {
193 const int nb_integration_points =
getGaussPts().size2();
201 for (
int gg = 0; gg != nb_integration_points; gg++) {
202 const double a = t_w * edge;
204 for (
int rr = 0; rr != nb_row_dofs; ++rr) {
208 for (
int cc = 0; cc != nb_col_dofs; cc++) {
256 const int nb_integration_points =
getGaussPts().size2();
266 for (
int gg = 0; gg != nb_integration_points; gg++) {
268 const double a = t_w * edge;
269 double boundary_term =
272 for (
int rr = 0; rr != nb_dofs; rr++) {
301 #endif //__POISSON2DLAGRANGEMULTIPLIER_HPP__
Data on single entity (This is passed as argument to DataOperator::doWork)
OpDomainRhsF(std::string field_name, ScalarFunc source_term_function, boost::shared_ptr< std::vector< unsigned char >> boundary_marker=nullptr)
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
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 transLocBoundaryLhs
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
@ OPROWCOL
operator doWork is executed on FE rows &columns
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
EntitiesFieldData::EntData EntData
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
OpBoundaryRhsG(std::string field_name, ScalarFunc boundary_function)
#define CHKERR
Inline error check.
friend class UserDataOperator
FTensor::Index< 'i', 2 > i
OpDomainLhsK(std::string row_field_name, std::string col_field_name, boost::shared_ptr< std::vector< unsigned char >> boundary_marker=nullptr)
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
default operator for TRI element
const VectorInt & getIndices() const
Get global indices of dofs on entity.
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
VectorDouble locBoundaryRhs
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
friend class UserDataOperator
constexpr auto field_name
MatrixDouble locBoundaryLhs
default operator for EDGE element
ScalarFunc sourceTermFunc
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
UBlasVector< double > VectorDouble
boost::function< double(const double, const double, const double)> ScalarFunc
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.
bool sYmm
If true assume that matrix is symmetric structure.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
OpBoundaryLhsC(std::string row_field_name, std::string col_field_name)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
@ OPROW
operator doWork function is executed on FE rows