xref: /petsc/src/mat/impls/htool/htool.cxx (revision dd39110ba674af879e763a5939b18b2b35b4c87f)
1c7a4214aSPierre Jolivet #include <../src/mat/impls/htool/htool.hpp> /*I "petscmat.h" I*/
2c7a4214aSPierre Jolivet #include <petscblaslapack.h>
3c7a4214aSPierre Jolivet #include <set>
4c7a4214aSPierre Jolivet 
5c7a4214aSPierre Jolivet const char *const MatHtoolCompressorTypes[] = { "sympartialACA", "fullACA", "SVD" };
698e73e17SPierre Jolivet const char *const MatHtoolClusteringTypes[] = { "PCARegular", "PCAGeometric", "BoundingBox1Regular", "BoundingBox1Geometric" };
7c7a4214aSPierre Jolivet const char HtoolCitation[] = "@article{marchand2020two,\n"
8c7a4214aSPierre Jolivet "  Author = {Marchand, Pierre and Claeys, Xavier and Jolivet, Pierre and Nataf, Fr\\'ed\\'eric and Tournier, Pierre-Henri},\n"
9c7a4214aSPierre Jolivet "  Title = {Two-level preconditioning for $h$-version boundary element approximation of hypersingular operator with {GenEO}},\n"
10c7a4214aSPierre Jolivet "  Year = {2020},\n"
11c7a4214aSPierre Jolivet "  Publisher = {Elsevier},\n"
12c7a4214aSPierre Jolivet "  Journal = {Numerische Mathematik},\n"
13c7a4214aSPierre Jolivet "  Volume = {146},\n"
14c7a4214aSPierre Jolivet "  Pages = {597--628},\n"
15c7a4214aSPierre Jolivet "  Url = {https://github.com/htool-ddm/htool}\n"
16c7a4214aSPierre Jolivet "}\n";
17c7a4214aSPierre Jolivet static PetscBool HtoolCite = PETSC_FALSE;
18c7a4214aSPierre Jolivet 
19c7a4214aSPierre Jolivet static PetscErrorCode MatGetDiagonal_Htool(Mat A,Vec v)
20c7a4214aSPierre Jolivet {
21c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
22c7a4214aSPierre Jolivet   PetscScalar    *x;
23c7a4214aSPierre Jolivet   PetscBool      flg;
24c7a4214aSPierre Jolivet 
25c7a4214aSPierre Jolivet   PetscFunctionBegin;
269566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A,&flg));
275f80ce2aSJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only congruent layouts supported");
289566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(v,&x));
29c7a4214aSPierre Jolivet   a->hmatrix->copy_local_diagonal(x);
309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(v,&x));
319566063dSJacob Faibussowitsch   PetscCall(VecScale(v,a->s));
32c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
33c7a4214aSPierre Jolivet }
34c7a4214aSPierre Jolivet 
35c7a4214aSPierre Jolivet static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A,Mat *b)
36c7a4214aSPierre Jolivet {
37c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
38c7a4214aSPierre Jolivet   Mat            B;
39c7a4214aSPierre Jolivet   PetscScalar    *ptr;
40c7a4214aSPierre Jolivet   PetscBool      flg;
41c7a4214aSPierre Jolivet 
42c7a4214aSPierre Jolivet   PetscFunctionBegin;
439566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A,&flg));
445f80ce2aSJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only congruent layouts supported");
459566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B)); /* same logic as in MatGetDiagonalBlock_MPIDense() */
46c7a4214aSPierre Jolivet   if (!B) {
479566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PETSC_COMM_SELF,A->rmap->n,A->rmap->n,A->rmap->n,A->rmap->n,NULL,&B));
489566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite(B,&ptr));
49c7a4214aSPierre Jolivet     a->hmatrix->copy_local_diagonal_block(ptr);
509566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite(B,&ptr));
519566063dSJacob Faibussowitsch     PetscCall(MatPropagateSymmetryOptions(A,B));
529566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
539566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
549566063dSJacob Faibussowitsch     PetscCall(MatScale(B,a->s));
559566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B));
56c7a4214aSPierre Jolivet     *b   = B;
579566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
58c7a4214aSPierre Jolivet   } else *b = B;
59c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
60c7a4214aSPierre Jolivet }
61c7a4214aSPierre Jolivet 
62c7a4214aSPierre Jolivet static PetscErrorCode MatMult_Htool(Mat A,Vec x,Vec y)
63c7a4214aSPierre Jolivet {
64c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
65c7a4214aSPierre Jolivet   const PetscScalar *in;
66c7a4214aSPierre Jolivet   PetscScalar       *out;
67c7a4214aSPierre Jolivet 
68c7a4214aSPierre Jolivet   PetscFunctionBegin;
699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&in));
709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y,&out));
71c7a4214aSPierre Jolivet   a->hmatrix->mvprod_local_to_local(in,out);
729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&in));
739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y,&out));
749566063dSJacob Faibussowitsch   PetscCall(VecScale(y,a->s));
75c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
76c7a4214aSPierre Jolivet }
77c7a4214aSPierre Jolivet 
78c7a4214aSPierre Jolivet /* naive implementation of MatMultAdd() needed for FEM-BEM coupling via MATNEST */
79c7a4214aSPierre Jolivet static PetscErrorCode MatMultAdd_Htool(Mat A,Vec v1,Vec v2,Vec v3)
80c7a4214aSPierre Jolivet {
81c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
82c7a4214aSPierre Jolivet   Vec               tmp;
83c7a4214aSPierre Jolivet   const PetscScalar scale = a->s;
84c7a4214aSPierre Jolivet 
85c7a4214aSPierre Jolivet   PetscFunctionBegin;
869566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(v2,&tmp));
879566063dSJacob Faibussowitsch   PetscCall(VecCopy(v2,v3)); /* no-op in MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]) since VecCopy() checks for x == y */
88c7a4214aSPierre Jolivet   a->s = 1.0; /* set s to 1.0 since VecAXPY() may be used to scale the MatMult() output Vec */
899566063dSJacob Faibussowitsch   PetscCall(MatMult_Htool(A,v1,tmp));
909566063dSJacob Faibussowitsch   PetscCall(VecAXPY(v3,scale,tmp));
919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp));
92c7a4214aSPierre Jolivet   a->s = scale; /* set s back to its original value */
93c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
94c7a4214aSPierre Jolivet }
95c7a4214aSPierre Jolivet 
968b8ddd11SPierre Jolivet static PetscErrorCode MatMultTranspose_Htool(Mat A,Vec x,Vec y)
978b8ddd11SPierre Jolivet {
988b8ddd11SPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
998b8ddd11SPierre Jolivet   const PetscScalar *in;
1008b8ddd11SPierre Jolivet   PetscScalar       *out;
1018b8ddd11SPierre Jolivet 
1028b8ddd11SPierre Jolivet   PetscFunctionBegin;
1039566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&in));
1049566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y,&out));
1058b8ddd11SPierre Jolivet   a->hmatrix->mvprod_transp_local_to_local(in,out);
1069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&in));
1079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y,&out));
1089566063dSJacob Faibussowitsch   PetscCall(VecScale(y,a->s));
1098b8ddd11SPierre Jolivet   PetscFunctionReturn(0);
1108b8ddd11SPierre Jolivet }
1118b8ddd11SPierre Jolivet 
112c7a4214aSPierre Jolivet static PetscErrorCode MatIncreaseOverlap_Htool(Mat A,PetscInt is_max,IS is[],PetscInt ov)
113c7a4214aSPierre Jolivet {
114c7a4214aSPierre Jolivet   std::set<PetscInt> set;
115c7a4214aSPierre Jolivet   const PetscInt     *idx;
1168f308287SPierre Jolivet   PetscInt           *oidx,size,bs[2];
117c7a4214aSPierre Jolivet   PetscMPIInt        csize;
118c7a4214aSPierre Jolivet 
119c7a4214aSPierre Jolivet   PetscFunctionBegin;
1209566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(A,bs,bs+1));
1218f308287SPierre Jolivet   if (bs[0] != bs[1]) bs[0] = 1;
122c7a4214aSPierre Jolivet   for (PetscInt i=0; i<is_max; ++i) {
123c7a4214aSPierre Jolivet     /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */
124c7a4214aSPierre Jolivet     /* needed to avoid subdomain matrices to replicate A since it is dense                           */
1259566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]),&csize));
1265f80ce2aSJacob Faibussowitsch     PetscCheck(csize == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported parallel IS");
1279566063dSJacob Faibussowitsch     PetscCall(ISGetSize(is[i],&size));
1289566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is[i],&idx));
129c7a4214aSPierre Jolivet     for (PetscInt j=0; j<size; ++j) {
130c7a4214aSPierre Jolivet       set.insert(idx[j]);
131c7a4214aSPierre Jolivet       for (PetscInt k=1; k<=ov; ++k) {               /* for each layer of overlap      */
132c7a4214aSPierre Jolivet         if (idx[j] - k >= 0) set.insert(idx[j] - k); /* do not insert negative indices */
133c7a4214aSPierre Jolivet         if (idx[j] + k < A->rmap->N && idx[j] + k < A->cmap->N) set.insert(idx[j] + k); /* do not insert indices greater than the dimension of A */
134c7a4214aSPierre Jolivet       }
135c7a4214aSPierre Jolivet     }
1369566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is[i],&idx));
1379566063dSJacob Faibussowitsch     PetscCall(ISDestroy(is+i));
1388f308287SPierre Jolivet     if (bs[0] > 1) {
1398f308287SPierre Jolivet       for (std::set<PetscInt>::iterator it=set.cbegin(); it!=set.cend(); it++) {
1408f308287SPierre Jolivet         std::vector<PetscInt> block(bs[0]);
1418f308287SPierre Jolivet         std::iota(block.begin(),block.end(),(*it/bs[0])*bs[0]);
1428f308287SPierre Jolivet         set.insert(block.cbegin(),block.cend());
1438f308287SPierre Jolivet       }
1448f308287SPierre Jolivet     }
145c7a4214aSPierre Jolivet     size = set.size(); /* size with overlap */
1469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size,&oidx));
147c7a4214aSPierre Jolivet     for (const PetscInt j : set) *oidx++ = j;
148c7a4214aSPierre Jolivet     oidx -= size;
1499566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,size,oidx,PETSC_OWN_POINTER,is+i));
150c7a4214aSPierre Jolivet   }
151c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
152c7a4214aSPierre Jolivet }
153c7a4214aSPierre Jolivet 
154c7a4214aSPierre Jolivet static PetscErrorCode MatCreateSubMatrices_Htool(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
155c7a4214aSPierre Jolivet {
156c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
157c7a4214aSPierre Jolivet   Mat               D,B,BT;
158c7a4214aSPierre Jolivet   const PetscScalar *copy;
159c7a4214aSPierre Jolivet   PetscScalar       *ptr;
160c7a4214aSPierre Jolivet   const PetscInt    *idxr,*idxc,*it;
161c7a4214aSPierre Jolivet   PetscInt          nrow,m,i;
162c7a4214aSPierre Jolivet   PetscBool         flg;
163c7a4214aSPierre Jolivet 
164c7a4214aSPierre Jolivet   PetscFunctionBegin;
165c7a4214aSPierre Jolivet   if (scall != MAT_REUSE_MATRIX) {
1669566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(n,submat));
167c7a4214aSPierre Jolivet   }
168c7a4214aSPierre Jolivet   for (i=0; i<n; ++i) {
1699566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(irow[i],&nrow));
1709566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(icol[i],&m));
1719566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(irow[i],&idxr));
1729566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(icol[i],&idxc));
173c7a4214aSPierre Jolivet     if (scall != MAT_REUSE_MATRIX) {
1749566063dSJacob Faibussowitsch       PetscCall(MatCreateDense(PETSC_COMM_SELF,nrow,m,nrow,m,NULL,(*submat)+i));
175c7a4214aSPierre Jolivet     }
1769566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite((*submat)[i],&ptr));
177c7a4214aSPierre Jolivet     if (irow[i] == icol[i]) { /* same row and column IS? */
1789566063dSJacob Faibussowitsch       PetscCall(MatHasCongruentLayouts(A,&flg));
179c7a4214aSPierre Jolivet       if (flg) {
1809566063dSJacob Faibussowitsch         PetscCall(ISSorted(irow[i],&flg));
181c7a4214aSPierre Jolivet         if (flg) { /* sorted IS? */
182c7a4214aSPierre Jolivet           it = std::lower_bound(idxr,idxr+nrow,A->rmap->rstart);
183c7a4214aSPierre Jolivet           if (it != idxr+nrow && *it == A->rmap->rstart) { /* rmap->rstart in IS? */
184c7a4214aSPierre Jolivet             if (std::distance(idxr,it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */
185c7a4214aSPierre Jolivet               for (PetscInt j=0; j<A->rmap->n && flg; ++j) if (PetscUnlikely(it[j] != A->rmap->rstart+j)) flg = PETSC_FALSE;
186c7a4214aSPierre Jolivet               if (flg) { /* complete local diagonal block in IS? */
187c7a4214aSPierre Jolivet                 /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM
188c7a4214aSPierre Jolivet                  *      [   B   C   E   ]
189c7a4214aSPierre Jolivet                  *  A = [   B   D   E   ]
190c7a4214aSPierre Jolivet                  *      [   B   F   E   ]
191c7a4214aSPierre Jolivet                  */
192c7a4214aSPierre Jolivet                 m = std::distance(idxr,it); /* shift of the coefficient (0,0) of block D from above */
1939566063dSJacob Faibussowitsch                 PetscCall(MatGetDiagonalBlock_Htool(A,&D));
1949566063dSJacob Faibussowitsch                 PetscCall(MatDenseGetArrayRead(D,&copy));
195c7a4214aSPierre Jolivet                 for (PetscInt k=0; k<A->rmap->n; ++k) {
1969566063dSJacob Faibussowitsch                   PetscCall(PetscArraycpy(ptr+(m+k)*nrow+m,copy+k*A->rmap->n,A->rmap->n)); /* block D from above */
197c7a4214aSPierre Jolivet                 }
1989566063dSJacob Faibussowitsch                 PetscCall(MatDenseRestoreArrayRead(D,&copy));
199c7a4214aSPierre Jolivet                 if (m) {
200c7a4214aSPierre Jolivet                   a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* vertical block B from above */
201c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
202c7a4214aSPierre Jolivet                   if (A->symmetric || A->hermitian) {
2039566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF,A->rmap->n,m,A->rmap->n,m,ptr+m,&B));
2049566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B,nrow));
2059566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF,m,A->rmap->n,m,A->rmap->n,ptr+m*nrow,&BT));
2069566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT,nrow));
207c7a4214aSPierre Jolivet                     if (A->hermitian && PetscDefined(USE_COMPLEX)) {
2089566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT));
209c7a4214aSPierre Jolivet                     } else {
2109566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B,MAT_REUSE_MATRIX,&BT));
211c7a4214aSPierre Jolivet                     }
2129566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2139566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
214c7a4214aSPierre Jolivet                   } else {
215c7a4214aSPierre Jolivet                     for (PetscInt k=0; k<A->rmap->n; ++k) { /* block C from above */
216c7a4214aSPierre Jolivet                       a->wrapper->copy_submatrix(m,1,idxr,idxc+m+k,ptr+(m+k)*nrow);
217c7a4214aSPierre Jolivet                     }
218c7a4214aSPierre Jolivet                   }
219c7a4214aSPierre Jolivet                 }
220c7a4214aSPierre Jolivet                 if (m+A->rmap->n != nrow) {
221c7a4214aSPierre Jolivet                   a->wrapper->copy_submatrix(nrow,std::distance(it+A->rmap->n,idxr+nrow),idxr,idxc+m+A->rmap->n,ptr+(m+A->rmap->n)*nrow); /* vertical block E from above */
222c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
223c7a4214aSPierre Jolivet                   if (A->symmetric || A->hermitian) {
2249566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF,A->rmap->n,nrow-(m+A->rmap->n),A->rmap->n,nrow-(m+A->rmap->n),ptr+(m+A->rmap->n)*nrow+m,&B));
2259566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B,nrow));
2269566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF,nrow-(m+A->rmap->n),A->rmap->n,nrow-(m+A->rmap->n),A->rmap->n,ptr+m*nrow+m+A->rmap->n,&BT));
2279566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT,nrow));
228c7a4214aSPierre Jolivet                     if (A->hermitian && PetscDefined(USE_COMPLEX)) {
2299566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT));
230c7a4214aSPierre Jolivet                     } else {
2319566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B,MAT_REUSE_MATRIX,&BT));
232c7a4214aSPierre Jolivet                     }
2339566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2349566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
235c7a4214aSPierre Jolivet                   } else {
236c7a4214aSPierre Jolivet                     for (PetscInt k=0; k<A->rmap->n; ++k) { /* block F from above */
237c7a4214aSPierre Jolivet                       a->wrapper->copy_submatrix(std::distance(it+A->rmap->n,idxr+nrow),1,it+A->rmap->n,idxc+m+k,ptr+(m+k)*nrow+m+A->rmap->n);
238c7a4214aSPierre Jolivet                     }
239c7a4214aSPierre Jolivet                   }
240c7a4214aSPierre Jolivet                 }
241c7a4214aSPierre Jolivet               } /* complete local diagonal block not in IS */
242c7a4214aSPierre Jolivet             } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */
243c7a4214aSPierre Jolivet           } else flg = PETSC_FALSE; /* rmap->rstart not in IS */
244c7a4214aSPierre Jolivet         } /* unsorted IS */
245c7a4214aSPierre Jolivet       }
246c7a4214aSPierre Jolivet     } else flg = PETSC_FALSE; /* different row and column IS */
247c7a4214aSPierre Jolivet     if (!flg) a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* reassemble everything */
2489566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(irow[i],&idxr));
2499566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(icol[i],&idxc));
2509566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite((*submat)[i],&ptr));
2519566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin((*submat)[i],MAT_FINAL_ASSEMBLY));
2529566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd((*submat)[i],MAT_FINAL_ASSEMBLY));
2539566063dSJacob Faibussowitsch     PetscCall(MatScale((*submat)[i],a->s));
254c7a4214aSPierre Jolivet   }
255c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
256c7a4214aSPierre Jolivet }
257c7a4214aSPierre Jolivet 
258c7a4214aSPierre Jolivet static PetscErrorCode MatDestroy_Htool(Mat A)
259c7a4214aSPierre Jolivet {
260c7a4214aSPierre Jolivet   Mat_Htool               *a = (Mat_Htool*)A->data;
261c7a4214aSPierre Jolivet   PetscContainer          container;
262c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
263c7a4214aSPierre Jolivet 
264c7a4214aSPierre Jolivet   PetscFunctionBegin;
2659566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A,NULL));
2669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",NULL));
2679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",NULL));
2689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",NULL));
2699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",NULL));
2709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",NULL));
2719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",NULL));
2729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",NULL));
2739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",NULL));
2749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",NULL));
2759566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A,"KernelTranspose",(PetscObject*)&container));
276c7a4214aSPierre Jolivet   if (container) { /* created in MatTranspose_Htool() */
2779566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container,(void**)&kernelt));
2789566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
2799566063dSJacob Faibussowitsch     PetscCall(PetscFree(kernelt));
2809566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&container));
2819566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)A,"KernelTranspose",NULL));
282c7a4214aSPierre Jolivet   }
283c7a4214aSPierre Jolivet   if (a->gcoords_source != a->gcoords_target) {
2849566063dSJacob Faibussowitsch     PetscCall(PetscFree(a->gcoords_source));
285c7a4214aSPierre Jolivet   }
2869566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->gcoords_target));
2879566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a->work_source,a->work_target));
288c7a4214aSPierre Jolivet   delete a->wrapper;
289c7a4214aSPierre Jolivet   delete a->hmatrix;
2909566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
291c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
292c7a4214aSPierre Jolivet }
293c7a4214aSPierre Jolivet 
294c7a4214aSPierre Jolivet static PetscErrorCode MatView_Htool(Mat A,PetscViewer pv)
295c7a4214aSPierre Jolivet {
296c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
297c7a4214aSPierre Jolivet   PetscBool      flg;
298c7a4214aSPierre Jolivet 
299c7a4214aSPierre Jolivet   PetscFunctionBegin;
3009566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pv,PETSCVIEWERASCII,&flg));
301c7a4214aSPierre Jolivet   if (flg) {
3029566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"symmetry: %c\n",a->hmatrix->get_symmetry_type()));
303c7a4214aSPierre Jolivet     if (PetscAbsScalar(a->s-1.0) > PETSC_MACHINE_EPSILON) {
304c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX)
3059566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(pv,"scaling: %g+%gi\n",(double)PetscRealPart(a->s),(double)PetscImaginaryPart(a->s)));
306c7a4214aSPierre Jolivet #else
3079566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(pv,"scaling: %g\n",(double)a->s));
308c7a4214aSPierre Jolivet #endif
309c7a4214aSPierre Jolivet     }
3109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"minimum cluster size: %" PetscInt_FMT "\n",a->bs[0]));
3119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"maximum block size: %" PetscInt_FMT "\n",a->bs[1]));
3129566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"epsilon: %g\n",(double)a->epsilon));
3139566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"eta: %g\n",(double)a->eta));
3149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"minimum target depth: %" PetscInt_FMT "\n",a->depth[0]));
3159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"minimum source depth: %" PetscInt_FMT "\n",a->depth[1]));
3169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"compressor: %s\n",MatHtoolCompressorTypes[a->compressor]));
3179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"clustering: %s\n",MatHtoolClusteringTypes[a->clustering]));
3189566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"compression ratio: %s\n",a->hmatrix->get_infos("Compression_ratio").c_str()));
3199566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"space saving: %s\n",a->hmatrix->get_infos("Space_saving").c_str()));
3209566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"number of dense (resp. low rank) matrices: %s (resp. %s)\n",a->hmatrix->get_infos("Number_of_dmat").c_str(),a->hmatrix->get_infos("Number_of_lrmat").c_str()));
3219566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"(minimum, mean, maximum) dense block sizes: (%s, %s, %s)\n",a->hmatrix->get_infos("Dense_block_size_min").c_str(),a->hmatrix->get_infos("Dense_block_size_mean").c_str(),a->hmatrix->get_infos("Dense_block_size_max").c_str()));
3229566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"(minimum, mean, maximum) low rank block sizes: (%s, %s, %s)\n",a->hmatrix->get_infos("Low_rank_block_size_min").c_str(),a->hmatrix->get_infos("Low_rank_block_size_mean").c_str(),a->hmatrix->get_infos("Low_rank_block_size_max").c_str()));
3239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"(minimum, mean, maximum) ranks: (%s, %s, %s)\n",a->hmatrix->get_infos("Rank_min").c_str(),a->hmatrix->get_infos("Rank_mean").c_str(),a->hmatrix->get_infos("Rank_max").c_str()));
324c7a4214aSPierre Jolivet   }
325c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
326c7a4214aSPierre Jolivet }
327c7a4214aSPierre Jolivet 
328c7a4214aSPierre Jolivet static PetscErrorCode MatScale_Htool(Mat A,PetscScalar s)
329c7a4214aSPierre Jolivet {
330c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
331c7a4214aSPierre Jolivet 
332c7a4214aSPierre Jolivet   PetscFunctionBegin;
333c7a4214aSPierre Jolivet   a->s *= s;
334c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
335c7a4214aSPierre Jolivet }
336c7a4214aSPierre Jolivet 
337c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */
338c7a4214aSPierre Jolivet static PetscErrorCode MatGetRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
339c7a4214aSPierre Jolivet {
340c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
341c7a4214aSPierre Jolivet   PetscInt       *idxc;
342c7a4214aSPierre Jolivet   PetscBLASInt   one = 1,bn;
343c7a4214aSPierre Jolivet 
344c7a4214aSPierre Jolivet   PetscFunctionBegin;
345c7a4214aSPierre Jolivet   if (nz) *nz = A->cmap->N;
346c7a4214aSPierre Jolivet   if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */
3479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N,&idxc));
348c7a4214aSPierre Jolivet     for (PetscInt i=0; i<A->cmap->N; ++i) idxc[i] = i;
349c7a4214aSPierre Jolivet   }
350c7a4214aSPierre Jolivet   if (idx) *idx = idxc;
351c7a4214aSPierre Jolivet   if (v) {
3529566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N,v));
353c7a4214aSPierre Jolivet     if (a->wrapper) a->wrapper->copy_submatrix(1,A->cmap->N,&row,idxc,*v);
354cab00e0dSPierre Jolivet     else reinterpret_cast<htool::VirtualGenerator<PetscScalar>*>(a->kernelctx)->copy_submatrix(1,A->cmap->N,&row,idxc,*v);
3559566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(A->cmap->N,&bn));
356c7a4214aSPierre Jolivet     PetscStackCallBLAS("BLASscal",BLASscal_(&bn,&a->s,*v,&one));
357c7a4214aSPierre Jolivet   }
358c7a4214aSPierre Jolivet   if (!idx) {
3599566063dSJacob Faibussowitsch     PetscCall(PetscFree(idxc));
360c7a4214aSPierre Jolivet   }
361c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
362c7a4214aSPierre Jolivet }
363c7a4214aSPierre Jolivet 
364c7a4214aSPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
365c7a4214aSPierre Jolivet {
366c7a4214aSPierre Jolivet   PetscFunctionBegin;
367c7a4214aSPierre Jolivet   if (nz) *nz = 0;
368c7a4214aSPierre Jolivet   if (idx) {
3699566063dSJacob Faibussowitsch     PetscCall(PetscFree(*idx));
370c7a4214aSPierre Jolivet   }
371c7a4214aSPierre Jolivet   if (v) {
3729566063dSJacob Faibussowitsch     PetscCall(PetscFree(*v));
373c7a4214aSPierre Jolivet   }
374c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
375c7a4214aSPierre Jolivet }
376c7a4214aSPierre Jolivet 
377c7a4214aSPierre Jolivet static PetscErrorCode MatSetFromOptions_Htool(PetscOptionItems *PetscOptionsObject,Mat A)
378c7a4214aSPierre Jolivet {
379c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
380c7a4214aSPierre Jolivet   PetscInt       n;
381c7a4214aSPierre Jolivet   PetscBool      flg;
382c7a4214aSPierre Jolivet 
383c7a4214aSPierre Jolivet   PetscFunctionBegin;
3849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"Htool options"));
3859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_min_cluster_size","Minimal leaf size in cluster tree",NULL,a->bs[0],a->bs,NULL));
3869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_max_block_size","Maximal number of coefficients in a dense block",NULL,a->bs[1],a->bs + 1,NULL));
3879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_htool_epsilon","Relative error in Frobenius norm when approximating a block",NULL,a->epsilon,&a->epsilon,NULL));
3889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_htool_eta","Admissibility condition tolerance",NULL,a->eta,&a->eta,NULL));
3899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_min_target_depth","Minimal cluster tree depth associated with the rows",NULL,a->depth[0],a->depth,NULL));
3909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_min_source_depth","Minimal cluster tree depth associated with the columns",NULL,a->depth[1],a->depth + 1,NULL));
391c7a4214aSPierre Jolivet   n = 0;
392*dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_htool_compressor","Type of compression","MatHtoolCompressorType",MatHtoolCompressorTypes,PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes),MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA],&n,&flg));
393c7a4214aSPierre Jolivet   if (flg) a->compressor = MatHtoolCompressorType(n);
39498e73e17SPierre Jolivet   n = 0;
395*dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_htool_clustering","Type of clustering","MatHtoolClusteringType",MatHtoolClusteringTypes,PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes),MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR],&n,&flg));
39698e73e17SPierre Jolivet   if (flg) a->clustering = MatHtoolClusteringType(n);
3979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
398c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
399c7a4214aSPierre Jolivet }
400c7a4214aSPierre Jolivet 
401c7a4214aSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A,MatAssemblyType type)
402c7a4214aSPierre Jolivet {
403c7a4214aSPierre Jolivet   Mat_Htool                                                    *a = (Mat_Htool*)A->data;
404c7a4214aSPierre Jolivet   const PetscInt                                               *ranges;
405c7a4214aSPierre Jolivet   PetscInt                                                     *offset;
406c7a4214aSPierre Jolivet   PetscMPIInt                                                  size;
407c7a4214aSPierre Jolivet   char                                                         S = PetscDefined(USE_COMPLEX) && A->hermitian ? 'H' : (A->symmetric ? 'S' : 'N'),uplo = S == 'N' ? 'N' : 'U';
408cab00e0dSPierre Jolivet   htool::VirtualGenerator<PetscScalar>                         *generator = nullptr;
40998e73e17SPierre Jolivet   std::shared_ptr<htool::VirtualCluster>                       t,s = nullptr;
4103b9338faSPierre Jolivet   std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor = nullptr;
411c7a4214aSPierre Jolivet 
412c7a4214aSPierre Jolivet   PetscFunctionBegin;
4139566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(HtoolCitation,&HtoolCite));
414c7a4214aSPierre Jolivet   delete a->wrapper;
415c7a4214aSPierre Jolivet   delete a->hmatrix;
4169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
4179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*size,&offset));
4189566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(A,&ranges));
419c7a4214aSPierre Jolivet   for (PetscInt i=0; i<size; ++i) {
420c7a4214aSPierre Jolivet     offset[2*i] = ranges[i];
421c7a4214aSPierre Jolivet     offset[2*i+1] = ranges[i+1] - ranges[i];
422c7a4214aSPierre Jolivet   }
42398e73e17SPierre Jolivet   switch (a->clustering) {
42498e73e17SPierre Jolivet   case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
42598e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
42698e73e17SPierre Jolivet     break;
42798e73e17SPierre Jolivet   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
42898e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
42998e73e17SPierre Jolivet     break;
43098e73e17SPierre Jolivet   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
43198e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim);
43298e73e17SPierre Jolivet     break;
43398e73e17SPierre Jolivet   default:
43498e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim);
43598e73e17SPierre Jolivet   }
436c7a4214aSPierre Jolivet   t->set_minclustersize(a->bs[0]);
43798e73e17SPierre Jolivet   t->build(A->rmap->N,a->gcoords_target,offset);
438c7a4214aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx);
439c7a4214aSPierre Jolivet   else {
440c7a4214aSPierre Jolivet     a->wrapper = NULL;
441cab00e0dSPierre Jolivet     generator = reinterpret_cast<htool::VirtualGenerator<PetscScalar>*>(a->kernelctx);
442c7a4214aSPierre Jolivet   }
443c7a4214aSPierre Jolivet   if (a->gcoords_target != a->gcoords_source) {
4449566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangesColumn(A,&ranges));
445c7a4214aSPierre Jolivet     for (PetscInt i=0; i<size; ++i) {
446c7a4214aSPierre Jolivet       offset[2*i] = ranges[i];
447c7a4214aSPierre Jolivet       offset[2*i+1] = ranges[i+1] - ranges[i];
448c7a4214aSPierre Jolivet     }
44998e73e17SPierre Jolivet     switch (a->clustering) {
45098e73e17SPierre Jolivet     case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
45198e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
45298e73e17SPierre Jolivet       break;
45398e73e17SPierre Jolivet     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
45498e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
45598e73e17SPierre Jolivet       break;
45698e73e17SPierre Jolivet     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
45798e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim);
45898e73e17SPierre Jolivet       break;
45998e73e17SPierre Jolivet     default:
46098e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim);
46198e73e17SPierre Jolivet     }
462c7a4214aSPierre Jolivet     s->set_minclustersize(a->bs[0]);
46398e73e17SPierre Jolivet     s->build(A->cmap->N,a->gcoords_source,offset);
464c7a4214aSPierre Jolivet     S = uplo = 'N';
465c7a4214aSPierre Jolivet   }
4669566063dSJacob Faibussowitsch   PetscCall(PetscFree(offset));
467c7a4214aSPierre Jolivet   switch (a->compressor) {
468c7a4214aSPierre Jolivet   case MAT_HTOOL_COMPRESSOR_FULL_ACA:
4693b9338faSPierre Jolivet     compressor = std::make_shared<htool::fullACA<PetscScalar>>();
470c7a4214aSPierre Jolivet     break;
471c7a4214aSPierre Jolivet   case MAT_HTOOL_COMPRESSOR_SVD:
4723b9338faSPierre Jolivet     compressor = std::make_shared<htool::SVD<PetscScalar>>();
473c7a4214aSPierre Jolivet     break;
474c7a4214aSPierre Jolivet   default:
4753b9338faSPierre Jolivet     compressor = std::make_shared<htool::sympartialACA<PetscScalar>>();
476c7a4214aSPierre Jolivet   }
4773b9338faSPierre Jolivet   a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar>*>(new htool::HMatrix<PetscScalar>(t,s ? s : t,a->epsilon,a->eta,S,uplo));
4783b9338faSPierre Jolivet   a->hmatrix->set_compression(compressor);
479c7a4214aSPierre Jolivet   a->hmatrix->set_maxblocksize(a->bs[1]);
480c7a4214aSPierre Jolivet   a->hmatrix->set_mintargetdepth(a->depth[0]);
481c7a4214aSPierre Jolivet   a->hmatrix->set_minsourcedepth(a->depth[1]);
4823b9338faSPierre Jolivet   if (s) a->hmatrix->build(a->wrapper ? *a->wrapper : *generator,a->gcoords_target,a->gcoords_source);
4833b9338faSPierre Jolivet   else   a->hmatrix->build(a->wrapper ? *a->wrapper : *generator,a->gcoords_target);
484c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
485c7a4214aSPierre Jolivet }
486c7a4214aSPierre Jolivet 
487c7a4214aSPierre Jolivet static PetscErrorCode MatProductNumeric_Htool(Mat C)
488c7a4214aSPierre Jolivet {
489c7a4214aSPierre Jolivet   Mat_Product       *product = C->product;
490c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)product->A->data;
491c7a4214aSPierre Jolivet   const PetscScalar *in;
492c7a4214aSPierre Jolivet   PetscScalar       *out;
4938b8ddd11SPierre Jolivet   PetscInt          N,lda;
494c7a4214aSPierre Jolivet 
495c7a4214aSPierre Jolivet   PetscFunctionBegin;
496c7a4214aSPierre Jolivet   MatCheckProduct(C,1);
4979566063dSJacob Faibussowitsch   PetscCall(MatGetSize(C,NULL,&N));
4989566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(C,&lda));
4995f80ce2aSJacob Faibussowitsch   PetscCheck(lda == C->rmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")",lda,C->rmap->n);
5009566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(product->B,&in));
5019566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C,&out));
5028b8ddd11SPierre Jolivet   switch (product->type) {
5038b8ddd11SPierre Jolivet   case MATPRODUCT_AB:
504c7a4214aSPierre Jolivet     a->hmatrix->mvprod_local_to_local(in,out,N);
5058b8ddd11SPierre Jolivet     break;
5068b8ddd11SPierre Jolivet   case MATPRODUCT_AtB:
5078b8ddd11SPierre Jolivet     a->hmatrix->mvprod_transp_local_to_local(in,out,N);
508c7a4214aSPierre Jolivet     break;
509c7a4214aSPierre Jolivet   default:
51098921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType %s is not supported",MatProductTypes[product->type]);
511c7a4214aSPierre Jolivet   }
5129566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C,&out));
5139566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(product->B,&in));
5149566063dSJacob Faibussowitsch   PetscCall(MatScale(C,a->s));
515c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
516c7a4214aSPierre Jolivet }
517c7a4214aSPierre Jolivet 
518c7a4214aSPierre Jolivet static PetscErrorCode MatProductSymbolic_Htool(Mat C)
519c7a4214aSPierre Jolivet {
520c7a4214aSPierre Jolivet   Mat_Product    *product = C->product;
521c7a4214aSPierre Jolivet   Mat            A,B;
522c7a4214aSPierre Jolivet   PetscBool      flg;
523c7a4214aSPierre Jolivet 
524c7a4214aSPierre Jolivet   PetscFunctionBegin;
525c7a4214aSPierre Jolivet   MatCheckProduct(C,1);
526c7a4214aSPierre Jolivet   A = product->A;
527c7a4214aSPierre Jolivet   B = product->B;
5289566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,""));
5295f80ce2aSJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MatProduct_AB not supported for %s",((PetscObject)product->B)->type_name);
530c7a4214aSPierre Jolivet   switch (product->type) {
531c7a4214aSPierre Jolivet   case MATPRODUCT_AB:
532c7a4214aSPierre Jolivet     if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) {
5339566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(C,A->rmap->n,B->cmap->n,A->rmap->N,B->cmap->N));
534c7a4214aSPierre Jolivet     }
5358b8ddd11SPierre Jolivet     break;
5368b8ddd11SPierre Jolivet   case MATPRODUCT_AtB:
5378b8ddd11SPierre Jolivet     if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) {
5389566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(C,A->cmap->n,B->cmap->n,A->cmap->N,B->cmap->N));
5398b8ddd11SPierre Jolivet     }
5408b8ddd11SPierre Jolivet     break;
5418b8ddd11SPierre Jolivet   default:
54298921bdaSJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"ProductType %s is not supported",MatProductTypes[product->type]);
5438b8ddd11SPierre Jolivet   }
5449566063dSJacob Faibussowitsch   PetscCall(MatSetType(C,MATDENSE));
5459566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
5469566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE));
5479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
5489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
549c7a4214aSPierre Jolivet   C->ops->productsymbolic = NULL;
550c7a4214aSPierre Jolivet   C->ops->productnumeric = MatProductNumeric_Htool;
551c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
552c7a4214aSPierre Jolivet }
553c7a4214aSPierre Jolivet 
554c7a4214aSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Htool(Mat C)
555c7a4214aSPierre Jolivet {
556c7a4214aSPierre Jolivet   PetscFunctionBegin;
557c7a4214aSPierre Jolivet   MatCheckProduct(C,1);
5588b8ddd11SPierre Jolivet   if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool;
559c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
560c7a4214aSPierre Jolivet }
561c7a4214aSPierre Jolivet 
56298e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix)
563c7a4214aSPierre Jolivet {
564c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
565c7a4214aSPierre Jolivet 
566c7a4214aSPierre Jolivet   PetscFunctionBegin;
567c7a4214aSPierre Jolivet   *hmatrix = a->hmatrix;
568c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
569c7a4214aSPierre Jolivet }
570c7a4214aSPierre Jolivet 
571c7a4214aSPierre Jolivet /*@C
572c7a4214aSPierre Jolivet      MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a MATHTOOL.
573c7a4214aSPierre Jolivet 
574c7a4214aSPierre Jolivet    Input Parameter:
575c7a4214aSPierre Jolivet .     A - hierarchical matrix
576c7a4214aSPierre Jolivet 
577c7a4214aSPierre Jolivet    Output Parameter:
578c7a4214aSPierre Jolivet .     hmatrix - opaque pointer to a Htool virtual matrix
579c7a4214aSPierre Jolivet 
580c7a4214aSPierre Jolivet    Level: advanced
581c7a4214aSPierre Jolivet 
582c7a4214aSPierre Jolivet .seealso:  MATHTOOL
583c7a4214aSPierre Jolivet @*/
58498e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix)
585c7a4214aSPierre Jolivet {
586c7a4214aSPierre Jolivet   PetscFunctionBegin;
587c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
588c7a4214aSPierre Jolivet   PetscValidPointer(hmatrix,2);
589cac4c232SBarry Smith   PetscTryMethod(A,"MatHtoolGetHierarchicalMat_C",(Mat,const htool::VirtualHMatrix<PetscScalar>**),(A,hmatrix));
590c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
591c7a4214aSPierre Jolivet }
592c7a4214aSPierre Jolivet 
593c7a4214aSPierre Jolivet static PetscErrorCode MatHtoolSetKernel_Htool(Mat A,MatHtoolKernel kernel,void *kernelctx)
594c7a4214aSPierre Jolivet {
595c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
596c7a4214aSPierre Jolivet 
597c7a4214aSPierre Jolivet   PetscFunctionBegin;
598c7a4214aSPierre Jolivet   a->kernel    = kernel;
599c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
600c7a4214aSPierre Jolivet   delete a->wrapper;
601c7a4214aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx);
602c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
603c7a4214aSPierre Jolivet }
604c7a4214aSPierre Jolivet 
605c7a4214aSPierre Jolivet /*@C
60698e73e17SPierre Jolivet      MatHtoolSetKernel - Sets the kernel and context used for the assembly of a MATHTOOL.
607c7a4214aSPierre Jolivet 
608c7a4214aSPierre Jolivet    Input Parameters:
609c7a4214aSPierre Jolivet +     A - hierarchical matrix
610c7a4214aSPierre Jolivet .     kernel - computational kernel (or NULL)
611cab00e0dSPierre Jolivet -     kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
612c7a4214aSPierre Jolivet 
613c7a4214aSPierre Jolivet    Level: advanced
614c7a4214aSPierre Jolivet 
615c7a4214aSPierre Jolivet .seealso:  MATHTOOL, MatCreateHtoolFromKernel()
616c7a4214aSPierre Jolivet @*/
617c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat A,MatHtoolKernel kernel,void *kernelctx)
618c7a4214aSPierre Jolivet {
619c7a4214aSPierre Jolivet   PetscFunctionBegin;
620c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
621c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel,2);
622c7a4214aSPierre Jolivet   if (!kernel)    PetscValidPointer(kernelctx,3);
623cac4c232SBarry Smith   PetscTryMethod(A,"MatHtoolSetKernel_C",(Mat,MatHtoolKernel,void*),(A,kernel,kernelctx));
624c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
625c7a4214aSPierre Jolivet }
626c7a4214aSPierre Jolivet 
62798e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A,IS* is)
62898e73e17SPierre Jolivet {
62998e73e17SPierre Jolivet   Mat_Htool             *a = (Mat_Htool*)A->data;
63098e73e17SPierre Jolivet   std::vector<PetscInt> source;
63198e73e17SPierre Jolivet 
63298e73e17SPierre Jolivet   PetscFunctionBegin;
633a9918087SPierre Jolivet   source = a->hmatrix->get_source_cluster()->get_local_perm();
6349566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A),source.size(),source.data(),PETSC_COPY_VALUES,is));
6359566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
63698e73e17SPierre Jolivet   PetscFunctionReturn(0);
63798e73e17SPierre Jolivet }
63898e73e17SPierre Jolivet 
63998e73e17SPierre Jolivet /*@C
64098e73e17SPierre Jolivet      MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster.
64198e73e17SPierre Jolivet 
64298e73e17SPierre Jolivet    Input Parameter:
64398e73e17SPierre Jolivet .     A - hierarchical matrix
64498e73e17SPierre Jolivet 
64598e73e17SPierre Jolivet    Output Parameter:
64698e73e17SPierre Jolivet .     is - permutation
64798e73e17SPierre Jolivet 
64898e73e17SPierre Jolivet    Level: advanced
64998e73e17SPierre Jolivet 
65098e73e17SPierre Jolivet .seealso:  MATHTOOL, MatHtoolGetPermutationTarget(), MatHtoolUsePermutation()
65198e73e17SPierre Jolivet @*/
65298e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat A,IS* is)
65398e73e17SPierre Jolivet {
65498e73e17SPierre Jolivet   PetscFunctionBegin;
65598e73e17SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
65698e73e17SPierre Jolivet   if (!is) PetscValidPointer(is,2);
657cac4c232SBarry Smith   PetscTryMethod(A,"MatHtoolGetPermutationSource_C",(Mat,IS*),(A,is));
65898e73e17SPierre Jolivet   PetscFunctionReturn(0);
65998e73e17SPierre Jolivet }
66098e73e17SPierre Jolivet 
66198e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A,IS* is)
66298e73e17SPierre Jolivet {
66398e73e17SPierre Jolivet   Mat_Htool             *a = (Mat_Htool*)A->data;
66498e73e17SPierre Jolivet   std::vector<PetscInt> target;
66598e73e17SPierre Jolivet 
66698e73e17SPierre Jolivet   PetscFunctionBegin;
667a9918087SPierre Jolivet   target = a->hmatrix->get_target_cluster()->get_local_perm();
6689566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A),target.size(),target.data(),PETSC_COPY_VALUES,is));
6699566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
67098e73e17SPierre Jolivet   PetscFunctionReturn(0);
67198e73e17SPierre Jolivet }
67298e73e17SPierre Jolivet 
67398e73e17SPierre Jolivet /*@C
67498e73e17SPierre Jolivet      MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster.
67598e73e17SPierre Jolivet 
67698e73e17SPierre Jolivet    Input Parameter:
67798e73e17SPierre Jolivet .     A - hierarchical matrix
67898e73e17SPierre Jolivet 
67998e73e17SPierre Jolivet    Output Parameter:
68098e73e17SPierre Jolivet .     is - permutation
68198e73e17SPierre Jolivet 
68298e73e17SPierre Jolivet    Level: advanced
68398e73e17SPierre Jolivet 
68498e73e17SPierre Jolivet .seealso:  MATHTOOL, MatHtoolGetPermutationSource(), MatHtoolUsePermutation()
68598e73e17SPierre Jolivet @*/
68698e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat A,IS* is)
68798e73e17SPierre Jolivet {
68898e73e17SPierre Jolivet   PetscFunctionBegin;
68998e73e17SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
69098e73e17SPierre Jolivet   if (!is) PetscValidPointer(is,2);
691cac4c232SBarry Smith   PetscTryMethod(A,"MatHtoolGetPermutationTarget_C",(Mat,IS*),(A,is));
69298e73e17SPierre Jolivet   PetscFunctionReturn(0);
69398e73e17SPierre Jolivet }
69498e73e17SPierre Jolivet 
69598e73e17SPierre Jolivet static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A,PetscBool use)
69698e73e17SPierre Jolivet {
69798e73e17SPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
69898e73e17SPierre Jolivet 
69998e73e17SPierre Jolivet   PetscFunctionBegin;
70098e73e17SPierre Jolivet   a->hmatrix->set_use_permutation(use);
70198e73e17SPierre Jolivet   PetscFunctionReturn(0);
70298e73e17SPierre Jolivet }
70398e73e17SPierre Jolivet 
70498e73e17SPierre Jolivet /*@C
70598e73e17SPierre Jolivet      MatHtoolUsePermutation - Sets whether MATHTOOL should permute input (resp. output) vectors following its internal source (resp. target) permutation.
70698e73e17SPierre Jolivet 
70798e73e17SPierre Jolivet    Input Parameters:
70898e73e17SPierre Jolivet +     A - hierarchical matrix
70998e73e17SPierre Jolivet -     use - Boolean value
71098e73e17SPierre Jolivet 
71198e73e17SPierre Jolivet    Level: advanced
71298e73e17SPierre Jolivet 
71398e73e17SPierre Jolivet .seealso:  MATHTOOL, MatHtoolGetPermutationSource(), MatHtoolGetPermutationTarget()
71498e73e17SPierre Jolivet @*/
71598e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat A,PetscBool use)
71698e73e17SPierre Jolivet {
71798e73e17SPierre Jolivet   PetscFunctionBegin;
71898e73e17SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
71998e73e17SPierre Jolivet   PetscValidLogicalCollectiveBool(A,use,2);
720cac4c232SBarry Smith   PetscTryMethod(A,"MatHtoolUsePermutation_C",(Mat,PetscBool),(A,use));
72198e73e17SPierre Jolivet   PetscFunctionReturn(0);
72298e73e17SPierre Jolivet }
72398e73e17SPierre Jolivet 
724c7a4214aSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B)
725c7a4214aSPierre Jolivet {
726c7a4214aSPierre Jolivet   Mat            C;
727c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
728c7a4214aSPierre Jolivet   PetscInt       lda;
729c7a4214aSPierre Jolivet   PetscScalar    *array;
730c7a4214aSPierre Jolivet 
731c7a4214aSPierre Jolivet   PetscFunctionBegin;
732c7a4214aSPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
733c7a4214aSPierre Jolivet     C = *B;
7345f80ce2aSJacob Faibussowitsch     PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible dimensions");
7359566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(C,&lda));
7365f80ce2aSJacob Faibussowitsch     PetscCheck(lda == C->rmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")",lda,C->rmap->n);
737c7a4214aSPierre Jolivet   } else {
7389566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C));
7399566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N));
7409566063dSJacob Faibussowitsch     PetscCall(MatSetType(C,MATDENSE));
7419566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
7429566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE));
743c7a4214aSPierre Jolivet   }
7449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
7459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
7469566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C,&array));
747c7a4214aSPierre Jolivet   a->hmatrix->copy_local_dense_perm(array);
7489566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C,&array));
7499566063dSJacob Faibussowitsch   PetscCall(MatScale(C,a->s));
750c7a4214aSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
7519566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&C));
752c7a4214aSPierre Jolivet   } else *B = C;
753c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
754c7a4214aSPierre Jolivet }
755c7a4214aSPierre Jolivet 
75698e73e17SPierre Jolivet static PetscErrorCode GenEntriesTranspose(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *rows,const PetscInt *cols,PetscScalar *ptr,void *ctx)
757c7a4214aSPierre Jolivet {
758c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose*)ctx;
75998e73e17SPierre Jolivet   PetscScalar             *tmp;
76098e73e17SPierre Jolivet 
76198e73e17SPierre Jolivet   PetscFunctionBegin;
76298e73e17SPierre Jolivet   generator->kernel(sdim,N,M,cols,rows,ptr,generator->kernelctx);
7639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M*N,&tmp));
7649566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tmp,ptr,M*N));
76598e73e17SPierre Jolivet   for (PetscInt i=0; i<M; ++i) {
76698e73e17SPierre Jolivet     for (PetscInt j=0; j<N; ++j) ptr[i+j*M] = tmp[j+i*N];
76798e73e17SPierre Jolivet   }
7689566063dSJacob Faibussowitsch   PetscCall(PetscFree(tmp));
76998e73e17SPierre Jolivet   PetscFunctionReturn(0);
770c7a4214aSPierre Jolivet }
771c7a4214aSPierre Jolivet 
772c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */
773c7a4214aSPierre Jolivet static PetscErrorCode MatTranspose_Htool(Mat A,MatReuse reuse,Mat *B)
774c7a4214aSPierre Jolivet {
775c7a4214aSPierre Jolivet   Mat                     C;
776c7a4214aSPierre Jolivet   Mat_Htool               *a = (Mat_Htool*)A->data,*c;
777c7a4214aSPierre Jolivet   PetscInt                M = A->rmap->N,N = A->cmap->N,m = A->rmap->n,n = A->cmap->n;
778c7a4214aSPierre Jolivet   PetscContainer          container;
779c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
780c7a4214aSPierre Jolivet 
781c7a4214aSPierre Jolivet   PetscFunctionBegin;
7825f80ce2aSJacob Faibussowitsch   PetscCheck(reuse != MAT_INPLACE_MATRIX,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatTranspose() with MAT_INPLACE_MATRIX not supported");
783c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
7849566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C));
7859566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C,n,m,N,M));
7869566063dSJacob Faibussowitsch     PetscCall(MatSetType(C,((PetscObject)A)->type_name));
7879566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
7889566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)C),&container));
7899566063dSJacob Faibussowitsch     PetscCall(PetscNew(&kernelt));
7909566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(container,kernelt));
7919566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)C,"KernelTranspose",(PetscObject)container));
792c7a4214aSPierre Jolivet   } else {
793c7a4214aSPierre Jolivet     C = *B;
7949566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)C,"KernelTranspose",(PetscObject*)&container));
7955f80ce2aSJacob Faibussowitsch     PetscCheck(container,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatTranspose() with MAT_INITIAL_MATRIX first");
7969566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container,(void**)&kernelt));
797c7a4214aSPierre Jolivet   }
798c7a4214aSPierre Jolivet   c                  = (Mat_Htool*)C->data;
799c7a4214aSPierre Jolivet   c->dim             = a->dim;
800c7a4214aSPierre Jolivet   c->s               = a->s;
80198e73e17SPierre Jolivet   c->kernel          = GenEntriesTranspose;
802c7a4214aSPierre Jolivet   if (kernelt->A != A) {
8039566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
804c7a4214aSPierre Jolivet     kernelt->A       = A;
8059566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
806c7a4214aSPierre Jolivet   }
807c7a4214aSPierre Jolivet   kernelt->kernel    = a->kernel;
808c7a4214aSPierre Jolivet   kernelt->kernelctx = a->kernelctx;
809c7a4214aSPierre Jolivet   c->kernelctx       = kernelt;
810c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
8119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(N*c->dim,&c->gcoords_target));
8129566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(c->gcoords_target,a->gcoords_source,N*c->dim));
813c7a4214aSPierre Jolivet     if (a->gcoords_target != a->gcoords_source) {
8149566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(M*c->dim,&c->gcoords_source));
8159566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->gcoords_source,a->gcoords_target,M*c->dim));
816c7a4214aSPierre Jolivet     } else c->gcoords_source = c->gcoords_target;
8179566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(M,&c->work_source,N,&c->work_target));
818c7a4214aSPierre Jolivet   }
8199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
8209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
821c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) *B = C;
822c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
823c7a4214aSPierre Jolivet }
824c7a4214aSPierre Jolivet 
825c7a4214aSPierre Jolivet /*@C
826c7a4214aSPierre Jolivet      MatCreateHtoolFromKernel - Creates a MATHTOOL from a user-supplied kernel.
827c7a4214aSPierre Jolivet 
828c7a4214aSPierre Jolivet    Input Parameters:
829c7a4214aSPierre Jolivet +     comm - MPI communicator
830c7a4214aSPierre Jolivet .     m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
831c7a4214aSPierre Jolivet .     n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
832c7a4214aSPierre Jolivet .     M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
833c7a4214aSPierre Jolivet .     N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
834c7a4214aSPierre Jolivet .     spacedim - dimension of the space coordinates
835c7a4214aSPierre Jolivet .     coords_target - coordinates of the target
836c7a4214aSPierre Jolivet .     coords_source - coordinates of the source
837c7a4214aSPierre Jolivet .     kernel - computational kernel (or NULL)
838cab00e0dSPierre Jolivet -     kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
839c7a4214aSPierre Jolivet 
840c7a4214aSPierre Jolivet    Output Parameter:
841c7a4214aSPierre Jolivet .     B - matrix
842c7a4214aSPierre Jolivet 
843c7a4214aSPierre Jolivet    Options Database Keys:
844c7a4214aSPierre Jolivet +     -mat_htool_min_cluster_size <PetscInt> - minimal leaf size in cluster tree
845c7a4214aSPierre Jolivet .     -mat_htool_max_block_size <PetscInt> - maximal number of coefficients in a dense block
846c7a4214aSPierre Jolivet .     -mat_htool_epsilon <PetscReal> - relative error in Frobenius norm when approximating a block
847c7a4214aSPierre Jolivet .     -mat_htool_eta <PetscReal> - admissibility condition tolerance
848c7a4214aSPierre Jolivet .     -mat_htool_min_target_depth <PetscInt> - minimal cluster tree depth associated with the rows
849c7a4214aSPierre Jolivet .     -mat_htool_min_source_depth <PetscInt> - minimal cluster tree depth associated with the columns
85098e73e17SPierre Jolivet .     -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression
85198e73e17SPierre Jolivet -     -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering
852c7a4214aSPierre Jolivet 
853c7a4214aSPierre Jolivet    Level: intermediate
854c7a4214aSPierre Jolivet 
85553022affSStefano Zampini .seealso:  MatCreate(), MATHTOOL, PCSetCoordinates(), MatHtoolSetKernel(), MatHtoolCompressorType, MATH2OPUS, MatCreateH2OpusFromKernel()
856c7a4214aSPierre Jolivet @*/
857c7a4214aSPierre Jolivet PetscErrorCode MatCreateHtoolFromKernel(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt spacedim,const PetscReal coords_target[],const PetscReal coords_source[],MatHtoolKernel kernel,void *kernelctx,Mat *B)
858c7a4214aSPierre Jolivet {
859c7a4214aSPierre Jolivet   Mat            A;
860c7a4214aSPierre Jolivet   Mat_Htool      *a;
861c7a4214aSPierre Jolivet 
862c7a4214aSPierre Jolivet   PetscFunctionBegin;
8639566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,&A));
864c7a4214aSPierre Jolivet   PetscValidLogicalCollectiveInt(A,spacedim,6);
865c7a4214aSPierre Jolivet   PetscValidRealPointer(coords_target,7);
866c7a4214aSPierre Jolivet   PetscValidRealPointer(coords_source,8);
867c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel,9);
868c7a4214aSPierre Jolivet   if (!kernel)    PetscValidPointer(kernelctx,10);
8699566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,m,n,M,N));
8709566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATHTOOL));
8719566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
872c7a4214aSPierre Jolivet   a            = (Mat_Htool*)A->data;
873c7a4214aSPierre Jolivet   a->dim       = spacedim;
874c7a4214aSPierre Jolivet   a->s         = 1.0;
875c7a4214aSPierre Jolivet   a->kernel    = kernel;
876c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
8779566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(A->rmap->N*spacedim,&a->gcoords_target));
8789566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->gcoords_target+A->rmap->rstart*spacedim,coords_target,A->rmap->n*spacedim));
8791c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(MPI_IN_PLACE,a->gcoords_target,A->rmap->N*spacedim,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)A))); /* global target coordinates */
880c7a4214aSPierre Jolivet   if (coords_target != coords_source) {
8819566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(A->cmap->N*spacedim,&a->gcoords_source));
8829566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->gcoords_source+A->cmap->rstart*spacedim,coords_source,A->cmap->n*spacedim));
8831c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(MPI_IN_PLACE,a->gcoords_source,A->cmap->N*spacedim,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)A))); /* global source coordinates */
884c7a4214aSPierre Jolivet   } else a->gcoords_source = a->gcoords_target;
8859566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(A->cmap->N,&a->work_source,A->rmap->N,&a->work_target));
886c7a4214aSPierre Jolivet   *B = A;
887c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
888c7a4214aSPierre Jolivet }
889c7a4214aSPierre Jolivet 
890c7a4214aSPierre Jolivet /*MC
891c7a4214aSPierre Jolivet      MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package.
892c7a4214aSPierre Jolivet 
893c7a4214aSPierre Jolivet   Use ./configure --download-htool to install PETSc to use Htool.
894c7a4214aSPierre Jolivet 
895c7a4214aSPierre Jolivet    Options Database Keys:
896e85a1b63SPatrick Sanan .     -mat_type htool - matrix type to "htool" during a call to MatSetFromOptions()
897c7a4214aSPierre Jolivet 
898c7a4214aSPierre Jolivet    Level: beginner
899c7a4214aSPierre Jolivet 
90053022affSStefano Zampini .seealso: MATH2OPUS, MATDENSE, MatCreateHtoolFromKernel(), MatHtoolSetKernel()
901c7a4214aSPierre Jolivet M*/
902c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A)
903c7a4214aSPierre Jolivet {
904c7a4214aSPierre Jolivet   Mat_Htool *a;
905c7a4214aSPierre Jolivet 
906c7a4214aSPierre Jolivet   PetscFunctionBegin;
9079566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(A,&a));
908c7a4214aSPierre Jolivet   A->data = (void*)a;
9099566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATHTOOL));
9109566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(A->ops,sizeof(struct _MatOps)));
911c7a4214aSPierre Jolivet   A->ops->getdiagonal       = MatGetDiagonal_Htool;
912c7a4214aSPierre Jolivet   A->ops->getdiagonalblock  = MatGetDiagonalBlock_Htool;
913c7a4214aSPierre Jolivet   A->ops->mult              = MatMult_Htool;
914c7a4214aSPierre Jolivet   A->ops->multadd           = MatMultAdd_Htool;
9158b8ddd11SPierre Jolivet   A->ops->multtranspose     = MatMultTranspose_Htool;
9168b8ddd11SPierre Jolivet   if (!PetscDefined(USE_COMPLEX)) A->ops->multhermitiantranspose = MatMultTranspose_Htool;
917c7a4214aSPierre Jolivet   A->ops->increaseoverlap   = MatIncreaseOverlap_Htool;
918c7a4214aSPierre Jolivet   A->ops->createsubmatrices = MatCreateSubMatrices_Htool;
919c7a4214aSPierre Jolivet   A->ops->transpose         = MatTranspose_Htool;
920c7a4214aSPierre Jolivet   A->ops->destroy           = MatDestroy_Htool;
921c7a4214aSPierre Jolivet   A->ops->view              = MatView_Htool;
922c7a4214aSPierre Jolivet   A->ops->setfromoptions    = MatSetFromOptions_Htool;
923c7a4214aSPierre Jolivet   A->ops->scale             = MatScale_Htool;
924c7a4214aSPierre Jolivet   A->ops->getrow            = MatGetRow_Htool;
925c7a4214aSPierre Jolivet   A->ops->restorerow        = MatRestoreRow_Htool;
926c7a4214aSPierre Jolivet   A->ops->assemblyend       = MatAssemblyEnd_Htool;
927c7a4214aSPierre Jolivet   a->dim                    = 0;
928c7a4214aSPierre Jolivet   a->gcoords_target         = NULL;
929c7a4214aSPierre Jolivet   a->gcoords_source         = NULL;
930c7a4214aSPierre Jolivet   a->s                      = 1.0;
931c7a4214aSPierre Jolivet   a->bs[0]                  = 10;
932c7a4214aSPierre Jolivet   a->bs[1]                  = 1000000;
933c7a4214aSPierre Jolivet   a->epsilon                = PetscSqrtReal(PETSC_SMALL);
934c7a4214aSPierre Jolivet   a->eta                    = 10.0;
935c7a4214aSPierre Jolivet   a->depth[0]               = 0;
936c7a4214aSPierre Jolivet   a->depth[1]               = 0;
937c7a4214aSPierre Jolivet   a->compressor             = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA;
9389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",MatProductSetFromOptions_Htool));
9399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",MatProductSetFromOptions_Htool));
9409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",MatConvert_Htool_Dense));
9419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",MatConvert_Htool_Dense));
9429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",MatHtoolGetHierarchicalMat_Htool));
9439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",MatHtoolSetKernel_Htool));
9449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",MatHtoolGetPermutationSource_Htool));
9459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",MatHtoolGetPermutationTarget_Htool));
9469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",MatHtoolUsePermutation_Htool));
947c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
948c7a4214aSPierre Jolivet }
949