1 /*
2   HMat-OSS (HMatrix library, open source software)
3 
4   Copyright (C) 2014-2015 Airbus Group SAS
5 
6   This program is free software; you can redistribute it and/or
7   modify it under the terms of the GNU General Public License
8   as published by the Free Software Foundation; either version 2
9   of the License, or (at your option) any later version.
10 
11   This program is distributed in the hope that it will be useful,
12   but WITHOUT ANY WARRANTY; without even the implied warranty of
13   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
14   GNU General Public License for more details.
15 
16   You should have received a copy of the GNU General Public License
17   along with this program; if not, write to the Free Software
18   Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301, USA.
19 
20   http://github.com/jeromerobert/hmat-oss
21 */
22 
23 /*! \file
24   \ingroup HMatrix
25   \brief Implementation of the HMatrix class.
26 */
27 #ifndef _H_MATRIX_HPP
28 #define _H_MATRIX_HPP
29 
30 #include "tree.hpp"
31 #include "assembly.hpp"
32 #include "cluster_tree.hpp"
33 #include "rk_matrix.hpp"
34 
35 namespace hmat {
36     /** Identify the current user level operation */
37     enum class MainOp {OTHER, SOLVELOWER, SOLVEUPPER, GEMM};
38 }
39 
40 #include "recursion.hpp"
41 #include <cassert>
42 #include <fstream>
43 #include <iostream>
44 #include <deque>
45 
46 
47 namespace hmat {
48 
49 template<typename T> class Vector;
50 class AdmissibilityCondition;
51 template<typename T> class FullMatrix;
52 
53 /** Flag used to describe the symmetry of a matrix.
54  */
55 enum SymmetryFlag {kNotSymmetric, kLowerSymmetric};
56 
57 /** Default rank value for blocks that dont have an actual computed rank
58    */
59 // TODO: the meaning/usage of UNINITIALIZED_BLOCK is not clear, it should be reworked
60 // or removed
61 enum DefaultRank {UNINITIALIZED_BLOCK = -3, NONLEAF_BLOCK = -2, FULL_BLOCK = -1};
62 
63 /** Settings global to a whole matrix */
64 struct MatrixSettings {
65 };
66 
67 /** Settings local to a matrix bloc */
68 struct LocalSettings {
69     const MatrixSettings * global;
LocalSettingshmat::LocalSettings70     LocalSettings(const MatrixSettings * s, double epsilon): global(s), epsilon_(epsilon) {}
71     /// epsilon used for SVD truncations
72     double epsilon_;
73 };
74 
75 /** Degrees of freedom permutation of a vector required in HMatrix context.
76 
77      In order that the subsets of rows and columns are
78      contiguous in HMatrix, we must reorder the elements of the vector. This
79      order is induced by the array of indices after the construction of
80      ClusterTree, which must be passed as a parameter.
81 
82      \param v Vector to reorder.
83      \param indices Array of indices after construction ClusterTree.
84 
85  */
86 template<typename T> void reorderVector(ScalarArray<T>* v, int* indices, int axis);
87 
88 /** Inverse permutation of a vector.
89 
90      See \a reorderVector () for more details.
91 
92      \param v Vector to reorder of the problem.
93      \param indices Array of indices after construction ClusterTree.
94 
95  */
96 template<typename T> void restoreVectorOrder(ScalarArray<T>* v, int *indices, int axis);
97 
98 template<typename T> class HMatrix;
99 
100 enum class Axis {ROW, COL};
101 
102 /** Precompute compatibility between children of a and b for GEMM.
103 
104      \param a first matrix
105      \param axisA check rows or columns of a
106      \param transA tells whether a is transposed or not
107      \param b second matrix
108      \param axisB check rows or columns of b
109      \param transB tells whether b is transposed or not
110      \return byte array
111  */
112 template<typename T> unsigned char * compatibilityGridForGEMM(const HMatrix<T>* a, Axis axisA, char transA, const HMatrix<T>* b, Axis axisB, char transB);
113 
114 /*! \brief The HMatrix class, representing a HMatrix.
115 
116   It is a tree of arity arity(ClusterTree)^2, 4 in most cases.
117   An HMatrix is a tree-like structure that is:
118     - a Leaf : in this case the node is either really a RkMatrix
119       (compressed block), or a small dense block.
120     - an internal node : in this case, it has 4 children that form a partition
121       of the HMatrix Dofs, and the node doesn't carry data itself.
122  */
123 template<typename T> class HMatrix : public Tree<HMatrix<T> >, public RecursionMatrix<T, HMatrix<T> > {
124   friend class RkMatrix<T>;
125 
126   /// Rows of this HMatrix block
127   const ClusterTree * rows_;
128   /// Columns of this HMatrix block
129   const ClusterTree * cols_;
130   union {
131    /// Compressed block, or NULL if the block is not a leaf or is full.
132    RkMatrix<T> * rk_;
133    /// Full block, or NULL if the block is not a leaf or is compressed.
134    FullMatrix<T> * full_;
135   };
136   /// rank of the block for Rk matrices, or: UNINITIALIZED_BLOCK=-3 for an uninitialized matrix, NONLEAF_BLOCK=-2 for non leaf, FULL_BLOCK=-1 for full a matrix
137   int rank_;
138   /// approximate rank of the block, or: UNINITIALIZED_BLOCK=-3 for an uninitialized matrix
139   int approximateRank_;
140   void uncompatibleGemm(char transA, char transB, T alpha, const HMatrix<T>* a, const HMatrix<T>*b);
141   void recursiveGemm(char transA, char transB, T alpha, const HMatrix<T>* a, const HMatrix<T>*b);
142   void leafGemm(char transA, char transB, T alpha, const HMatrix<T>* a, const HMatrix<T>*b);
143   HMatrix<T> * fullRkSubset(const IndexSet* subset, bool col) const;
144 
145   /** Only used by internalCopy */
146   HMatrix(const MatrixSettings * settings);
147   /** This <- This + alpha * b
148 
149       \param alpha
150       \param b
151    */
152   void axpy(T alpha, const RkMatrix<T>* b);
153   /** This <- This + alpha * b
154 
155       \param alpha
156       \param b
157       \param rows
158       \param cols
159    */
160   void axpy(T alpha, const FullMatrix<T>* b);
161 public:
162   /*! \brief Create a HMatrix based on a row and column ClusterTree.
163 
164     \param _rows The row cluster tree
165     \param _cols The column cluster tree
166    */
167   HMatrix(const ClusterTree* _rows, const ClusterTree* _cols, const MatrixSettings * settings,
168        int depth, SymmetryFlag symmetryFlag,
169        AdmissibilityCondition * admissibilityCondition);
170 
171   /*! \brief Create a copy of this matrix for internal use only.
172    * Only copy this node, not the whole tree. The created matrix
173    * is an uninitialized leaf with same rows and cols as this.
174    */
175   HMatrix<T> * internalCopy(bool temporary_ = false, bool withRowChild=false, bool withColChild=false) const;
176   ~HMatrix();
177 
178   /** Return a full null block at the given offset and of the give size */
179   HMatrix<T> * internalCopy(const ClusterTree * rows, const ClusterTree * cols) const;
180 
181   /**
182    * \brief Split this block according and admissibility condition
183    *
184    * \return true if the block was actually splitted
185   */
186   bool split(AdmissibilityCondition * admissibilityCondition, bool lowRank,
187              SymmetryFlag symmetryFlag = kNotSymmetric);
188 
189   /** \brief Add the list of leaves to a list */
190   void listAllLeaves(std::deque<const HMatrix<T> *> & out) const;
191   void listAllLeaves(std::deque<HMatrix<T> *> & out);
192 
193   /**
194    * Create a temporary block from a list of children.
195    * @param children The list of children ordered as insertChild and get
196    * methods expect.
197    */
198   HMatrix(const ClusterTree * rows, const ClusterTree * cols, std::vector<HMatrix*> & children);
199 
200   /*! \brief HMatrix coarsening.
201 
202      If all children are Rk leaves, then we try to merge them into a single Rk-leaf.
203      This is done if the memory of the resulting leaf is less than the sum of the initial
204      leaves. Note that this operation could be used hierarchically.
205      \param epsilon the truncate epsilon
206      \param upper the symmetric of 'this', when building a non-sym matrix with a sym content
207      \param force if true the block is kept coarsened even if it's larger
208      \return true if all leaves are rk (i.e. if coarsening was tryed, not if it succeded)
209    */
210   bool coarsen(double epsilon, HMatrix<T>* upper = NULL, bool force=false) ;
211   /*! \brief HMatrix assembly.
212    */
213   void assemble(Assembly<T>& f, const AllocationObserver & = AllocationObserver());
214   /*! \brief HMatrix assembly.
215 
216     \param f the assembly function
217     \param upper the upper part of the matrix. If NULL, it is assumed
218                  that upper=this (that is, the current block is on the diagonal)
219     \param onlyLower if true, only assemble the lower part of the matrix, ie don't copy.
220    */
221   void assembleSymmetric(Assembly<T>& f,
222      HMatrix<T>* upper=NULL, bool onlyLower=false,
223      const AllocationObserver & = AllocationObserver());
224   /*! \brief Evaluate the HMatrix, ie converts it to a full matrix.
225 
226     This conversion does the reorderng of the unknowns such that the resulting
227     matrix can directly be used or compared with a full matrix assembled
228     otherwise.
229 
230     \param result a FullMatrix that has to be preallocated at the same size than
231     this.
232    */
233   void eval(FullMatrix<T>* result, bool renumber = true) const;
234   /*! \brief Evaluate this as a subblock of the larger matrix result.
235 
236     _rows and _cols are the rows and columns of the result matrix. This has to
237     be a subset of _rows and _cols, and it is put at its place inside the result
238     matrix. This function does not do any reodering.
239 
240     \param result Result matrix, of size (_rows->n, _cols->n)
241     \param _rows Rows of the result matrix
242     \param _cols Columns of the result matrix
243    */
244   void evalPart(FullMatrix<T>* result, const IndexSet* _rows, const IndexSet* _cols) const;
245 
246   void info(hmat_info_t &);
247 
248   /** This *= alpha
249 
250       \param alpha scaling factor
251    */
252   void scale(T alpha);
253   /** Compute y <- alpha * op(this) * x + beta * y if side is Side::LEFT or
254               y <- alpha * x * op(this) + beta * y if side is Side::RIGHT
255 
256       The arguments are similar to BLAS GEMV.
257    */
258   void gemv(char trans, T alpha, const FullMatrix<T>* x, T beta, FullMatrix<T>* y, Side side = Side::LEFT) const;
259   /** Compute y <- alpha * op(this) * x + beta * y if side is Side::LEFT or
260               y <- alpha * x * op(this) + beta * y if side is Side::RIGHT
261 
262       The arguments are similar to BLAS GEMV.
263    */
264   void gemv(char trans, T alpha, const ScalarArray<T>* x, T beta, ScalarArray<T>* y, Side side = Side::LEFT) const;
265   /*! \brief this <- alpha * op(A) * op(B) + beta * this
266 
267     \param transA 'N' or 'T', as in BLAS
268     \param transB 'N' or 'T', as in BLAS
269     \param alpha alpha
270     \param a the matrix A
271     \param b the matrix B
272     \param beta beta
273    */
274   void gemm(char transA, char transB, T alpha, const HMatrix<T>* a, const HMatrix<T>*b, T beta, MainOp=MainOp::OTHER);
275   /*! \brief this <- this - M * D * M^T, where 'this' is symmetric (Lower stored),
276       D diagonal
277 
278       \warning D has to be reduced in ldlt form with d->ldltDecomposition() before
279 
280       \param m M
281       \param d D : only the diagonal of this matrix is considered
282    */
283   void mdmtProduct(const HMatrix<T> * m, const HMatrix<T> * d);
284 
285   /*! \brief this <- this - M * D * N^T with D diagonal
286 
287       \warning D has to be reduced in ldlt form with d->ldltDecomposition() before
288 
289       \param m M
290       \param d D : only the diagonal of this matrix is considered
291       \param n N
292    */
293   void mdntProduct(const HMatrix<T>* m, const HMatrix<T>* d, const HMatrix<T>* n);
294 
295   /** Create a matrix filled with 0s, with the same structure as H.
296 
297       \param h the model matrix,
298       \return a 0 matrix with the same structure as H.
299    */
300   static HMatrix<T>* Zero(const HMatrix<T>* h);
301 
302   /**
303    * Used internally for deserialization
304    * @see serialization.hpp
305    */
306   static HMatrix<T> * unmarshall(const MatrixSettings * settings, int rank, int rankApprox, char bitfield, double epsilon);
307 
308   /** Returns a copy of this (with all the structure and data)
309        */
310   HMatrix<T>* copy() const ;
311   /** this <- o (copy)
312 
313       \param o The HMatrix to copy. 'this' must be allready created and have the right structure.
314    */
315   void copy(const HMatrix<T>* o);
316   /** Copy the structure of an HMatrix without copying its content.
317 
318       \return an empty HMatrix (not even the Full leaves are
319       allocated) mirroring the structure of this.
320    */
321   HMatrix<T>* copyStructure() const;
322   /*! \brief Return square of the Frobenius norm of the matrix.
323    */
324   double normSqr() const;
325   /*! \brief Return the Frobenius norm of the matrix.
326    */
norm() const327   double norm() const {
328     return sqrt(normSqr());
329   }
330 
331   /*! \brief Return an approximation of the largest eigenvalue via the power method.
332    */
333   T approximateLargestEigenvalue(int max_iter, double epsilon) const;
334 
335   /*! \brief Return the approximated rank.
336    */
approximateRank() const337   int approximateRank() const {
338     return approximateRank_;
339   }
340 
approximateRank(int a)341   void approximateRank(int a)  {
342     approximateRank_ = a;
343   }
344 
345   /*! \brief Return low rank epsilon
346    */
lowRankEpsilon() const347   double lowRankEpsilon() const {
348     return localSettings.epsilon_;
349   }
350 
351   /** Recursively set low-rank epsilon member
352    */
353   void lowRankEpsilon(double epsilon, bool recursive = true);
354 
355   /** Set a matrix to 0.
356    */
357   void clear();
358   /** Inverse an HMatrix in place.
359 
360       \param tmp temporary HMatrix used in the inversion. If set, it must have
361       the same structure as this. Otherwise, it is allocated at the start of the
362       computation (and will be freed at the end).
363       \param depth The depth, used for pretty printing purposes
364    */
365   void inverse();
366   /*! \brief Transpose the H-matrix in place
367    */
368   void transpose();
369 
370   void conjugate();
371   /**
372    * Swap non diagonal blocks and cluster trees.
373    * Only used internally.
374    */
375   void transposeMeta(bool temporaryOnly=false);
376   /**
377    * Swap Rk or Full blocks around the diagonal
378    * Only used internally.
379    */
380   void transposeData();
381 
382   /*! \brief this <- o^t
383 
384     \param o
385    */
386   void copyAndTranspose(const HMatrix<T>* o);
387 
388   /*! \brief Truncate Rk matrices with respect to their respective epsilon_
389    */
390   void truncate();
391 
392   /*! \brief LU decomposition in place.
393 
394     \warning Do not use. Doesn't work
395    */
396   void luDecomposition(hmat_progress_t * progress);
397   /** \brief LDL^t decomposition in place
398      \warning this has to be created with the flag lower
399      \warning this has to be assembled with assembleSymmetric with onlyLower = true
400    */
401   void ldltDecomposition(hmat_progress_t * progress);
402   void lltDecomposition(hmat_progress_t * progress);
403 
404   /** This <- This + alpha * b
405 
406       \param alpha
407       \param b
408    */
409   void axpy(T alpha, const HMatrix<T>* b);
410   /** This <- This + alpha * Id
411 
412       \param alpha
413    */
414   void addIdentity(T alpha);
415   /** This <- This + A , norm(A) ~ epsilon*norm(this)
416 
417       \param epsilon
418    */
419   void addRand(double epsilon);
420   /*! Return true if this is a full block.
421    */
isFullMatrix() const422   inline bool isFullMatrix() const {
423     return rank_ == FULL_BLOCK && full_ != NULL;
424   }
425   /** Return the full matrix corresponding to the current leaf */
getFullMatrix() const426   FullMatrix<T>* getFullMatrix() const {
427     assert(isFullMatrix());
428     return full_;
429   }
430   /*! Return true if this is a compressed block.
431    */
isRkMatrix() const432   inline bool isRkMatrix() const {
433     return rank_ >= 0;
434   }
435   /*! \brief Multiplication de deux HMatrix dont au moins une est une RkMatrix.
436 
437       Le resultat est alors une RkMatrix.
438 
439       \param transA 'T' ou 'N' selon si A est transposee ou non
440       \param transB 'T' ou 'N' selon si B est transposee ou non
441       \param a A
442       \param b B
443    */
444   static RkMatrix<T>* multiplyRkMatrix(double epsilon, char transA, char transB, const HMatrix<T>* a, const HMatrix<T>* b);
445   /** Multiplication de deux HMatrix dont au moins une est une matrice pleine,
446       et aucune n'est une RkMatrix.
447 
448       Le resultat est alors une matrice pleine.
449   */
450   static FullMatrix<T>* multiplyFullMatrix(char transA, char transB, const HMatrix<T>* a, const HMatrix<T>* b);
451   /*! \brief B <- B*D where B is this
452 
453     \warning D must a have been decomposed by LDLt
454     \param d matrice D
455     \param left run B <- D*B instead of B <- B*D
456     \param inverse run B <- B * D^-1
457   */
458   void multiplyWithDiag(const HMatrix<T>* d, Side side = Side::RIGHT, bool inverse = false) const;
459   /*! \brief Resolution du systeme L X = B, avec this = L, et X = B.
460 
461     \param b la matrice B en entree, et X en sortie.
462    */
463   void solveLowerTriangularLeft(HMatrix<T>* b, Factorization algo, Diag diag, Uplo uplo, MainOp=MainOp::OTHER) const;
464   /*! \brief Resolution du systeme L x = x, avec this = L, et x = b vecteur.
465 
466     B est un vecteur a plusieurs colonnes, donc une FullMatrix.
467 
468     \param b Le vecteur b en entree, et x en sortie.
469    */
470   void solveLowerTriangularLeft(ScalarArray<T>* b, Factorization algo, Diag diag, Uplo uplo) const;
471   void solveLowerTriangularLeft(FullMatrix<T>* b, Factorization algo, Diag diag, Uplo uplo) const;
472   /*! Resolution de X U = B, avec U = this, et X = B.
473 
474     \param b la matrice B en entree, X en sortie
475    */
476   void solveUpperTriangularRight(HMatrix<T>* b, Factorization algo, Diag diag, Uplo uplo) const;
477   /*! Resolution de U X = B, avec U = this, et X = B.
478 
479     \param b la matrice B en entree, X en sortie
480    */
481   void solveUpperTriangularLeft(HMatrix<T>* b, Factorization algo, Diag diag, Uplo uplo, MainOp=MainOp::OTHER) const;
482   /*! Resolution de x U = b, avec U = this, et x = b.
483 
484     \warning b est un vecteur ligne et non colonne.
485 
486     \param b Le vecteur b en entree, x en sortie.
487    */
488   void solveUpperTriangularRight(ScalarArray<T>* b, Factorization algo, Diag diag, Uplo uplo) const;
489   void solveUpperTriangularRight(FullMatrix<T>* b, Factorization algo, Diag diag, Uplo uplo) const;
490   /*! Resolution de U x = b, avec U = this, et x = b.
491     U peut etre en fait L^T ou L est une matrice stockee inferieurement
492     en precisant uplo = true
493 
494     \param b Le vecteur b en entree, x en sortie.
495     \param indice les indices portes par le vecteur
496     \param uplo indique le stockage de la matrice U ou L^T
497   */
498   void solveUpperTriangularLeft(ScalarArray<T>* b, Factorization algo, Diag diag, Uplo uplo) const;
499   void solveUpperTriangularLeft(FullMatrix<T>* b, Factorization algo, Diag diag, Uplo uplo) const;
500   /*! Solve D x = b, in place with D a diagonal matrix.
501 
502      \param b Input: B, Output: X
503    */
504   void solveDiagonal(ScalarArray<T>* b) const;
505   void solveDiagonal(FullMatrix<T>* b) const;
506   /*! Resolution de This * x = b.
507 
508     \warning This doit etre factorisee avec \a HMatrix::luDecomposition() avant.
509    */
510   void solve(ScalarArray<T>* b) const;
511   void solve(FullMatrix<T>* b) const;
512   /*! Resolution de This * X = b.
513 
514     \warning This doit etre factorisee avec \a HMatrix::luDecomposition() avant.
515    */
516   void solve(HMatrix<T>* b, Factorization algo) const;
517 
518   void trsm( char side, char uplo, char trans, char diag, T alpha,
519 	     HMatrix<T>* b ) const;
520   void trsm( char side, char uplo, char trans, char diag, T alpha,
521 	     ScalarArray<T>* b ) const;
522 
523   /*! Resolution de This * x = b.
524 
525     \warning This doit etre factorisee avec \a HMatrix::ldltDecomposition() avant.
526    */
527   void solveLdlt(ScalarArray<T>* b) const ;
528   void solveLdlt(FullMatrix<T>* b) const ;
529   /*! Resolution de This * x = b.
530 
531     \warning This doit etre factorisee avec \a HMatrix::lltDecomposition() avant.
532    */
533   void solveLlt(ScalarArray<T>* b) const ;
534   void solveLlt(FullMatrix<T>* b) const ;
535   /*! Triggers an assertion if the HMatrix contains any NaN.
536    */
537   void checkNan() const;
538   /*! Triggers an assertion if children of an HMatrix are not contained within
539       this HMatrix.
540    */
541   void checkStructure() const;
542   /** Recursively set the isTriLower flag on this matrix diagonal blocks */
543   void setTriLower(bool value);
544   /** Recursively set the isLower flag on this matrix diagonal blocks */
545   void setLower(bool value);
546 
547   const ClusterData* rows() const;
548   const ClusterData* cols() const;
549 
550   /*! \brief Return the number of children in the row dimension.
551     */
nrChildRow() const552   inline int nrChildRow() const {
553     // if rows admissible, only one child = itself
554     return keepSameRows ? 1 : rows_->nrChild() ;
555   }
556 
557   /*! \brief Return the number of children in the column dimension.
558     */
nrChildCol() const559   inline int nrChildCol() const {
560     // if cols admissible, only one child = itself
561     return keepSameCols ? 1 : cols_->nrChild() ;
562   }
563   /*! \brief Destroy the HMatrix.
564     */
destroy()565   void destroy() {
566     delete this;
567   }
568 
569   /*! Return the child (i, j) of this.
570 
571     \warning do not use on a leaf !
572 
573     \param i row
574     \param j column
575     \return the (i,j) child of this.
576    */
get(int i,int j) const577   HMatrix<T>* get(int i, int j) const {
578     assert(i>=0 && i<nrChildRow());
579     assert(j>=0 && j<nrChildCol());
580     assert(i + j * nrChildRow() < this->nrChild());
581     return this->getChild(i + j * nrChildRow());
582   }
583 
584   /*! Set the child (i, j) of this.
585 
586     \warning do not use on a leaf ! <- how can I know, since I am inserting a child ???
587 
588     \param i row index
589     \param j column index
590     \param child the child (i, j) of this.
591    */
592   using Tree<HMatrix<T> >::insertChild;
insertChild(int i,int j,HMatrix<T> * child)593   void insertChild(int i, int j, HMatrix<T>* child) {
594     insertChild(i+j*nrChildRow(), child) ;
595   }
596 
597   /*! \brief Find the correct child when recursing in GEMM or GEMV
598 
599     This function returns the child (i,j) of op(this) where 'op' is 'T' or 'N'.
600     If the matrix is symmetric (upper or lower), and if the child required is in the
601     part of the matrix that is not stored, it returns the symmetric block and switches 'op'.
602    \param[in,out] t input is the transpose flag for this, ouput is the transpose flag for the returned matrix
603    \param[in] i row index of the child
604    \param[in] j col index of the child
605    \return Pointer on the child
606   */
607   const HMatrix<T> * getChildForGEMM(char & t, int i, int j) const;
608 
609   void setClusterTrees(const ClusterTree* rows, const ClusterTree* cols);
610   HMatrix<T> * subset(const IndexSet * rows, const IndexSet * cols) const;
611 
612   /* \brief Retrieve diagonal values.
613   */
614   void extractDiagonal(T* diag) const;
615 
616   /// Should try to coarsen the matrix at assembly
617   static bool coarsening;
618   /// Should recompress the matrix after assembly
619   static bool recompress;//TODO: remove
620   /// Validate the functions is_guaranteed_null_col/row() (user provided)
621   static bool validateNullRowCol;
622   /// Validate the rk-matrices after compression
623   static bool validateCompression;
624   /// For blocks above error threshold, re-run the compression algorithm
625   static bool validationReRun;
626   /// For blocks above error threshold, dump the faulty block to disk
627   static bool validationDump;
628   /// Error threshold for the compression validation
629   static double validationErrorThreshold;
630   short isUpper:1, isLower:1,       /// symmetric, upper or lower stored
631        isTriUpper:1, isTriLower:1, /// upper/lower triangular
632        keepSameRows:1, keepSameCols:1,
633        temporary_:1, ownRowsClusterTree_:1, ownColsClusterTree_:1;
634   LocalSettings localSettings;
635 
rank() const636   int rank() const {
637       assert(rank_ >= 0);
638       return rank_;
639   }
640 
641   /// Set the rank of an evicted rk block
642   void rank(int rank);
643 
rk() const644   RkMatrix<T> * rk() const {
645       assert(rank_ >= 0);
646       return rk_;
647   }
648 
649   /*! \brief Set 'this' as an Rk matrix using copy of a and b
650      */
651   void rk(const ScalarArray<T> *a, const ScalarArray<T> *b);
652 
rk(RkMatrix<T> * m)653   void rk(RkMatrix<T> * m) {
654       rk_ = m;
655       rank_ = m == NULL ? 0 : m->rank();
656   }
657 
full() const658   FullMatrix<T> * full() const {
659       assert(rank_ == FULL_BLOCK);
660       return full_;
661   }
662 
full(FullMatrix<T> * m)663   void full(FullMatrix<T> * m) {
664       full_ = m;
665       rank_ = FULL_BLOCK;
666   }
667 
isNull() const668   bool isNull() const {
669       assert(rank_ >= FULL_BLOCK);
670       return rank_ == 0 || (rank_ == FULL_BLOCK && full_ == NULL);
671   }
672 
673   bool isRecursivelyNull() const;
674   // TODO: the meaning/usage of UNINITIALIZED_BLOCK is not clear, it should be reworked
675   // or removed
isAssembled() const676   bool isAssembled() const {
677       return rank_ > UNINITIALIZED_BLOCK;
678   }
679 
isVoid() const680   bool isVoid() const {
681       return rows()->size() == 0 || cols()->size() == 0;
682   }
683   /**
684    * Tag a not leaf block as assembled.
685    * Must only be called when all leaves of this block have been
686    * assembled (no coherency check).
687    */
assembled()688   void assembled() {
689       assert(!this->isLeaf());
690       rank_ = NONLEAF_BLOCK;
691   }
692 
693   /**
694    * Tag an entire subtree (except the leaves) as assembled.
695    * (with recursion and coherency check: the leaves *must* already be tagged as assembled).
696    */
assembledRecurse()697   void assembledRecurse() {
698     if (!this->isLeaf()) {
699       for (int i=0 ; i<this->nrChild() ; i++)
700         if (this->getChild(i))
701           this->getChild(i)->assembledRecurse();
702       rank_ = NONLEAF_BLOCK;
703     }
704     assert(isAssembled());
705   }
706 
707   /** Set the entire subtree as temporary flag */
708   void temporary(bool b);
709 
rowsTree() const710   const ClusterTree * rowsTree() const {
711       return rows_;
712   }
713 
colsTree() const714   const ClusterTree * colsTree() const {
715       return cols_;
716   }
717 
ownClusterTrees(bool owns_row,bool owns_col)718   void ownClusterTrees(bool owns_row, bool owns_col) {
719       ownRowsClusterTree_ = owns_row;
720       ownColsClusterTree_ = owns_col;
721   }
722 
setColsTree(ClusterTree * clusterTree,bool ownClusterTree)723   void setColsTree(ClusterTree * clusterTree, bool ownClusterTree) {
724       ownColsClusterTree_ = ownClusterTree;
725       cols_ = clusterTree;
726   }
727 
728   /**
729    * Convert this HMatrix to a string for debug.
730    * This is better than overriding << because it allows to use printf.
731    */
732   std::string toString() const;
733 
734   /*! \brief Return a short string describing the content of this HMatrix for debug (like: "HMatrix [320, 452]x[760, 890] norm=13.23442" or "uninitialized" if needed)
735     */
description() const736   std::string description() const {
737     std::ostringstream convert;
738     convert << "HMatrix " << rows()->description() << "x" << cols()->description() ;
739     if (isAssembled())
740       convert << "norm=" << norm();
741     else
742       convert << "uninitialized";
743     return convert.str();
744   }
745 };
746 
747 }  // end namespace hmat
748 
749 #endif
750