MWNode¶
This is an introduction to the mwtree class. We write a small overarching summary of the class where we define the algorithm/equation/structure reasoning for having this class or where it fits with the rest of the code.
digraph { "MWNode" -> "FunctionNode" "MWNode" -> "OperatorNode" }-
template<int D>
class MWNode¶ Base class for Multiwavelet nodes.
A MWNode will contain the scaling and wavelet coefficients to represent functions or operators within a Multiwavelet framework. The nodes are in multidimensional. The dimensionality is set thoucgh the template parameter D=1,2,3. In addition to the coefficients the node contains metadata such as the scale, the translation index, the norm, pointers to parent node and child nodes, pointer to the corresponding MWTree etc… See memeber and data descriptions for details.
Subclassed by mrcpp::FunctionNode< D >
Public Functions
-
MWNode(const MWNode<D> &node, bool allocCoef = true)¶
MWNode copy constructor.
Creates loose nodes and optionally copy coefs. The node does not “belong” to the tree: it cannot be accessed by traversing the tree.
- Parameters:
node – [in] the original node
allocCoef – [in] if true MW coefficients are allocated and copied from the original node
-
bool hasCoord(const Coord<D> &r) const¶
Test if a given coordinate is within the boundaries of the node.
- Parameters:
r – [in] point coordinates
-
bool isCompatible(const MWNode<D> &node)¶
Testing if nodes are compatible wrt NodeIndex and Tree (order, rootScale, relPrec, etc).
-
bool isAncestor(const NodeIndex<D> &idx) const¶
Test if the node is decending from a given NodeIndex, that is, if they have overlapping support.
- Parameters:
idx – [in] the NodeIndex of the requested node
-
double getScalingNorm() const¶
Calculate and return the squared scaling norm.
-
virtual double getWaveletNorm() const¶
Calculate and return the squared wavelet norm.
-
void getCoefs(Eigen::VectorXd &c) const¶
wraps the MW coefficients into an eigen vector object
-
void printCoefs() const¶
Printout of node coefficients.
-
void getPrimitiveQuadPts(Eigen::MatrixXd &pts) const¶
Returns the quadrature points in a given node.
The original quadrature points are fetched and then dilated and translated. For each cartesian direction \( \alpha = x,y,z... \) the set of quadrature points becomes \( x^\alpha_i = 2^{-n} (x_i + l^\alpha \). By taking all possible \((k+1)^d\) combinations, they will then define a d-dimensional grid of quadrature points.
- Parameters:
pts – [inout] quadrature points in a \( d \times (k+1) \) matrix form.
-
void getPrimitiveChildPts(Eigen::MatrixXd &pts) const¶
Returns the quadrature points in a given node.
The original quadrature points are fetched and then dilated and translated to match the quadrature points in the children of the given node. For each cartesian direction \( \alpha = x,y,z... \) the set of quadrature points becomes \( x^\alpha_i = 2^{-n-1} (x_i + 2 l^\alpha + t^\alpha) \), where \( t^\alpha = 0,1 \). By taking all possible \((k+1)^d\combinations \), they will then define a d-dimensional grid of quadrature points for the child nodes.
- Parameters:
pts – [inout] quadrature points in a \( d \times (k+1) \) matrix form.
-
void getExpandedQuadPts(Eigen::MatrixXd &pts) const¶
Returns the quadrature points in a given node.
The primitive quadrature points are used to obtain a tensor-product representation collecting all \( (k+1)^d \) vectors of quadrature points.
- Parameters:
pts – [inout] expanded quadrature points in a \( d \times (k+1)^d \) matrix form.
-
void getExpandedChildPts(Eigen::MatrixXd &pts) const¶
Returns the quadrature points in a given node.
The primitive quadrature points of the children are used to obtain a tensor-product representation collecting all \( 2^d (k+1)^d \) vectors of quadrature points.
- Parameters:
pts – [inout] expanded quadrature points in a \( d \times 2^d(k+1)^d \) matrix form.
-
void zeroCoefs()¶
sets all MW coefficients and the norms to zero
-
void setCoefBlock(int block, int block_size, const double *c)¶
assigns values to a block of coefficients
a block is typically containing one kind of coefficients (given scaling/wavelet in each direction). Its size is then \( (k+1)^D \) and the index is between 0 and \( 2^D-1 \).
- Parameters:
c – [in] the input coefficients
block – [in] the block index
block_size – [in] size of the block
-
void addCoefBlock(int block, int block_size, const double *c)¶
adds values to a block of coefficients
a block is typically containing one kind of coefficients (given scaling/wavelet in each direction). Its size is then \( (k+1)^D \) and the index is between 0 and \( 2^D-1 \).
- Parameters:
c – [in] the input coefficients
block – [in] the block index
block_size – [in] size of the block
-
void zeroCoefBlock(int block, int block_size)¶
sets values of a block of coefficients to zero
a block is typically containing one kind of coefficients (given scaling/wavelet in each direction). Its size is then \( (k+1)^D \) and the index is between 0 and \( 2^D-1 \).
- Parameters:
block – [in] the block index
block_size – [in] size of the block
-
void attachCoefs(double *coefs)¶
Attach a set of coefs to this node. Only used locally (the tree is not aware of this).
-
void calcNorms()¶
Calculate and store square norm and component norms, if allocated.
-
void zeroNorms()¶
Set all norms to zero.
-
void clearNorms()¶
Set all norms to Undefined.
-
virtual void deleteChildren()¶
Recursive deallocation of children and all their descendants.
Leaves node as LeafNode and children[] as null pointer.
-
virtual void deleteParent()¶
Recursive deallocation of parent and all their forefathers.
-
virtual void cvTransform(int kind)¶
Coefficient-Value transform.
This routine transforms the scaling coefficients of the node to the function values in the corresponding quadrature roots (of its children).
NOTE: this routine assumes a 0/1 (scaling on child 0 and 1) representation, instead of s/d (scaling and wavelet).
- Parameters:
operation – [in] forward (coef->value) or backward (value->coef).
-
virtual void mwTransform(int kind)¶
Multiwavelet transform.
Application of the filters on one node to pass from a 0/1 (scaling on child 0 and 1) representation to an s/d (scaling and wavelet) representation. Bit manipulation is used in order to determine the correct filters and whether to apply them or just pass to the next couple of indexes. The starting coefficients are preserved until the application is terminated, then they are overwritten. With minor modifications this code can also be used for the inverse mw transform (just use the transpose filters) or for the application of an operator (using A, B, C and T parts of an operator instead of G1, G0, H1, H0). This is the version where the three directions are operated one after the other. Although this is formally faster than the other algorithm, the separation of the three dimensions prevent the possibility to use the norm of the operator in order to discard a priori negligible contributions.
- Parameters:
operation – [in] compression (s0,s1->s,d) or reconstruction (s,d->s0,s1).
-
double getNodeNorm(const NodeIndex<D> &idx) const¶
Gives the norm (absolute value) of the node at the given NodeIndex.
Recursive routine to find the node with a given NodeIndex. When an EndNode is found, do not generate any new node, but rather give the value of the norm assuming the function is uniformly distributed within the node.
- Parameters:
idx – [in] the NodeIndex of the requested node
Protected Functions
-
MWNode()¶
MWNode default constructor.
Should be used only by NodeAllocator to obtain virtual table pointers for the derived classes.
-
MWNode(MWTree<D> *tree, int rIdx)¶
MWNode constructor.
Constructor for root nodes. It requires the corresponding MWTree and an integer to fetch the right NodeIndex
- Parameters:
tree – [in] the MWTree the root node belongs to
rIdx – [in] the integer specifying the corresponding root node
-
MWNode(MWTree<D> *tree, const NodeIndex<D> &idx)¶
MWNode constructor.
Constructor for an empty node, given the corresponding MWTree and NodeIndex
- Parameters:
tree – [in] the MWTree the root node belongs to
idx – [in] the NodeIndex defining scale and translation of the node
-
MWNode(MWNode<D> *parent, int cIdx)¶
MWNode constructor.
Constructor for leaf nodes. It requires the corresponding parent and an integer to identify the correct child.
- Parameters:
parent – [in] parent node
cIdx – [in] child index of the current node
-
virtual void dealloc()¶
Dummy deallocation of MWNode coefficients.
This is just to make sure this method never really gets called (derived classes must implement their own version). This was to avoid having pure virtual methods in the base class.
-
bool crop(double prec, double splitFac, bool absPrec)¶
Recurse down until an EndNode is found, and then crop children below the given precision threshold.
- Parameters:
prec – [in] precision required
splitFac – [in] factor used in the split check (larger factor means tighter threshold for finer nodes)
absPrec – [in] flag to switch from relative (false) to absolute (true) precision.
-
virtual void allocCoefs(int n_blocks, int block_size)¶
Allocate the coefs vector.
This is only used by loose nodes, because the loose nodes are not treated by the NodeAllocator class.
-
virtual void freeCoefs()¶
Deallocate the coefs vector.
This is only used by loose nodes, because the loose nodes are not treated by the NodeAllocator class.
-
void setMaxSquareNorm()¶
recursively set maxSquaredNorm and maxWSquareNorm of parent and descendants
normalization is such that a constant function gives constant value, i.e. not same normalization as a squareNorm
-
void resetMaxSquareNorm()¶
recursively reset maxSquaredNorm and maxWSquareNorm of parent and descendants to value -1
-
virtual double calcComponentNorm(int i) const¶
Calculate the norm of one component (NOT the squared norm!).
-
virtual void reCompress()¶
Update the coefficients of the node by a mw transform of the scaling coefficients of the children.
-
virtual void giveChildrenCoefs(bool overwrite = true)¶
forward MW transform from this node to its children
it performs forward MW transform inserting the result directly in the right place for each child node. The children must already be present and its memory allocated for this to work properly.
- Parameters:
overwrite – [in] if true the coefficients of the children are overwritten. If false the values are summed to the already present ones.
-
virtual void giveChildCoefs(int cIdx, bool overwrite = true)¶
forward MW transform to compute scaling coefficients of a single child
it performs forward MW transform in place on a loose node. The scaling coefficients of the selected child are then copied/summed in the correct child node.
- Parameters:
cIdx – [in] child index
overwrite – [in] if true the coefficients of the children are overwritten. If false the values are summed to the already present ones.
-
virtual void giveParentCoefs(bool overwrite = true)¶
backward MW transform to compute scaling/wavelet coefficients of a parent
Takes a MWParent and generates coefficients, reverse operation from giveChildrenCoefs
Warning
This routine is only used in connection with Periodic Boundary Conditions
-
virtual void copyCoefsFromChildren()¶
Copy scaling coefficients from children to parent.
Takes the scaling coefficients of the children and stores them consecutively in the corresponding block of the parent, following the usual bitwise notation.
-
int getChildIndex(const NodeIndex<D> &nIdx) const¶
Routine to find the path along the tree.
Given the translation indices at the final scale, computes the child m to be followed at the current scale in oder to get to the requested node at the final scale. The result is the index of the child needed. The index is obtained by bit manipulation of of the translation indices.
- Parameters:
nIdx – [in] the sought after node through its NodeIndex
-
int getChildIndex(const Coord<D> &r) const¶
Routine to find the path along the tree.
@detailsGiven a point in space, determines which child should be followed to get to the corresponding terminal node.
- Parameters:
r – [in] the sought after node through the coordinates of a point in space
-
MWNode<D> *retrieveNode(const Coord<D> &r, int depth)¶
Node retriever that ALWAYS returns the requested node.
Recursive routine to find and return the node with a given NodeIndex. This routine always returns the appropriate node, and will generate nodes that does not exist. Recursion starts at this node and ASSUMES the requested node is in fact decending from this node.
- Parameters:
r – [in] the coordinates of a point in the node
depth – [in] the depth which one needs to descend
-
MWNode<D> *retrieveNode(const NodeIndex<D> &idx)¶
Node retriever that ALWAYS returns the requested node, possibly without coefs.
Recursive routine to find and return the node with a given NodeIndex. This routine always returns the appropriate node, and will generate nodes that does not exist. Recursion starts at this node and ASSUMES the requested node is in fact descending from this node.
- Parameters:
idx – [in] the NodeIndex of the requested node
-
MWNode<D> *retrieveParent(const NodeIndex<D> &idx)¶
Node retriever that ALWAYS returns the requested node.
WARNING: This routine is NOT thread safe! Must be used within omp critical.
Recursive routine to find and return the node with a given NodeIndex. This routine always returns the appropriate node, and will generate nodes that does not exist. Recursion starts at this node and ASSUMES the requested node is in fact related to this node.
- Parameters:
idx – [in] the NodeIndex of the requested node
-
const MWNode<D> *retrieveNodeNoGen(const NodeIndex<D> &idx) const¶
Const version of node retriever that NEVER generates.
Recursive routine to find and return the node with a given NodeIndex. This routine returns the appropriate Node, or a NULL pointer if the node does not exist, or if it is a GenNode. Recursion starts at at this node and ASSUMES the requested node is in fact decending from this node.
- Parameters:
idx – [in] the requested NodeIndex
-
MWNode<D> *retrieveNodeNoGen(const NodeIndex<D> &idx)¶
Node retriever that NEVER generates.
Recursive routine to find and return the node with a given NodeIndex. This routine returns the appropriate Node, or a NULL pointer if the node does not exist, or if it is a GenNode. Recursion starts at at this node and ASSUMES the requested node is in fact decending from this node.
- Parameters:
idx – [in] the requested NodeIndex
-
const MWNode<D> *retrieveNodeOrEndNode(const Coord<D> &r, int depth) const¶
Node retriever that returns requested Node or EndNode (const version).
Recursive routine to find and return the node given the coordinates of a point in space. This routine returns the appropriate Node, or the EndNode on the path to the requested node, and will never create or return GenNodes. Recursion starts at at this node and ASSUMES the requested node is in fact decending from this node.
- Parameters:
r – [in] the coordinates of a point in the node
depth – [in] the depth which one needs to descend
-
MWNode<D> *retrieveNodeOrEndNode(const Coord<D> &r, int depth)¶
Node retriever that returns requested Node or EndNode.
Recursive routine to find and return the node given the coordinates of a point in space. This routine returns the appropriate Node, or the EndNode on the path to the requested node, and will never create or return GenNodes. Recursion starts at at this node and ASSUMES the requested node is in fact decending from this node.
- Parameters:
r – [in] the coordinates of a point in the node
depth – [in] the depth which one needs to descend
-
const MWNode<D> *retrieveNodeOrEndNode(const NodeIndex<D> &idx) const¶
Node retriever that returns requested Node or EndNode (const version).
Recursive routine to find and return the node given the coordinates of a point in space. This routine returns the appropriate Node, or the EndNode on the path to the requested node, and will never create or return GenNodes. Recursion starts at at this node and ASSUMES the requested node is in fact decending from this node.
- Parameters:
idx – [in] the NodeIndex of the requested node
-
MWNode<D> *retrieveNodeOrEndNode(const NodeIndex<D> &idx)¶
Node retriever that returns requested Node or EndNode.
Recursive routine to find and return the node given the coordinates of a point in space. This routine returns the appropriate Node, or the EndNode on the path to the requested node, and will never create or return GenNodes. Recursion starts at at this node and ASSUMES the requested node is in fact decending from this node.
- Parameters:
idx – [in] the NodeIndex of the requested node
-
void threadSafeGenChildren()¶
Generates scaling cofficients of children.
If the node is a leafNode, it takes the scaling&wavelet coefficients of the parent and it generates the scaling coefficients for the children
-
void deleteGenerated()¶
Deallocation of all generated nodes .
-
virtual std::ostream &print(std::ostream &o) const¶
printout ofm the node content.
- Parameters:
o – [in] the output stream
Protected Attributes
-
double squareNorm = {-1.0}¶
Squared norm of all 2^D (k+1)^D coefficients.
-
double maxSquareNorm = {-1.0}¶
Largest squared norm among itself and descendants.
-
double maxWSquareNorm = {-1.0}¶
Largest wavelet squared norm among itself and descendants. NB: must be set before used.
-
double *coefs = {nullptr}¶
the 2^D (k+1)^D MW coefficients For example, in case of a one dimensional function \( f \) this array equals \( s_0, \ldots, s_k, d_0, \ldots, d_k \), where scaling coefficients \( s_j = s_{jl}^n(f) \) and wavelet coefficients \( d_j = d_{jl}^n(f) \). Here \( n, l \) are unique for every node.
-
int serialIx = {-1}¶
index in serial Tree
-
int parentSerialIx = {-1}¶
index of parent in serial Tree, or -1 for roots
-
int childSerialIx = {-1}¶
index of first child in serial Tree, or -1 for leafnodes/endnodes
-
MWNode(const MWNode<D> &node, bool allocCoef = true)¶