Medusa  1.1
Coordinate Free Mehless Method implementation
Operators

Implements operators for intuitive equation notation. More...

Classes

class  mm::ExplicitOperators< shape_storage_type >
 A class for evaluating typical operators needed in spatial discretization. More...
 
class  mm::ExplicitVectorOperators< shape_storage_type >
 A class for evaluating typical operators needed in spatial discretization. More...
 
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. More...
 
class  mm::ImplicitVectorOperators< shape_storage_type, matrix_type, rhs_type >
 This class represents implicit vector 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. More...
 
class  mm::RaggedShapeStorage< vec_t, OpFamilies >
 Efficiently stores shape functions of different lengths. More...
 
class  mm::ShapeStorage< Derived, vec_t, OpFamilies >
 Shape storage base class. More...
 
class  mm::UniformShapeStorage< vec_t, OpFamilies >
 Efficiently stores shape functions of uniform length. More...
 

Functions

template<typename approx_t , typename shape_storage_t , typename ... OpFamilies>
void mm::computeShapes (const DomainDiscretization< typename approx_t::vector_t > &domain, approx_t approx, const indexes_t &indexes, const std::tuple< OpFamilies... > &operators, shape_storage_t *storage)
 Computes shape functions (stencil weights) for given nodes using support domains from domain and approximations provided by approx. More...
 

Detailed Description

Implements operators for intuitive equation notation.

Shapes defined by the approximation can be computed with computeShapes functions, which has built in parallelization.

Computed stencil weights / shape functions are stored in ShapeStorage and can be used to construct operators. These operators enable an intuitive and expressive syntax which directly corresponds to the mathematics notation.

For implicit operators, one simply writes the equations:

for (int i : domain.interior()) {
op.lap(i) = 0;
}
for (int i : domain.boundary()) {
op.value(i) = 1;
}

Explicit operators can be applied directly to scalar or vector fields:

Vec2d gd = op.graddiv(field, i);

Function Documentation

◆ computeShapes()

template<typename approx_t , typename shape_storage_t , typename ... OpFamilies>
void mm::computeShapes ( const DomainDiscretization< typename approx_t::vector_t > &  domain,
approx_t  approx,
const indexes_t indexes,
const std::tuple< OpFamilies... > &  operators,
shape_storage_t *  storage 
)

Computes shape functions (stencil weights) for given nodes using support domains from domain and approximations provided by approx.

The computed shapes are stored in storage.

Template Parameters
approx_tType of the approximation engine to use. Must satisfy the Approximation engine concept.
shape_storage_tType of the shape storage. Must satisfy the Shape storage concept.
OpFamiliesA list of operator families for which to compute the shapes. The families must satisfy the Operator family concept. The family types must be unique.
Parameters
[in]domainDomain discretization for which to compute the shape functions.
[in]approxApproximation engine specifying the shape computation procedure.
[in]indexesA set of indexes for which to compute the shape functions. Common values are domain.all(), domain.interior(), domain.boundary().
[in]operatorsObject representing operator families.
[out]storageAn object for storing the shape functions. The object is filled during computation, but must be appropriately reserved before computing shapes, by using ShapeStorage::resize.
Note
This function supports parallel execution if OpenMP is enabled in the compiler, e.g.\ for GCC: -fopenmp, for ICC: -openmp.

Usage example:

BoxShape<Vec2d> box(0.0, 1.0);
DomainDiscretization<Vec2d> d = box.discretizeWithStep(0.1);
d.findSupport(FindClosest(9));
WLS<Monomials<Vec2d>> approx(2);
auto shapes = d.computeShapes<sh::lap|sh::d1>(approx); // implicit call via domain hook
// direct call which overwrites shapes for Laplacian for internal nodes in given storage
computeShapes(d, approx, d.interior(), operators, &shapes);

Usage with (dummy) custom operator:

// Dummy zero operator
struct Zero : Operator<Zero> {
template <typename basis_t> typename basis_t::scalar_t applyAt0(
const basis_t&, int, const std::vector<typename basis_t::vector_t>&,
typename basis_t::scalar_t) const { return 0.0; }
};
std::tuple<Lap<2>, Der1s<2>, Zero> ops;
auto shapes = d.computeShapes<decltype(ops)>(ops, approx);
shapes.get<Zero>(3, 2); // coefficient 4 of shape for zero in 3rd node.
shapes.get<Der1s<2>>(0, 3, 4); // coefficient 4 of shape for dx in 3rd node.
Examples
test/end2end/diffusion_explicit.cpp, test/operators/computeShapes_test.cpp, test/operators/ExplicitOperators_test.cpp, and test/operators/ExplicitVectorOperators_test.cpp.
mm::sh::d1
static const shape_flags d1
Indicates to calculate d1 shapes.
Definition: shape_flags.hpp:23
scalar_t
Scalar scalar_t
Type of the elements, alias of Scalar.
Definition: MatrixBaseAddons.hpp:16
mm::sh::lap
static const shape_flags lap
Indicates to calculate laplace shapes.
Definition: shape_flags.hpp:24
mm::computeShapes
void computeShapes(const DomainDiscretization< typename approx_t::vector_t > &domain, approx_t approx, const indexes_t &indexes, const std::tuple< OpFamilies... > &operators, shape_storage_t *storage)
Computes shape functions (stencil weights) for given nodes using support domains from domain and appr...
mm::sh::operator_tuple::type
void type
The type of the operator tuple.
Definition: shape_flags.hpp:43
mm::Vec2d
Vec< double, 2 > Vec2d
Convenience typedef for 2d vector of doubles.
Definition: Vec_fwd.hpp:34