Medusa  1.1
Coordinate Free Mehless Method implementation
mm::DomainDiscretization< vec_t > Class Template Reference

#include <DomainDiscretization_fwd.hpp>

Detailed Description

template<class vec_t>
class mm::DomainDiscretization< vec_t >

Class representing domain discretization along with an associated shape.

Domain discretization consists of points in the interior and on the boundary. Each point has an associated type, with types of internal points being positive and types of boundary points being negative. Boundary points also have an outer unit normal normal associated with them.

Domains can be added and subtracted, nodes can be added and removed (removal currently at O(N) cost). Support / stencil indices can be stored and accessed.

Additionally, domain interiors can be filled, regularized and used when computing shapes.

Usage example:

BoxShape<Vec2d> b({0.0, 0.3}, {0.4, 2.7});
DomainDiscretization<Vec2d> d = b.discretizeWithStep(0.1);
std::cout << d.dim << std::endl; // dimension of the domain
d.size(); // number of points
int i = 0;
d.pos(i); // position of node i
if (d.type(i) < 0) {
d.normal(i); // boundary nodes have normals
}
d.boundary(); // indices of boundary nodes
d.interior(); // indices of interior nodes
d.shape(); // get geometric shape of the domain
d.addBoundaryNode({0.001, 0.3}, -2, {-1, 0}); // add new boundary point
d.addInternalNode({0.25, 0.75}, 2); // add new internal point
BallShape<Vec2d> c(0.0, 1.0);
auto d2 = c.discretizeBoundaryWithStep(0.1);
d -= d2; // subtract or add discretizations
std::cout << d << std::endl;
See also
DomainShape, FindClosest
Examples
test/end2end/poisson_explicit.cpp.

Definition at line 46 of file DomainDiscretization_fwd.hpp.

+ Collaboration diagram for mm::DomainDiscretization< vec_t >:

Public Member Functions

 DomainDiscretization (const DomainShape< vec_t > &shape)
 Construct an empty discretization for a given shape. More...
 
const Range< vec_t > & positions () const
 Returns positions of all nodes. More...
 
const vec_t & pos (int i) const
 Returns the position of i-th node. More...
 
vec_t & pos (int i)
 Returns writeable position of i-th node. More...
 
scalar_t pos (int i, int j) const
 Returns j-th coordinate of the position of i-th node. More...
 
const Range< Range< int > > & supports () const
 Returns support indices for all nodes. More...
 
const Range< int > & support (int i) const
 Returns support indices for i-th node. More...
 
Range< int > & support (int i)
 Returns writeable array of support indices for i-th node. More...
 
int support (int i, int j) const
 Returns j-th support node of i-th node. More...
 
Range< vec_t > supportNodes (int i) const
 Returns positions of support nodes of i-th node. More...
 
vec_t supportNode (int i, int j) const
 Returns position of j-th support node of i-th node. More...
 
scalar_t dr (int i) const
 Returns Euclidean distance to the second support node. More...
 
int supportSize (int i) const
 Returns size of i-th node support. More...
 
Range< int > supportSizes () const
 Returns a vector of support sizes for each node. More...
 
const Range< int > & types () const
 Returns types of all nodes. More...
 
Range< int > & types ()
 Returns mutable types of all nodes. More...
 
int type (int i) const
 Returns type of i-th node. More...
 
int & type (int i)
 Returns writeable type of i-th node. More...
 
const DomainShape< vec_t > & shape () const
 Returns geometric shape of the underlying domain. More...
 
const Range< int > & bmap () const
 Returns boundary map. More...
 
int bmap (int node) const
 Returns index of node node among only boundary nodes. More...
 
const Range< vec_t > & normals () const
 Returns normals of all boundary nodes. More...
 
const vec_t & normal (int i) const
 Returns outside unit normal of i-th node. More...
 
vec_t & normal (int i)
 Returns writable outside unit normal of i-th node. More...
 
indexes_t boundary () const
 Returns indexes of all boundary nodes. More...
 
indexes_t interior () const
 Returns indexes of all internal nodes. More...
 
indexes_t all () const
 Returns indexes of all nodes, i.e. {0, 1, ..., N-1}. More...
 
int size () const
 Returns N, the number of nodes in this discretization. More...
 
int addInternalNode (const vec_t &point, int type)
 Adds a single interior node with specified type to this discretization. More...
 
int addBoundaryNode (const vec_t &point, int type, const vec_t &normal)
 Adds a boundary node with given type and normal to the domain. More...
 
indexes_t addNodes (const DomainDiscretization< vec_t > &d)
 Add nodes from another discretization to this discretization. More...
 
void changeToBoundary (int i, const vec_t &point, int type, const vec_t &normal)
 Changes node i to boundary point with given type and normal. More...
 
void changeToInterior (int i, const vec_t &point, int type)
 Changes node i to interior point with given type. More...
 
Range< int > addGhostNodes (scalar_t h, int type=0)
 Overload of addGhostNodes with constant h and adds ghost nodes to all boundary nodes. More...
 
template<typename func_t >
Range< int > addGhostNodes (func_t h, int type=0)
 Overload of addGhostNodes which adds ghost nodes to all boundary nodes. More...
 
Range< int > addGhostNodes (scalar_t h, int type, const indexes_t &indexes)
 Overload of addGhostNodes with constant h. More...
 
template<typename func_t >
Range< int > addGhostNodes (func_t h, int type, const indexes_t &indexes)
 Add ghost or fictitious nodes to some boundary nodes. More...
 
void shuffleNodes (const indexes_t &permutation)
 Shuffles node order according to given permutation. More...
 
template<typename compare_t = std::less<std::pair<vec_t, int>>>
Range< int > reorderNodes (const compare_t &cmp=compare_t())
 Reorders nodes according to the compare function. More...
 
DomainDiscretizationadd (const DomainDiscretization &d)
 Merges given discretization to the current discretization. More...
 
DomainDiscretizationoperator+= (const DomainDiscretization &d)
 Operator version of DomainDiscretization::add. More...
 
DomainDiscretizationsubtract (const DomainDiscretization &d)
 Subtracts given discretized domain from *this. More...
 
DomainDiscretizationoperator-= (const DomainDiscretization &d)
 Operator version of DomainDiscretization::subtract. More...
 
DomainDiscretizationtranslate (const vec_t &a)
 Translate the domain by a given vector a. More...
 
DomainDiscretizationrotate (const Eigen::Matrix< scalar_t, dim, dim > &Q)
 Transform the domain by given orthogonal matrix Q. More...
 
DomainDiscretizationrotate (scalar_t angle)
 2D version of rotate accepting an angle. More...
 
void assert_is_valid () const
 Checks if domain is in valid state. More...
 
void clear ()
 Clears all data about this discretization. More...
 
void removeNodes (const Range< int > &to_remove)
 Remove nodes with given indices. This removes their types and potential normals as well. More...
 
void removeInternalNodes ()
 Removes all internal nodes. More...
 
void removeBoundaryNodes ()
 Removes all boundary nodes. More...
 
bool contains (const vec_t &point) const
 Returns true if point is inside the domain. More...
 
template<typename search_structure_t >
bool discreteContains (const vec_t &point, const search_structure_t &search) const
 A discrete version of contains using given boundary discretization. More...
 
template<typename search_structure_t >
bool contains (const vec_t &point, const search_structure_t &search) const
 Returns true if point is inside the domain. More...
 
template<typename search_structure_t >
void makeDiscreteContainsStructure (search_structure_t &search) const
 Fills the search structure search with boundary points for use with contains() or discreteContains(). More...
 
template<typename callable_t , typename... Args>
void findSupport (callable_t &callable, Args &&... args)
 Enables more readable calls to support engines. More...
 
template<typename callable_t , typename... Args>
void findSupport (const callable_t &callable, Args &&... args)
 Const version of DomainDiscretization::findSupport. More...
 
template<typename callable_t , typename... Args>
void fill (callable_t &callable, Args &&... args)
 Enables more readable calls to fill engines. More...
 
template<typename callable_t , typename... Args>
void fill (const callable_t &callable, Args &&... args)
 Const version of DomainDiscretization::fill. More...
 
template<typename callable_t , typename... Args>
void relax (callable_t &callable, Args &&... args)
 Enables more readable calls to relax engines. More...
 
template<typename callable_t , typename... Args>
void relax (const callable_t &callable, Args &&... args)
 Const version of DomainDiscretization::relax. More...
 
template<sh::shape_flags mask = sh::all, typename approx_t >
RaggedShapeStorage< vec_t, typename sh::operator_tuple< mask, vec_t::dim >::typecomputeShapes (approx_t approx, const indexes_t &indexes={}) const
 Compute shapes, specified with shape flags, for this domain with given approximation for given indexes. More...
 
template<typename OperatorTuple , typename approx_t >
RaggedShapeStorage< vec_t, OperatorTuple > computeShapes (approx_t approx, const indexes_t &indexes={}) const
 Compute shapes for default constructable families. More...
 
template<typename OperatorTuple , typename approx_t >
RaggedShapeStorage< vec_t, OperatorTuple > computeShapes (OperatorTuple operators, approx_t approx, const indexes_t &indexes={}) const
 Compute shapes for non-default constructable families. More...
 

Static Public Member Functions

template<typename hdf5_file >
static DomainDiscretization< vec_t > load (hdf5_file &file, const std::string &name)
 Load discretization from a HD5 file. More...
 

Public Types

enum  { dim = vec_t::dim }
 Store dimension of the domain. More...
 
typedef vec_t vector_t
 Vector data type used in computations. More...
 
typedef vec_t::Scalar scalar_t
 Scalar data type used in computation. More...
 

Friends

template<class V >
std::ostream & operator<< (std::ostream &os, const DomainDiscretization< V > &d)
 Output basic info about given domain. More...
 

Protected Attributes

Range< vec_t > positions_
 Positions of internal discretization points. More...
 
Range< Range< int > > support_
 Attribute used to store support points of each node. More...
 
Range< int > types_
 Storing types of nodes aligned with positions. More...
 
Range< int > boundary_map_
 Mapping index of a boundary node among all nodes, to its index among boundary nodes. More...
 
Range< vec_t > normals_
 List of normals of boundary nodes. More...
 
deep_copy_unique_ptr< DomainShape< vec_t > > shape_
 Geometric shape of the domain. More...
 

Private Member Functions

int addNode (const vec_t &point, int type)
 Helper for adding points to the domain, which keeps the object consistent. More...
 
void chopOff (const DomainDiscretization< vec_t > &d, const Range< int > &relevant)
 Remove a portion of this domain which is contained in d or is only dx away. More...
 
std::ostream & output_report (std::ostream &os=std::cout) const
 Outputs a simple report about out domain, like number of nodes, support size. More...
 

Member Enumeration Documentation

◆ anonymous enum

template<class vec_t >
anonymous enum

Store dimension of the domain.

Enumerator
dim 

Dimensionality of the domain.

Definition at line 51 of file DomainDiscretization_fwd.hpp.

Constructor & Destructor Documentation

◆ DomainDiscretization()

template<class vec_t >
mm::DomainDiscretization< vec_t >::DomainDiscretization ( const DomainShape< vec_t > &  shape)
explicit

Construct an empty discretization for a given shape.

Definition at line 239 of file DomainDiscretization.hpp.

Member Function Documentation

◆ add()

template<class vec_t >
DomainDiscretization< vec_t > & mm::DomainDiscretization< vec_t >::add ( const DomainDiscretization< vec_t > &  d)

Merges given discretization to the current discretization.

The shape of current domain changes to become the union of shape() and d.shape(). The discretization d is added to *this. Types of nodes are preserved. This makes it possible to set custom types before the += operation and refer to those types in the union as well.

Example:

double h = 0.1;
DomainDiscretization<Vec2d> rect = BoxShape<Vec2d>({-1, -1}, {1, 1}).discretizeWithStep(h);
rect.types()[rect.interior()] = 2;
DomainDiscretization<Vec2d> circ = BallShape<Vec2d>({0, 0}, 0.7).discretizeWithStep(h);
circ.types()[circ.boundary()] = -5; // -1 to -4 taken by rectangle sides
circ += rect;
std::cout << (circ.types() == 2).size() << std::endl; // nonzero number is printed

More specifically, when performing d1 += d2, first, boundary nodes in d2 that are contained in d1 are removed. Then, all nodes from d1 that are included in d2 are removed. Additionally, nodes from d1 are removed if they are too close to any of the nodes in d2. The threshold dx for point p in d1 is computed by finding the closest neighbour to p among the remaining nodes of d1 and then finding its closest neighbour among those same nodes. A factor of the distance between these last two nodes is taken as dx. The remaining nodes in d2 are then added to d1.

See the images below for visualization of the add operation.

See also
ShapeUnion, subtract

Definition at line 274 of file DomainDiscretization.hpp.

◆ addBoundaryNode()

template<class vec_t >
int mm::DomainDiscretization< vec_t >::addBoundaryNode ( const vec_t &  point,
int  type,
const vec_t &  normal 
)

Adds a boundary node with given type and normal to the domain.

Parameters
pointCoordinates of the node to add.
typeType of the point, must be negative.
normalOutside unit normal to the boundary at point point.
Returns
The index of the new node.
See also
addInternalNode

Definition at line 56 of file DomainDiscretization.hpp.

◆ addGhostNodes() [1/4]

template<class vec_t >
template<typename func_t >
Range< int > mm::DomainDiscretization< vec_t >::addGhostNodes ( func_t  h,
int  type,
const indexes_t indexes 
)

Add ghost or fictitious nodes to some boundary nodes.

Parameters
hNodal spacing function. The ghost node associated with some boundary node p and with normal n is added at position p + h(p)*n.
typeType of the added ghost node. The usual choice is type 0, making them neither internal nor boundary. This means that shapes in these nodes will not be computed, nor will their supports be found, however, they can appear as support nodes of other domain nodes. Note that this means that their support size is 0, which has to be taken into account when preallocating matrix sizes for implicit solving. If a positive number is chosen for type, the nodes are added as internal nodes and appear in the interior index list. Similarly, if type is negative, they are added as boundary nodes with normal n and appear in the boundary index list.
indexesA list of indexes for which to add boundary nodes. If not given, all boundary indices as given by boundary are assumed.
Returns
A list of indices of the size of the original domain, mapping each node index to the index of its corresponding ghost node, or -1, if it has no corresponding ghost node. Precisely the slots with indices in indexes have an entry not equal to -1.

This method can be useful when solving PDEs implicitly with Neumann boundary conditions, see the Ghost nodes example and the snippet below:

// 3D Poisson problem with mixed BC
double step = 0.1;
BallShape<Vec3d> sh(0, 1);
DomainDiscretization<Vec3d> domain = sh.discretizeBoundaryWithStep(step);
auto dir = domain.positions().filter([](const Vec3d& p) { return p[0] > 0; });
auto neu = domain.positions().filter([](const Vec3d& p) { return p[0] <= 0; });
GeneralFill<Vec3d> fill; fill.seed(20);
domain.fill(fill, step);
// Add ghost nodes and remember the mapping from boundary to ghost nodes.
Range<int> gh = domain.addGhostNodes(step, 0, neu);
int N = domain.size();
// By default, support is found for non-zero nodes, searching among all nodes
int n = 25;
domain.findSupport(FindClosest(n));
RBFFD<Polyharmonic<double, 3>, Vec3d, ScaleToClosest> approx({}, 2);
// Approximations are computed for non-zero nodes.
auto storage = domain.computeShapes(approx);
// Matrix size is equal to the count of all nodes.
Eigen::SparseMatrix<double, Eigen::RowMajor> M(N, N);
// Warning: support size vector contains zeros in ghost nodes, but we want matrix to allocate
// storage for ghost nodes as well.
M.reserve(Range<int>(N, n));
Eigen::VectorXd rhs(N); rhs.setZero();
auto op = storage.implicitOperators(M, rhs);
for (int i : domain.interior()) {
op.lap(i) = 1;
}
for (int i : dir) {
op.value(i) = 0;
}
for (int i : neu) {
op.neumann(i, domain.normal(i)) = 0;
op.lap(i, gh[i]) = 1;
}
Eigen::SparseLU<decltype(M)> solver;
solver.compute(M);
ScalarFieldd u = solver.solve(rhs);
Warning
Support nodes are not computed for ghost nodes, but ghost nodes are included in other nodes' supports. This means that support size in ghost nodes will be zero, also when returned by ShapeStorage::supportSizes, however, the matrix rows for ghost nodes must still be preallocated.
See also
addInternalNode, addBoundaryNode

Definition at line 155 of file DomainDiscretization.hpp.

◆ addGhostNodes() [2/4]

template<class vec_t >
template<typename func_t >
Range< int > mm::DomainDiscretization< vec_t >::addGhostNodes ( func_t  h,
int  type = 0 
)

Overload of addGhostNodes which adds ghost nodes to all boundary nodes.

Definition at line 149 of file DomainDiscretization.hpp.

◆ addGhostNodes() [3/4]

template<class vec_t >
Range< int > mm::DomainDiscretization< vec_t >::addGhostNodes ( scalar_t  h,
int  type,
const indexes_t indexes 
)

Overload of addGhostNodes with constant h.

Definition at line 142 of file DomainDiscretization.hpp.

◆ addGhostNodes() [4/4]

template<class vec_t >
Range< int > mm::DomainDiscretization< vec_t >::addGhostNodes ( scalar_t  h,
int  type = 0 
)

Overload of addGhostNodes with constant h and adds ghost nodes to all boundary nodes.

Definition at line 137 of file DomainDiscretization.hpp.

◆ addInternalNode()

template<class vec_t >
int mm::DomainDiscretization< vec_t >::addInternalNode ( const vec_t &  point,
int  type 
)

Adds a single interior node with specified type to this discretization.

Parameters
pointCoordinates of the node to add.
typeType of the node to add. Must be positive.
Returns
The index of the new node.
See also
addBoundaryNode

Definition at line 49 of file DomainDiscretization.hpp.

◆ addNode()

template<class vec_t >
int mm::DomainDiscretization< vec_t >::addNode ( const vec_t &  point,
int  type 
)
private

Helper for adding points to the domain, which keeps the object consistent.

Definition at line 66 of file DomainDiscretization.hpp.

◆ addNodes()

template<class vec_t >
indexes_t mm::DomainDiscretization< vec_t >::addNodes ( const DomainDiscretization< vec_t > &  d)

Add nodes from another discretization to this discretization.

The shape of the domain is not changed. Returns indices of newly added nodes.

Definition at line 104 of file DomainDiscretization.hpp.

◆ all()

template<class vec_t >
indexes_t mm::DomainDiscretization< vec_t >::all ( ) const
inline

Returns indexes of all nodes, i.e. {0, 1, ..., N-1}.

Definition at line 186 of file DomainDiscretization_fwd.hpp.

◆ assert_is_valid()

template<class vec_t >
void mm::DomainDiscretization< vec_t >::assert_is_valid

Checks if domain is in valid state.

This includes that all points are contained in the domain, that all boundary nodes have a normal, etc...

Definition at line 179 of file DomainDiscretization.hpp.

◆ bmap() [1/2]

template<class vec_t >
const Range<int>& mm::DomainDiscretization< vec_t >::bmap ( ) const
inline

Returns boundary map.

See also
boundary_map_

Definition at line 165 of file DomainDiscretization_fwd.hpp.

◆ bmap() [2/2]

template<class vec_t >
int mm::DomainDiscretization< vec_t >::bmap ( int  node) const
inline

Returns index of node node among only boundary nodes.

The returned index is in range [0, boundary().size()) if node is a boundary node, and -1 otherwise.

Definition at line 170 of file DomainDiscretization_fwd.hpp.

◆ boundary()

template<class vec_t >
indexes_t mm::DomainDiscretization< vec_t >::boundary ( ) const
inline

Returns indexes of all boundary nodes.

Definition at line 182 of file DomainDiscretization_fwd.hpp.

◆ changeToBoundary()

template<class vec_t >
void mm::DomainDiscretization< vec_t >::changeToBoundary ( int  i,
const vec_t &  point,
int  type,
const vec_t &  normal 
)

Changes node i to boundary point with given type and normal.

Definition at line 75 of file DomainDiscretization.hpp.

◆ changeToInterior()

template<class vec_t >
void mm::DomainDiscretization< vec_t >::changeToInterior ( int  i,
const vec_t &  point,
int  type 
)

Changes node i to interior point with given type.

Definition at line 88 of file DomainDiscretization.hpp.

◆ chopOff()

template<class vec_t >
void mm::DomainDiscretization< vec_t >::chopOff ( const DomainDiscretization< vec_t > &  d,
const Range< int > &  relevant 
)
private

Remove a portion of this domain which is contained in d or is only dx away.

The dx is computed locally from the discretization of d. No nodes are added.

Helper for add and subtract methods.

Parameters
dThe discretization of the domain which specifies which part to chop off.
relevantOnly the nodes with indices in relevant are used for computing the local dx. This is mostly a performance improvement which allows to construct a smaller k-d tree.

Definition at line 243 of file DomainDiscretization.hpp.

◆ clear()

template<class vec_t >
void mm::DomainDiscretization< vec_t >::clear

Clears all data about this discretization.

The state of the object is as if it were newly constructed using shape().

Definition at line 170 of file DomainDiscretization.hpp.

◆ computeShapes() [1/3]

template<class vec_t >
template<sh::shape_flags mask = sh::all, typename approx_t >
RaggedShapeStorage<vec_t, typename sh::operator_tuple<mask, vec_t::dim>::type> mm::DomainDiscretization< vec_t >::computeShapes ( approx_t  approx,
const indexes_t indexes = {} 
) const

Compute shapes, specified with shape flags, for this domain with given approximation for given indexes.

For complete documentation, se mm::computeShapes function.

◆ computeShapes() [2/3]

template<class vec_t >
template<typename OperatorTuple , typename approx_t >
RaggedShapeStorage<vec_t, OperatorTuple> mm::DomainDiscretization< vec_t >::computeShapes ( approx_t  approx,
const indexes_t indexes = {} 
) const

Compute shapes for default constructable families.

See also
computeShapes

◆ computeShapes() [3/3]

template<class vec_t >
template<typename OperatorTuple , typename approx_t >
RaggedShapeStorage<vec_t, OperatorTuple> mm::DomainDiscretization< vec_t >::computeShapes ( OperatorTuple  operators,
approx_t  approx,
const indexes_t indexes = {} 
) const

Compute shapes for non-default constructable families.

See also
computeShapes

◆ contains() [1/2]

template<class vec_t >
bool mm::DomainDiscretization< vec_t >::contains ( const vec_t &  point) const

Returns true if point is inside the domain.

Exceptions
Assertionfails if domain shape does not have a contains method. In that case, discreteContains might be more appropriate.
See also
DomainShape::hasContains, DomainShape::contains, DomainDiscretization::discreteContains

Definition at line 315 of file DomainDiscretization.hpp.

◆ contains() [2/2]

template<class vec_t >
template<typename search_structure_t >
bool mm::DomainDiscretization< vec_t >::contains ( const vec_t &  point,
const search_structure_t &  search 
) const

Returns true if point is inside the domain.

If possible, this function uses analytical DomainShape::contains, otherwise, it falls back to discreteContains. For parameter explanations see discreteContains.

See also
contains, discreteContains, DomainShape::hasContains, DomainShape::contains

Definition at line 335 of file DomainDiscretization.hpp.

◆ discreteContains()

template<class vec_t >
template<typename search_structure_t >
bool mm::DomainDiscretization< vec_t >::discreteContains ( const vec_t &  point,
const search_structure_t &  search 
) const

A discrete version of contains using given boundary discretization.

A point p is defined to be contained in the domain, iff the dot product (p-q).n is negative, where q is the nearest boundary point and n is its outer normal. This works sufficiently well, but can fail if boundary discretization is too sparse.

Parameters
pointThe point being tested.
searchNearest neighbour search structure containing boundary points in the same order as defined by boundary_map_. Such suitable structure is returned by makeDiscreteContainsStructure.

Usage example:

BoxShape<Vec3d> box({0, 0, 0}, {2, 2, 2});
BallShape<Vec3d> circ({0, 0, 0}, 1);
auto shape = box - circ;
auto discr = shape.discretizeBoundaryWithStep(0.1);
KDTree<Vec3d> search;
discr.makeDiscreteContainsStructure(search);
EXPECT_TRUE(discr.discreteContains({1.2, 1.2, 1.2}, search));
EXPECT_FALSE(discr.discreteContains({0, 0, 0}, search));
// discreteContains might fail near the shape boundary
EXPECT_FALSE(discr.discreteContains({2, 2, 2}, search));
See also
contains, makeDiscreteContainsStructure

Definition at line 323 of file DomainDiscretization.hpp.

◆ dr()

template<class vec_t >
scalar_t mm::DomainDiscretization< vec_t >::dr ( int  i) const
inline

Returns Euclidean distance to the second support node.

Definition at line 148 of file DomainDiscretization_fwd.hpp.

◆ fill() [1/2]

template<class vec_t >
template<typename callable_t , typename... Args>
void mm::DomainDiscretization< vec_t >::fill ( callable_t &  callable,
Args &&...  args 
)
inline

Enables more readable calls to fill engines.

Definition at line 459 of file DomainDiscretization_fwd.hpp.

◆ fill() [2/2]

template<class vec_t >
template<typename callable_t , typename... Args>
void mm::DomainDiscretization< vec_t >::fill ( const callable_t &  callable,
Args &&...  args 
)
inline

Const version of DomainDiscretization::fill.

Definition at line 461 of file DomainDiscretization_fwd.hpp.

◆ findSupport() [1/2]

template<class vec_t >
template<typename callable_t , typename... Args>
void mm::DomainDiscretization< vec_t >::findSupport ( callable_t &  callable,
Args &&...  args 
)
inline

Enables more readable calls to support engines.

Definition at line 455 of file DomainDiscretization_fwd.hpp.

◆ findSupport() [2/2]

template<class vec_t >
template<typename callable_t , typename... Args>
void mm::DomainDiscretization< vec_t >::findSupport ( const callable_t &  callable,
Args &&...  args 
)
inline

Const version of DomainDiscretization::findSupport.

Definition at line 457 of file DomainDiscretization_fwd.hpp.

◆ interior()

template<class vec_t >
indexes_t mm::DomainDiscretization< vec_t >::interior ( ) const
inline

Returns indexes of all internal nodes.

Definition at line 184 of file DomainDiscretization_fwd.hpp.

◆ load()

template<class vec_t >
template<typename hdf5_file >
DomainDiscretization< vec_t > mm::DomainDiscretization< vec_t >::load ( hdf5_file &  file,
const std::string &  name 
)
static

Load discretization from a HD5 file.

This function looks for datasets named pos, types, normals and bmap. pos should be a dim x N dataset of doubles types and bmap should be a N datasets of ints normals should be a dim x B dataset of doubles, where B is the number of boundary nodes. This format is produced by the HDF::writeDomain function.

Usage example:

// Load domain from file `test_domain_load.h5`. It is stored in the `domain` group.
// Shape of `d` is UnknownShape.
HDF file("test/testdata/test_domain_load.h5", HDF::READONLY);
auto d = DomainDiscretization<Vec2d>::load(file, "domain");
Parameters
fileHDF object containing the discretization.
nameName of the folder where discretization is stored.
See also
HDF::writeDomain
Examples
test/domains/DomainDiscretization_test.cpp.

Definition at line 354 of file DomainDiscretization.hpp.

◆ makeDiscreteContainsStructure()

template<class vec_t >
template<typename search_structure_t >
void mm::DomainDiscretization< vec_t >::makeDiscreteContainsStructure ( search_structure_t &  search) const

Fills the search structure search with boundary points for use with contains() or discreteContains().

Parameters
searchA search structure to be used. Its contents are overwritten.
See also
discreteContains

Definition at line 342 of file DomainDiscretization.hpp.

◆ normal() [1/2]

template<class vec_t >
vec_t & mm::DomainDiscretization< vec_t >::normal ( int  i)

Returns writable outside unit normal of i-th node.

See also
normal

Definition at line 31 of file DomainDiscretization.hpp.

◆ normal() [2/2]

template<class vec_t >
const vec_t & mm::DomainDiscretization< vec_t >::normal ( int  i) const

Returns outside unit normal of i-th node.

The node must be a boundary node.

Exceptions
Assertionfails if the noe is not a boundary node, i.e. type(i) < 0 must hold.

Definition at line 40 of file DomainDiscretization.hpp.

◆ normals()

template<class vec_t >
const Range<vec_t>& mm::DomainDiscretization< vec_t >::normals ( ) const
inline

Returns normals of all boundary nodes.

Definition at line 172 of file DomainDiscretization_fwd.hpp.

◆ operator+=()

template<class vec_t >
DomainDiscretization& mm::DomainDiscretization< vec_t >::operator+= ( const DomainDiscretization< vec_t > &  d)
inline

Operator version of DomainDiscretization::add.

See also
add

Definition at line 331 of file DomainDiscretization_fwd.hpp.

◆ operator-=()

template<class vec_t >
DomainDiscretization& mm::DomainDiscretization< vec_t >::operator-= ( const DomainDiscretization< vec_t > &  d)
inline

Operator version of DomainDiscretization::subtract.

See also
subtract

Definition at line 361 of file DomainDiscretization_fwd.hpp.

◆ output_report()

template<class vec_t >
std::ostream & mm::DomainDiscretization< vec_t >::output_report ( std::ostream &  os = std::cout) const
private

Outputs a simple report about out domain, like number of nodes, support size.

Definition at line 209 of file DomainDiscretization.hpp.

◆ pos() [1/3]

template<class vec_t >
vec_t& mm::DomainDiscretization< vec_t >::pos ( int  i)
inline

Returns writeable position of i-th node.

Definition at line 131 of file DomainDiscretization_fwd.hpp.

◆ pos() [2/3]

template<class vec_t >
const vec_t& mm::DomainDiscretization< vec_t >::pos ( int  i) const
inline

Returns the position of i-th node.

Definition at line 129 of file DomainDiscretization_fwd.hpp.

◆ pos() [3/3]

template<class vec_t >
scalar_t mm::DomainDiscretization< vec_t >::pos ( int  i,
int  j 
) const
inline

Returns j-th coordinate of the position of i-th node.

Definition at line 133 of file DomainDiscretization_fwd.hpp.

◆ positions()

template<class vec_t >
const Range<vec_t>& mm::DomainDiscretization< vec_t >::positions ( ) const
inline

Returns positions of all nodes.

Definition at line 127 of file DomainDiscretization_fwd.hpp.

◆ relax() [1/2]

template<class vec_t >
template<typename callable_t , typename... Args>
void mm::DomainDiscretization< vec_t >::relax ( callable_t &  callable,
Args &&...  args 
)
inline

Enables more readable calls to relax engines.

Definition at line 463 of file DomainDiscretization_fwd.hpp.

◆ relax() [2/2]

template<class vec_t >
template<typename callable_t , typename... Args>
void mm::DomainDiscretization< vec_t >::relax ( const callable_t &  callable,
Args &&...  args 
)
inline

Const version of DomainDiscretization::relax.

Definition at line 465 of file DomainDiscretization_fwd.hpp.

◆ removeBoundaryNodes()

template<class vec_t >
void mm::DomainDiscretization< vec_t >::removeBoundaryNodes ( )
inline

Removes all boundary nodes.

Definition at line 391 of file DomainDiscretization_fwd.hpp.

◆ removeInternalNodes()

template<class vec_t >
void mm::DomainDiscretization< vec_t >::removeInternalNodes ( )
inline

Removes all internal nodes.

Definition at line 388 of file DomainDiscretization_fwd.hpp.

◆ removeNodes()

template<class vec_t >
void mm::DomainDiscretization< vec_t >::removeNodes ( const Range< int > &  to_remove)

Remove nodes with given indices. This removes their types and potential normals as well.

Definition at line 117 of file DomainDiscretization.hpp.

◆ reorderNodes()

template<class vec_t >
template<typename compare_t >
template mm::Range< int > mm::DomainDiscretization< vec_t >::reorderNodes ( const compare_t &  cmp = compare_t())

Reorders nodes according to the compare function.

This can be useful to make the matrix of the implicit system have a nicer structure. Sorts lexicographically (i.e. by x-coordinate first) by default.

Parameters
cmpA compare function for elements of type std::pair<vec_t, int> as passed to std::sort. The first value in the pair is the coordinate of the point and the second is its index.
Returns
The permutation used to shuffle the nodes.
See also
shuffleNodes

Definition at line 441 of file DomainDiscretization.hpp.

◆ rotate() [1/2]

template<class vec_t >
DomainDiscretization< vec_t > & mm::DomainDiscretization< vec_t >::rotate ( const Eigen::Matrix< scalar_t, dim, dim > &  Q)

Transform the domain by given orthogonal matrix Q.

Definition at line 420 of file DomainDiscretization.hpp.

◆ rotate() [2/2]

template<class vec_t >
DomainDiscretization< vec_t > & mm::DomainDiscretization< vec_t >::rotate ( scalar_t  angle)

2D version of rotate accepting an angle.

Definition at line 431 of file DomainDiscretization.hpp.

◆ shape()

template<class vec_t >
const DomainShape<vec_t>& mm::DomainDiscretization< vec_t >::shape ( ) const
inline

Returns geometric shape of the underlying domain.

Definition at line 162 of file DomainDiscretization_fwd.hpp.

◆ shuffleNodes()

template<class vec_t >
void mm::DomainDiscretization< vec_t >::shuffleNodes ( const indexes_t permutation)

Shuffles node order according to given permutation.

This also modifies the support indices and all other necessary data. For given list I = {4, 3, 2, 0, 1} and original positions P = {a, b, c, d, e}, the new list of positions would be {e, d, c, a, b}, i.e. P = P(I). This means that old_domain.pos(p[i]) == new_domain.pos(i) holds for all nodes.

Parameters
permutationA list of n integers containing each number from 0 to n-1 exactly once.
Exceptions
Assertionfails if permutation is not a permutation.

Definition at line 453 of file DomainDiscretization.hpp.

◆ size()

template<class vec_t >
int mm::DomainDiscretization< vec_t >::size ( ) const
inline

Returns N, the number of nodes in this discretization.

Definition at line 189 of file DomainDiscretization_fwd.hpp.

◆ subtract()

template<class vec_t >
DomainDiscretization< vec_t > & mm::DomainDiscretization< vec_t >::subtract ( const DomainDiscretization< vec_t > &  d)

Subtracts given discretized domain from *this.

The shape of current domain changes to become the difference of shape() and d.shape(). The discretization nodes that are now outside of domain are removed and the appropriate boundary nodes from d are added. Types of nodes are preserved. This makes it possible to set custom types before the -= operation and refer to those types in the difference as well.

Example:

double h = 0.1;
DomainDiscretization<Vec2d> rect = BoxShape<Vec2d>({-1, -1}, {1, 1}).discretizeWithStep(h);
DomainDiscretization<Vec2d> circ = BallShape<Vec2d>({0, 0}, 0.7).discretizeBoundaryWithStep(h);
circ.types()[circ.boundary()] = -5; // -1 to -4 taken by rectangle sides
rect -= circ;
std::cout << (rect.types() == -5).size() << std::endl; // nonzero number is printed

More specifically, when performing d1 -= d2, nodes in d1 that are contained in d2 or less than dx away from d2 are removed. The local spacing dx for a point p in d1 is computed by first finding the nearest boundary node in d2 and then taking a fraction of the distance to its nearest boundary node. Boundary nodes from d2 that are contained in d1 are also added, with inverted normals.

See the images below for visualization of the subtract operation.

See also
ShapeDifference, add

Definition at line 298 of file DomainDiscretization.hpp.

◆ support() [1/3]

template<class vec_t >
Range<int>& mm::DomainDiscretization< vec_t >::support ( int  i)
inline

Returns writeable array of support indices for i-th node.

Definition at line 140 of file DomainDiscretization_fwd.hpp.

◆ support() [2/3]

template<class vec_t >
const Range<int>& mm::DomainDiscretization< vec_t >::support ( int  i) const
inline

Returns support indices for i-th node.

Indices are ordered by distance to support nodes, with the first being the closest, usually the node itself.

Definition at line 138 of file DomainDiscretization_fwd.hpp.

◆ support() [3/3]

template<class vec_t >
int mm::DomainDiscretization< vec_t >::support ( int  i,
int  j 
) const
inline

Returns j-th support node of i-th node.

Definition at line 142 of file DomainDiscretization_fwd.hpp.

◆ supportNode()

template<class vec_t >
vec_t mm::DomainDiscretization< vec_t >::supportNode ( int  i,
int  j 
) const
inline

Returns position of j-th support node of i-th node.

Definition at line 146 of file DomainDiscretization_fwd.hpp.

◆ supportNodes()

template<class vec_t >
Range<vec_t> mm::DomainDiscretization< vec_t >::supportNodes ( int  i) const
inline

Returns positions of support nodes of i-th node.

Definition at line 144 of file DomainDiscretization_fwd.hpp.

◆ supports()

template<class vec_t >
const Range<Range<int> >& mm::DomainDiscretization< vec_t >::supports ( ) const
inline

Returns support indices for all nodes.

Definition at line 135 of file DomainDiscretization_fwd.hpp.

◆ supportSize()

template<class vec_t >
int mm::DomainDiscretization< vec_t >::supportSize ( int  i) const
inline

Returns size of i-th node support.

Definition at line 150 of file DomainDiscretization_fwd.hpp.

◆ supportSizes()

template<class vec_t >
Range< int > mm::DomainDiscretization< vec_t >::supportSizes

Returns a vector of support sizes for each node.

Definition at line 21 of file DomainDiscretization.hpp.

◆ translate()

template<class vec_t >
DomainDiscretization< vec_t > & mm::DomainDiscretization< vec_t >::translate ( const vec_t &  a)

Translate the domain by a given vector a.

Definition at line 411 of file DomainDiscretization.hpp.

◆ type() [1/2]

template<class vec_t >
int& mm::DomainDiscretization< vec_t >::type ( int  i)
inline

Returns writeable type of i-th node.

Definition at line 160 of file DomainDiscretization_fwd.hpp.

◆ type() [2/2]

template<class vec_t >
int mm::DomainDiscretization< vec_t >::type ( int  i) const
inline

Returns type of i-th node.

Definition at line 158 of file DomainDiscretization_fwd.hpp.

◆ types() [1/2]

template<class vec_t >
Range<int>& mm::DomainDiscretization< vec_t >::types ( )
inline

Returns mutable types of all nodes.

Definition at line 156 of file DomainDiscretization_fwd.hpp.

◆ types() [2/2]

template<class vec_t >
const Range<int>& mm::DomainDiscretization< vec_t >::types ( ) const
inline

Returns types of all nodes.

Definition at line 154 of file DomainDiscretization_fwd.hpp.

Friends And Related Function Documentation

◆ operator<<

template<class vec_t >
template<class V >
std::ostream& operator<< ( std::ostream &  os,
const DomainDiscretization< V > &  d 
)
friend

Output basic info about given domain.

See also
computeShapes

Member Data Documentation

◆ boundary_map_

template<class vec_t >
Range<int> mm::DomainDiscretization< vec_t >::boundary_map_
protected

Mapping index of a boundary node among all nodes, to its index among boundary nodes.

Indices of internal nodes are set to -1.

positions_ = {p1, p2, p3, p4, p5};
types_ = { 2, -1, 2, 1, -1};
boundary_map_ = {-1, 0, -1, -1, 1};

Nodes p2 and p5 are the first and the second boundary nodes. There are no guarantees on the order of nodes.

Definition at line 85 of file DomainDiscretization_fwd.hpp.

◆ normals_

template<class vec_t >
Range<vec_t> mm::DomainDiscretization< vec_t >::normals_
protected

List of normals of boundary nodes.

positions_ = {p1, p2, p3, p4, p5};
types_ = { 2, -1, 2, 1, -1};
boundary_map_ = {-1, 0, -1, -1, 1};
normals = {n1, n2};

Normal n1 is the normal of node p2 and normal n2 is the normal of p5.

Definition at line 97 of file DomainDiscretization_fwd.hpp.

◆ positions_

template<class vec_t >
Range<vec_t> mm::DomainDiscretization< vec_t >::positions_
protected

Positions of internal discretization points.

Definition at line 54 of file DomainDiscretization_fwd.hpp.

◆ shape_

template<class vec_t >
deep_copy_unique_ptr<DomainShape<vec_t> > mm::DomainDiscretization< vec_t >::shape_
protected

Geometric shape of the domain.

Definition at line 100 of file DomainDiscretization_fwd.hpp.

◆ support_

template<class vec_t >
Range<Range<int> > mm::DomainDiscretization< vec_t >::support_
protected

Attribute used to store support points of each node.

For each node i, support[i] is list of indices of i's support ordered by distance from i. The first element in the i-th list is usually i followed by one or more others.

Definition at line 60 of file DomainDiscretization_fwd.hpp.

◆ types_

template<class vec_t >
Range<int> mm::DomainDiscretization< vec_t >::types_
protected

Storing types of nodes aligned with positions.

Negative types are reserved for boundary nodes, positive for interior. Zero type is not used. Example:

positions_ = {p1, p2, p3, p4, p5};
types_ = { 2, -1, 2, 1, -1};

Nodes p2 and p5 are of type -1, p1 and p3 of type 2 and so on.

Definition at line 72 of file DomainDiscretization_fwd.hpp.


The documentation for this class was generated from the following files:
mm::DomainDiscretization::boundary_map_
Range< int > boundary_map_
Mapping index of a boundary node among all nodes, to its index among boundary nodes.
Definition: DomainDiscretization_fwd.hpp:85
mm::DomainDiscretization::types_
Range< int > types_
Storing types of nodes aligned with positions.
Definition: DomainDiscretization_fwd.hpp:72
mm::DomainDiscretization::normals
const Range< vec_t > & normals() const
Returns normals of all boundary nodes.
Definition: DomainDiscretization_fwd.hpp:172
mm::HDF::READONLY
@ READONLY
Read only open, the file must exist.
Definition: HDF_fwd.hpp:87
mm::sh::d2
static const shape_flags d2
Indicates to calculate d2 shapes.
Definition: shape_flags.hpp:25
mm::DomainDiscretization::size
int size() const
Returns N, the number of nodes in this discretization.
Definition: DomainDiscretization_fwd.hpp:189
mm::DomainDiscretization::shape
const DomainShape< vec_t > & shape() const
Returns geometric shape of the underlying domain.
Definition: DomainDiscretization_fwd.hpp:162
mm::ScalarFieldd
ScalarField< double > ScalarFieldd
Convenience typedef for ScalarField of doubles.
Definition: ScalarField_fwd.hpp:72
mm::Vec3d
Vec< double, 3 > Vec3d
Convenience typedef for 3d vector of doubles.
Definition: Vec_fwd.hpp:35
mm::DomainDiscretization::fill
void fill(callable_t &callable, Args &&... args)
Enables more readable calls to fill engines.
Definition: DomainDiscretization_fwd.hpp:459
mm::DomainDiscretization::load
static DomainDiscretization< vec_t > load(hdf5_file &file, const std::string &name)
Load discretization from a HD5 file.
Definition: DomainDiscretization.hpp:354
mm::DomainDiscretization::positions_
Range< vec_t > positions_
Positions of internal discretization points.
Definition: DomainDiscretization_fwd.hpp:54