OperatorTree¶
This is an introduction to the OperatorTree 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 { "MWTree" -> "FunctionTree" "MWTree" -> "OperatorTree" }-
class OperatorTree : public mrcpp::MWTree<2>¶
Subclassed by mrcpp::CornerOperatorTree
Public Functions
-
virtual void calcBandWidth(double prec = -1.0)¶
Calculates band widths of the non-standard form matrices.
It is starting from \( l = 0 \) and updating the band width value each time we encounter considerable value while keeping increasing \( l \), that stands for the distance to the diagonal.
- Parameters:
prec – [in] Precision used for thresholding
-
virtual bool isOutsideBand(int oTransl, int o_depth, int idx)¶
Checks if the distance to diagonal is bigger than the operator band width.
- Parameters:
oTransl – [in] distance to diagonal
o_depth – [in] scaling order
idx – [in] index corresponding to one of the matrices \( A, B, C \) or \( T \).
- Returns:
True if oTransl is outside of the band and False otherwise.
-
void removeRoughScaleNoise(int trust_scale = 10)¶
Cleans up end nodes.
Traverses the tree and rewrites end nodes having branch node twins, i. e. identical with respect to scale and translation. This method is very handy, when an adaptive operator construction can make a significunt noise at low scaling depth. Its need comes from the fact that mwTransform up cannot override rubbish that can potentially stick to end nodes at a particular level, and as a result spread further up to the root with mwTransform.
- Parameters:
trust_scale – [in] there is no cleaning down below trust_scale (it speeds up operator building).
-
void setupOperNodeCache()¶
Make 1D lists, adressable from [-l, l] scale by scale, of operator node pointers for fast operator retrieval. This method is not thread safe, since it projects missing operator nodes on the fly. Hence, it must NEVER be called within a parallel region, or all hell will break loose. This is not really a problem, but you have been warned.
-
inline OperatorNode &getNode(int n, int l)¶
TODO: It has to be specified more. l is distance to the diagonal.
-
virtual void mwTransformDown(bool overwrite) override¶
Regenerate all scaling coeffs by MW transformation of existing s/w-coeffs on coarser scales, starting at the rootNodes. Option to overwrite or add up existing scaling coefficients (can be used after operator application). Reimplementation of MWTree::mwTransform() without OMP, as calculation of OperatorNorm is done using random vectors, which is non-deterministic in parallel. FunctionTrees should be fine.
-
virtual void mwTransformUp() override¶
Regenerate all s/d-coeffs by backtransformation, starting at the bottom and thus purifying all coefficients. Option to overwrite or add up existing coefficients of BranchNodes (can be used after operator application). Reimplementation of MWTree::mwTransform() without OMP, as calculation of OperatorNorm is done using random vectors, which is non-deterministic in parallel. FunctionTrees should be fine.
-
void setZero()¶
Set the MW coefficients to zero, keeping the same tree structure.
Keeps the node structure of the tree, even though the zero function is representable at depth zero. One should then use cropTree to remove unnecessary nodes.
-
void clear()¶
Remove all nodes in the tree.
Leaves the tree in the same state as after construction, i.e. undefined tree structure containing only root nodes without coefficients. The assigned memory, including branch and leaf nodes, (nodeChunks in NodeAllocator) is NOT released, but is immediately available to the new function.
-
inline double getSquareNorm() const¶
- Returns:
Squared L2 norm of the function
-
void calcSquareNorm(bool deep = false)¶
Calculate the squared norm \( ||f||^2_{\ldots} \) of a function represented as a tree.
The norm is calculated using endNodes only. The specific type of norm which is computed will depend on the derived class
-
inline int getNNodes() const¶
- Returns:
the total number of nodes in the tree
-
int getNNodesAtDepth(int i) const¶
- Parameters:
depth – [in] Tree depth (0 depth is the coarsest scale) to count.
- Returns:
Total number of nodes in the tree, at given depth (not in use)
-
int getSizeNodes() const¶
- Returns:
Size of all MW coefs in the tree, in kB
-
inline NodeBox<D, T> &getRootBox()¶
- Returns:
-
void mwTransform(int type, bool overwrite = true)¶
Full Multiwavelet transform of the tree in either directions.
It performs a Multiwavlet transform of the whole tree. The input parameters will specify the direction (upwards or downwards) and whether the result is added to the coefficients or it overwrites them. See the documentation for the mwTransformUp and mwTransformDown for details.
\[\begin{split} \pmatrix{ s_{nl}\\ d_{nl} } \rightleftarrows \pmatrix{ s_{n+1,2l}\\ s_{n+1,2l+1} } \end{split}\]- Parameters:
type – [in] TopDown (from roots to leaves) or BottomUp (from leaves to roots) which specifies the direction of the MW transform
overwrite – [in] if true, the result will overwrite preexisting coefficients.
-
MWNode<D, T> *findNode(NodeIndex<D> nIdx)¶
Finds and returns the node pointer with the given NodeIndex.
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 the appropriate rootNode.
-
const MWNode<D, T> *findNode(NodeIndex<D> nIdx) const¶
Finds and returns the node pointer with the given NodeIndex, const version.
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 the appropriate rootNode.
-
MWNode<D, T> &getNode(NodeIndex<D> nIdx, bool create = false)¶
Finds and returns the node reference with the given NodeIndex.
This routine ALWAYS returns the node you ask for. If the node does not exist, it will be generated by MW transform. Recursion starts at the appropriate rootNode and descends from this. The nodes are permanently added to the tree if create = true
-
MWNode<D, T> &getNode(Coord<D> r, int depth = -1)¶
Finds and returns the node at a given depth that contains a given coordinate.
This routine ALWAYS returns the node you ask for, and will generate nodes that do not exist. Recursion starts at the appropriate rootNode and decends from this.
- Parameters:
depth – [in] requested node depth from root scale.
r – [in] coordinates of an arbitrary point in space
-
MWNode<D, T> &getNodeOrEndNode(NodeIndex<D> nIdx)¶
Finds and returns the node with the given NodeIndex.
This routine returns the Node you ask for, or the EndNode on the path to the requested node, if the requested one is deeper than the leaf node ancestor. It will never create or return GenNodes. Recursion starts at the appropriate rootNode and decends from this.
-
const MWNode<D, T> &getNodeOrEndNode(NodeIndex<D> nIdx) const¶
Finds and returns the node reference with the given NodeIndex. Const version.
This routine ALWAYS returns the node you ask for. If the node does not exist, it will be generated by MW transform. Recursion starts at the appropriate rootNode and decends from this.
-
MWNode<D, T> &getNodeOrEndNode(Coord<D> r, int depth = -1)¶
Finds and returns the node at a given depth that contains a given coordinate.
This routine returns the Node you ask for, or the EndNode on the path to the requested node, and will never create or return GenNodes. Recursion starts at the appropriate rootNode and decends from this.
- Parameters:
depth – [in] requested node depth from root scale.
r – [in] coordinates of an arbitrary point in space
-
const MWNode<D, T> &getNodeOrEndNode(Coord<D> r, int depth = -1) const¶
Finds and returns the node at a given depth that contains a given coordinate. Const version.
This routine returns the Node you ask for, or the EndNode on the path to the requested node, and will never create or return GenNodes. Recursion starts at the appropriate rootNode and decends from this.
- Parameters:
depth – [in] requested node depth from root scale.
r – [in] coordinates of an arbitrary point in space
-
MWNodeVector<D, T> *copyEndNodeTable()¶
Returns the list of all EndNodes.
copies the list of all EndNode pointers into a new vector and returns it.
-
void deleteRootNodes()¶
Deletes all the nodes in the tree.
This method will recursively delete all the nodes, including the root nodes. Derived classes will call this method when the object is deleted.
-
void resetEndNodeTable()¶
Recreate the endNodeTable.
the endNodeTable is first deleted and then rebuilt from scratch. It makes use of the TreeIterator to traverse the tree.
-
int getIx(NodeIndex<D> nIdx)¶
gives serialIx of a node from its NodeIndex
gives a unique integer for each nodes corresponding to the position of the node in the serialized representation
-
void makeMaxSquareNorms()¶
sets values for maxSquareNorm in all nodes
it defines the upper bound of the squared norm \( ||f||^2_{\ldots} \) in this node or its descendents
Public Members
-
MWNodeVector<D, T> endNodeTable¶
Final projected nodes.
Protected Functions
-
virtual std::ostream &print(std::ostream &o) const override¶
Prints a summary of the tree structure on the output file.
-
void incrementNodeCount(int scale)¶
Increments node counter by one for non-GenNodes.
TO BE DOCUMENTED
Warning
: This routine is not thread safe, and must NEVER be called outside a critical region in parallel. It’s way. way too expensive to lock the tree, so don’t even think about it.
-
void decrementNodeCount(int scale)¶
Decrements node counter by one for non-GenNodes.
TO BE DOCUMENTED
Warning
: This routine is not thread safe, and must NEVER be called outside a critical region in parallel. It’s way. way too expensive to lock the tree, so don’t even think about it.
Protected Attributes
-
OperatorNode ***nodePtrStore¶
Avoids tree lookups.
-
OperatorNode ***nodePtrAccess¶
Center (l=0) of node list.
-
NodeBox<D, T> rootBox¶
The actual container of nodes.
-
std::vector<int> nodesAtDepth¶
Node counter.
-
std::vector<int> nodesAtNegativeDepth¶
Node counter.
-
virtual void calcBandWidth(double prec = -1.0)¶