Medusa  1.1
Coordinate Free Mehless Method implementation
mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type > Class Template Reference

#include <ImplicitOperators_fwd.hpp>

Detailed Description

template<class shape_storage_type, class matrix_type, class rhs_type>
class mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >

This class represents implicit operators that fill given matrix M and right hand side rhs with appropriate coefficients approximating differential operators with shape functions from given shape storage ss.

Warning
The operators grad(), value() etc. are lazy and will not write to the matrix until summed with another operator or assigned to. If eager evaluation is required, use the eval() method on each operation.
Each evaluated operation is only added to the matrix and rhs. Setting the same equation twice thus has the effect of multiplying it by 2. Therefore, the matrix and the rhs should be initialized (usually to zero) before setting any equations.
Template Parameters
shape_storage_typeAny shape storage class satisfying the Shape storage concept.
matrix_typeUsually and Eigen sparse matrix.
rhs_typeUsually an Eigen VectorXd.
See also
ImplicitVectorOperators

Usage example:

Eigen::SparseMatrix<double, Eigen::RowMajor> M(N, N);
M.reserve(shapes.supportSizes());
Eigen::VectorXd rhs(N); rhs.setZero();
auto op = shapes.implicitOperators(M, rhs);
std::cout << op << std::endl;
double alpha = 0.25;
for (int i : domain.interior()) {
double x = domain.pos(i, 0), y = domain.pos(i, 1);
op.grad(i, {y, x}) + alpha * op.lap(i) =
4*(alpha + x*y) + std::exp(x) * (2*x*std::cos(2*y) - (3*alpha - y)*std::sin(2*y));
}
// Set boundary conditions
for (int i : domain.types() == -1) {
double y = domain.pos(i, 1);
op.value(i) = y*y + std::sin(2*y);
}
for (int i : domain.types() == -2) {
double y = domain.pos(i, 1);
op.value(i) = 1 + y*y + std::exp(1) * std::sin(2*y);
}
for (int i : domain.types() == -3) {
double x = domain.pos(i, 0);
op.value(i) = x*x;
}
for (int i : domain.types() == -4) {
double x = domain.pos(i, 0);
op.value(i) = 1 + x*x + std::sin(2) * std::exp(x);
}

Definition at line 40 of file ImplicitOperators_fwd.hpp.

Public Member Functions

 ImplicitOperators ()
 Default constructor sets offset to 0 and pointers to nullptr. More...
 
 ImplicitOperators (const shape_storage_t &ss)
 Set only shape storage, the rest as default constructor. More...
 
 ImplicitOperators (const shape_storage_t &ss, matrix_t &M, rhs_t &rhs, int row_offset=0, int col_offset=0)
 Set shape storage and problem matrix and rhs. More...
 
void setProblem (matrix_t &M, rhs_t &rhs, int row_offset=0, int col_offset=0)
 Sets current matrix and right hand side. More...
 
void setRowOffset (int row_offset)
 Sets row offset for given matrix, treating it as is the first row had index row_offset. More...
 
void setColOffset (int col_offset)
 Sets col offset for given matrix, treating it as is the first column had index col_offset. More...
 
bool hasShapes () const
 Returns true if operators have a non-null pointer to storage. More...
 
bool hasMatrix () const
 Returns true if operators have a non-null pointer to problem matrix. More...
 
bool hasRhs () const
 Returns true if operators have a non-null pointer to problem right hand side. More...
 
ValueOp value (int node, int row)
 Sets implicit equation that value of a field is equal to some other value. More...
 
ValueOp value (int node)
 Same as value with row equal to current node. More...
 
BasicOp< Lap< dim > > lap (int node, int row)
 Sets implicit equation that Laplacian of a field \( u\) at node-th point is equal to some value. More...
 
BasicOp< Lap< dim > > lap (int node)
 Same as lap with row equal to current node. More...
 
GradOp grad (int node, const vector_t &v, int row)
 Sets implicit equation that gradient of a field \( u\) at node-th point along v is equal to some value. More...
 
GradOp grad (int node, const vector_t &v)
 Same as grad with row equal to current node. More...
 
GradOp neumann (int node, const vector_t &unit_normal, int row)
 Sets neumann boundary conditions in node node. More...
 
GradOp neumann (int node, const vector_t &n)
 Same as neumann with row equal to current node. More...
 
BasicOp< Der1s< dim > > der1 (int node, int var, int row)
 Sets implicit equation that derivative of a field \( u\) at node-th point with respect to var is equal to some value. More...
 
BasicOp< Der1s< dim > > der1 (int node, int var)
 Same as der1 with row equal to current node. More...
 
BasicOp< Der2s< dim > > der2 (int node, int varmin, int varmax, int row)
 Sets implicit equation that second derivative of a field \( u\) at node-th point with respect to varmin and varmax is equal to some value. More...
 
BasicOp< Der2s< dim > > der2 (int node, int varmin, int varmax)
 Same as der2 with row equal to current node. More...
 
template<typename op_family_t >
BasicOp< op_family_t > apply (int node, typename op_family_t::operator_t o, int row)
 Add the weights for operator o to the appropriate elements in the matrix. More...
 
template<typename op_family_t >
BasicOp< op_family_t > apply (int node, typename op_family_t::operator_t o)
 Same as apply with row equal to current node. More...
 
template<typename op_family_t >
BasicOp< op_family_t > apply (int node)
 Overload for default-constructible operator. More...
 

Public Types

enum  { dim = shape_storage_t::dim }
 Store dimension of the domain. More...
 
typedef shape_storage_type shape_storage_t
 Type of shape storage. More...
 
typedef matrix_type matrix_t
 Matrix type. More...
 
typedef rhs_type rhs_t
 Right hand side type. More...
 
typedef shape_storage_t::vector_t vector_t
 Vector type. More...
 
typedef matrix_t::Scalar scalar_t
 Scalar type of matrix elements. More...
 

Classes

class  BasicOp
 Class representing Basic operators i.e. Der1, Der2, Lap or user defined custom operators. More...
 
class  GradOp
 Class representing the directional derivative (gradient) operation. More...
 
class  OpBase
 Base class for all elementary implicit operations. More...
 
class  RowOp
 Class representing an operation on a specific row of the matrix. More...
 
class  ValueOp
 Class representing the "evaluate" operation, i.e. the zero-th derivative. More...
 

Friends

template<typename S , typename M , typename R >
std::ostream & operator<< (std::ostream &os, const ImplicitOperators< S, M, R > &op)
 Output basic info about given operators. More...
 

Private Member Functions

void addToM (int i, int j, scalar_t v)
 Add v to M(i, j). More...
 
void addToRhs (int i, scalar_t v)
 Adds v to rhs[i]. More...
 

Private Attributes

const shape_storage_tss
 Shape storage, but name is shortened for readability. More...
 
matrix_tM
 Pointer to problem matrix. More...
 
rhs_trhs
 Pointer to right hand side. More...
 
int row_offset
 Row offset to be used when accessing matrix or rhs coefficients. More...
 
int col_offset
 Column offset to be used when accessing matrix coefficients. More...
 

Member Enumeration Documentation

◆ anonymous enum

template<class shape_storage_type , class matrix_type , class rhs_type >
anonymous enum

Store dimension of the domain.

Enumerator
dim 

Dimensionality of the function domain.

Definition at line 50 of file ImplicitOperators_fwd.hpp.

Constructor & Destructor Documentation

◆ ImplicitOperators() [1/3]

template<class shape_storage_type , class matrix_type , class rhs_type >
mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::ImplicitOperators ( )
inline

Default constructor sets offset to 0 and pointers to nullptr.

Definition at line 225 of file ImplicitOperators_fwd.hpp.

◆ ImplicitOperators() [2/3]

template<class shape_storage_type , class matrix_type , class rhs_type >
mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::ImplicitOperators ( const shape_storage_t ss)
inlineexplicit

Set only shape storage, the rest as default constructor.

Definition at line 227 of file ImplicitOperators_fwd.hpp.

◆ ImplicitOperators() [3/3]

template<class shape_storage_type , class matrix_type , class rhs_type >
mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::ImplicitOperators ( const shape_storage_t ss,
matrix_t M,
rhs_t rhs,
int  row_offset = 0,
int  col_offset = 0 
)

Set shape storage and problem matrix and rhs.

Parameters
ssClass storing all computed shapes.
MProblem matrix. Must have at least ss->size() rows and cols.
rhsProblem right hand side. Must have at least ss->size() rows.
row_offsetInstead of counting rows from 0, count them from row row_offset.
col_offsetInstead of counting columns from 0, count them from row col_offset.
Warning
Matrices M and rhs have values only added to them and should be zero initialized by the user. Since this is a common mistake, a warning is printed if this is not the case when in debug mode.

This class is usually constructed directly from shape storage using the implicitOperators() member function.

Definition at line 15 of file ImplicitOperators.hpp.

Member Function Documentation

◆ addToM()

template<class S , class matrix_type , class rhs_type >
void mm::ImplicitOperators< S, matrix_type, rhs_type >::addToM ( int  i,
int  j,
scalar_t  v 
)
private

Add v to M(i, j).

Exceptions
Assertionfails if v is none or i or j are out of range.

Definition at line 41 of file ImplicitOperators.hpp.

◆ addToRhs()

template<class shape_storage_type , class matrix_type , class rhs_type >
void mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::addToRhs ( int  i,
scalar_t  v 
)
private

Adds v to rhs[i].

Exceptions
Assertionfails if v is none or i is out of range.

Definition at line 54 of file ImplicitOperators.hpp.

◆ apply() [1/3]

template<class shape_storage_type , class matrix_type , class rhs_type >
template<typename op_family_t >
BasicOp<op_family_t> mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::apply ( int  node)
inline

Overload for default-constructible operator.

Definition at line 439 of file ImplicitOperators_fwd.hpp.

◆ apply() [2/3]

template<class shape_storage_type , class matrix_type , class rhs_type >
template<typename op_family_t >
BasicOp<op_family_t> mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::apply ( int  node,
typename op_family_t::operator_t  o 
)
inline

Same as apply with row equal to current node.

Definition at line 433 of file ImplicitOperators_fwd.hpp.

◆ apply() [3/3]

template<class shape_storage_type , class matrix_type , class rhs_type >
template<typename op_family_t >
BasicOp<op_family_t> mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::apply ( int  node,
typename op_family_t::operator_t  o,
int  row 
)
inline

Add the weights for operator o to the appropriate elements in the matrix.

Parameters
node
oOperator family
rowWrite in this matrix row.
Returns
This function is lazy and returns an "operation" which will write the weights in the matrix, when evaluated.

Definition at line 427 of file ImplicitOperators_fwd.hpp.

◆ der1() [1/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
BasicOp<Der1s<dim> > mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::der1 ( int  node,
int  var 
)
inline

Same as der1 with row equal to current node.

Definition at line 386 of file ImplicitOperators_fwd.hpp.

◆ der1() [2/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
BasicOp<Der1s<dim> > mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::der1 ( int  node,
int  var,
int  row 
)
inline

Sets implicit equation that derivative of a field \( u\) at node-th point with respect to var is equal to some value.

The code alpha*op.der1(node, 0) = v fills the node-th row of matrix M so that node-th row of the equation \( M u = v \) is a good approximation of the equation

\[ \alpha \dpar{u}{x}(p) = v(p), \]

where \(p\) is the node-th point and \(x\) is the 0-th variable.

Parameters
nodeIndex of a node from 0 to N for which to write the equation for. This means only node-th row of the matrix is changed. User must make sure that the matrix is large enough.
varVariable with respect to which to derive.
rowWrite equation in this specific row. Row with index node is chosen by default.
Returns
A class representing this operation for lazy evaluation purposes.

Definition at line 382 of file ImplicitOperators_fwd.hpp.

◆ der2() [1/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
BasicOp<Der2s<dim> > mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::der2 ( int  node,
int  varmin,
int  varmax 
)
inline

Same as der2 with row equal to current node.

Definition at line 413 of file ImplicitOperators_fwd.hpp.

◆ der2() [2/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
BasicOp<Der2s<dim> > mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::der2 ( int  node,
int  varmin,
int  varmax,
int  row 
)
inline

Sets implicit equation that second derivative of a field \( u\) at node-th point with respect to varmin and varmax is equal to some value.

The code alpha*op.der2(node, 0, 1) = v fills the node-th row of matrix M so that node-th row of the equation \( M u = v \) is a good approximation of the equation

\[ \alpha \dpar{^2u}{x\partial y}(p) = v(p), \]

where \(p\) is the node-th point and \(x\) and \(y\) are the 0-th and the 1-st variables.

Parameters
nodeIndex of a node from 0 to N for which to write the equation for. This means only node-th row of the matrix is changed. User must make sure that the matrix is large enough.
varminSmaller of the two variables with respect to which to derive.
varmaxGrater of the two variables with respect to which to derive.
rowWrite equation in this specific row. Row with index node is chosen by default.
Returns
A class representing this operation for lazy evaluation purposes.

Definition at line 409 of file ImplicitOperators_fwd.hpp.

◆ grad() [1/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
GradOp mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::grad ( int  node,
const vector_t v 
)
inline

Same as grad with row equal to current node.

Definition at line 337 of file ImplicitOperators_fwd.hpp.

◆ grad() [2/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
GradOp mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::grad ( int  node,
const vector_t v,
int  row 
)
inline

Sets implicit equation that gradient of a field \( u\) at node-th point along v is equal to some value.

The code alpha*op.grad(node, v) = r fills the node-th row of matrix M so that node-th row of the equation \( M u = r \) is a good approximation of the equation

\[ \alpha \vec{v} \cdot (\nabla u)(p) = r(p), \]

where \(p\) is the node-th point and \(\vec{v}\) is the vector v. Alternatively, this can be viewed as setting the directional derivative along v to be equal to r.

Parameters
nodeIndex of a node from 0 to N for which to write the equation for. This means only node-th row of the matrix is changed. User must make sure that the matrix is large enough.
vVector to multiply the gradient with.
rowWrite equation in this specific row. Row with index node is chosen by default.
Returns
A class representing this operation for lazy evaluation purposes.

Example:

// Solve problem for S unknown scalar field:
// grad S . v + laplace S = phi
// grad S . (y, x) + alfa laplace S = 4 (alfa + xy) + e^x(2 x cos(2y) - (3 alfa - y) sin(2y))
// on [0, 1] x [0, 1] with BC:
// S(x, 0) = x^2
// S(x, 1) = 1 + x^2 + sin(2) e^x
// S(0, y) = y^2 + sin(2y)
// S(1, y) = 1 + y^2 + e * sin(2y)
// with solution:
// S(x, y) = e^x sin(2y) + x^2 + y^2
std::function<Vec2d(Vec2d)> analytical = [] (const Vec2d& p) {
return std::exp(p[0]) * std::sin(2*p[1]) + p[0]*p[0] + p[1]*p[1];
};
// Prepare domain
BoxShape<Vec2d> box(0.0, 1.0);
DomainDiscretization<Vec2d> domain = box.discretizeWithStep(0.05);
int N = domain.size();
domain.findSupport(FindClosest(9));
// Prepare operators and matrix
WLS<Monomials<Vec2d>, GaussianWeight<Vec2d>> mls(2, 15.0 / N);
auto shapes = domain.computeShapes<sh::lap|sh::grad>(mls, domain.interior());
Eigen::SparseMatrix<double, Eigen::RowMajor> M(N, N);
M.reserve(shapes.supportSizes());
Eigen::VectorXd rhs(N); rhs.setZero();
auto op = shapes.implicitOperators(M, rhs);
double alpha = 0.25;
for (int i : domain.interior()) {
double x = domain.pos(i, 0), y = domain.pos(i, 1);
op.grad(i, {y, x}) + alpha * op.lap(i) =
4*(alpha + x*y) + std::exp(x) * (2*x*std::cos(2*y) - (3*alpha - y)*std::sin(2*y));
}
// Set boundary conditions
for (int i : domain.types() == -1) {
double y = domain.pos(i, 1);
op.value(i) = y*y + std::sin(2*y);
}
for (int i : domain.types() == -2) {
double y = domain.pos(i, 1);
op.value(i) = 1 + y*y + std::exp(1) * std::sin(2*y);
}
for (int i : domain.types() == -3) {
double x = domain.pos(i, 0);
op.value(i) = x*x;
}
for (int i : domain.types() == -4) {
double x = domain.pos(i, 0);
op.value(i) = 1 + x*x + std::sin(2) * std::exp(x);
}
// Sparse solve
Eigen::BiCGSTAB<Eigen::SparseMatrix<double, Eigen::RowMajor>> solver;
solver.compute(M);
Eigen::VectorXd sol = solver.solve(rhs);
for (int i = 0; i < N; ++i) {
Vec2d correct = analytical(domain.pos(i));
ASSERT_NEAR(sol[i], correct[0], 5e-4);
}

Definition at line 335 of file ImplicitOperators_fwd.hpp.

◆ hasMatrix()

template<class shape_storage_type , class matrix_type , class rhs_type >
bool mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::hasMatrix ( ) const
inline

Returns true if operators have a non-null pointer to problem matrix.

Definition at line 268 of file ImplicitOperators_fwd.hpp.

◆ hasRhs()

template<class shape_storage_type , class matrix_type , class rhs_type >
bool mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::hasRhs ( ) const
inline

Returns true if operators have a non-null pointer to problem right hand side.

Definition at line 270 of file ImplicitOperators_fwd.hpp.

◆ hasShapes()

template<class shape_storage_type , class matrix_type , class rhs_type >
bool mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::hasShapes ( ) const
inline

Returns true if operators have a non-null pointer to storage.

Definition at line 266 of file ImplicitOperators_fwd.hpp.

◆ lap() [1/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
BasicOp<Lap<dim> > mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::lap ( int  node)
inline

Same as lap with row equal to current node.

Definition at line 311 of file ImplicitOperators_fwd.hpp.

◆ lap() [2/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
BasicOp<Lap<dim> > mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::lap ( int  node,
int  row 
)
inline

Sets implicit equation that Laplacian of a field \( u\) at node-th point is equal to some value.

The code alpha*op.lap(node) = v fills the node-th row of matrix M so that node-th row of the equation \( M u = v \) is a good approximation of the equation

\[ \alpha \nabla^2 u(p) = v(p), \]

where \(p\) is the node-th point.

Parameters
nodeIndex of a node from 0 to N for which to write the equation for. This means only the node-th row of the matrix is changed. User must make sure that the matrix is large enough.
rowWrite equation in this specific row. Row with index node is chosen by default.
Returns
A class representing this operation for lazy evaluation purposes.

Example:

// Solve 1-D boundary value problem
// u''(x) = x; u(0) = 1, u(1) = -2
// with solution u(x) = 1/6 (6 - 19 x + x^3)
BoxShape<Vec1d> box(0.0, 1.0);
DomainDiscretization<Vec1d> domain = box.discretizeWithStep(0.1);
WLS<Monomials<Vec1d>> wls(4);
int N = domain.size();
domain.findSupport(FindClosest(5));
auto shapes = domain.computeShapes<sh::lap>(wls, domain.interior());
Eigen::MatrixXd M = Eigen::MatrixXd::Zero(N, N);
Eigen::VectorXd rhs(N); rhs.setZero();
auto op = shapes.implicitOperators(M, rhs);
for (int i : domain.interior()) {
op.lap(i) = domain.pos(i, 0);
}
for (int i : domain.boundary()) {
op.value(i) = (domain.pos(i, 0) == 0) ? 1 : -2;
}
Eigen::VectorXd sol = M.lu().solve(rhs);
std::function<double(double)> analytic = [](double x) { return 1/6.0 * (6 - 19 * x + x*x*x); };
for (int i = 0; i < domain.size(); ++i) {
EXPECT_NEAR(sol[i], analytic(domain.pos(i, 0)), 1e-14);
}

Definition at line 309 of file ImplicitOperators_fwd.hpp.

◆ neumann() [1/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
GradOp mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::neumann ( int  node,
const vector_t n 
)
inline

Same as neumann with row equal to current node.

Definition at line 361 of file ImplicitOperators_fwd.hpp.

◆ neumann() [2/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
GradOp mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::neumann ( int  node,
const vector_t unit_normal,
int  row 
)
inline

Sets neumann boundary conditions in node node.

The code alpha*op.neumann(node, normal) = v fills the node-th row of matrix M so that node-th row of the equation \( M u = v \) is a good approximation of the equation

\[ \alpha \dpar{u}{\vec{n}}(p) = v(p), \]

where \(p\) is the node-th point and \(\vec{n}\) is the unit_normal.

This is the same as using grad(), but has additional semantic meaning of setting the boundary conditions.

Returns
A class representing this operation for lazy evaluation purposes.
See also
grad

Definition at line 355 of file ImplicitOperators_fwd.hpp.

◆ setColOffset()

template<class shape_storage_type , class matrix_type , class rhs_type >
void mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::setColOffset ( int  col_offset)
inline

Sets col offset for given matrix, treating it as is the first column had index col_offset.

Definition at line 260 of file ImplicitOperators_fwd.hpp.

◆ setProblem()

template<class shape_storage_type , class matrix_type , class rhs_type >
void mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::setProblem ( matrix_t M,
rhs_t rhs,
int  row_offset = 0,
int  col_offset = 0 
)
inline

Sets current matrix and right hand side.

Definition at line 248 of file ImplicitOperators_fwd.hpp.

◆ setRowOffset()

template<class shape_storage_type , class matrix_type , class rhs_type >
void mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::setRowOffset ( int  row_offset)
inline

Sets row offset for given matrix, treating it as is the first row had index row_offset.

Definition at line 255 of file ImplicitOperators_fwd.hpp.

◆ value() [1/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
ValueOp mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::value ( int  node)
inline

Same as value with row equal to current node.

Definition at line 289 of file ImplicitOperators_fwd.hpp.

◆ value() [2/2]

template<class shape_storage_type , class matrix_type , class rhs_type >
ValueOp mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::value ( int  node,
int  row 
)
inline

Sets implicit equation that value of a field is equal to some other value.

The code alpha*op.value(node) = v fills the node-th row of matrix M so that node-th row of the equation \( M u = v \) is a good approximation of the equation

\[ \alpha u(p) = v(p), \]

where \(p\) is the node-th point.

Parameters
nodeIndex of a node from 0 to N for which to write the equation for. This means only the node-th row of the matrix is changed. User must make sure that the matrix is large enough.
rowWrite equation in this specific row. Row with index node is chosen by default.
Returns
A class representing this operation for lazy evaluation purposes.

Definition at line 287 of file ImplicitOperators_fwd.hpp.

Friends And Related Function Documentation

◆ operator<<

template<class shape_storage_type , class matrix_type , class rhs_type >
template<typename S , typename M , typename R >
std::ostream& operator<< ( std::ostream &  os,
const ImplicitOperators< S, M, R > &  op 
)
friend

Output basic info about given operators.

Definition at line 65 of file ImplicitOperators.hpp.

Member Data Documentation

◆ col_offset

template<class shape_storage_type , class matrix_type , class rhs_type >
int mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::col_offset
private

Column offset to be used when accessing matrix coefficients.

Definition at line 58 of file ImplicitOperators_fwd.hpp.

◆ M

template<class shape_storage_type , class matrix_type , class rhs_type >
matrix_t* mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::M
private

Pointer to problem matrix.

Definition at line 55 of file ImplicitOperators_fwd.hpp.

◆ rhs

template<class shape_storage_type , class matrix_type , class rhs_type >
rhs_t* mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::rhs
private

Pointer to right hand side.

Definition at line 56 of file ImplicitOperators_fwd.hpp.

◆ row_offset

template<class shape_storage_type , class matrix_type , class rhs_type >
int mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::row_offset
private

Row offset to be used when accessing matrix or rhs coefficients.

Definition at line 57 of file ImplicitOperators_fwd.hpp.

◆ ss

template<class shape_storage_type , class matrix_type , class rhs_type >
const shape_storage_t* mm::ImplicitOperators< shape_storage_type, matrix_type, rhs_type >::ss
private

Shape storage, but name is shortened for readability.

Definition at line 54 of file ImplicitOperators_fwd.hpp.


The documentation for this class was generated from the following files:
mm::sh::lap
static const shape_flags lap
Indicates to calculate laplace shapes.
Definition: shape_flags.hpp:24
mm::ImplicitOperators::rhs
rhs_t * rhs
Pointer to right hand side.
Definition: ImplicitOperators_fwd.hpp:56
mm::ImplicitOperators::M
matrix_t * M
Pointer to problem matrix.
Definition: ImplicitOperators_fwd.hpp:55
mm::sh::grad
static const shape_flags grad
Indicates to prepare all shapes needed for grad.
Definition: shape_flags.hpp:27
mm::Vec2d
Vec< double, 2 > Vec2d
Convenience typedef for 2d vector of doubles.
Definition: Vec_fwd.hpp:34