xref: /petsc/src/mat/impls/htool/htool.cxx (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
1c7a4214aSPierre Jolivet #include <../src/mat/impls/htool/htool.hpp> /*I "petscmat.h" I*/
2c7a4214aSPierre Jolivet #include <petscblaslapack.h>
3c7a4214aSPierre Jolivet #include <set>
4c7a4214aSPierre Jolivet 
598e73e17SPierre Jolivet #define ALEN(a) (sizeof(a)/sizeof((a)[0]))
6c7a4214aSPierre Jolivet const char *const MatHtoolCompressorTypes[] = { "sympartialACA", "fullACA", "SVD" };
798e73e17SPierre Jolivet const char *const MatHtoolClusteringTypes[] = { "PCARegular", "PCAGeometric", "BoundingBox1Regular", "BoundingBox1Geometric" };
8c7a4214aSPierre Jolivet const char HtoolCitation[] = "@article{marchand2020two,\n"
9c7a4214aSPierre Jolivet "  Author = {Marchand, Pierre and Claeys, Xavier and Jolivet, Pierre and Nataf, Fr\\'ed\\'eric and Tournier, Pierre-Henri},\n"
10c7a4214aSPierre Jolivet "  Title = {Two-level preconditioning for $h$-version boundary element approximation of hypersingular operator with {GenEO}},\n"
11c7a4214aSPierre Jolivet "  Year = {2020},\n"
12c7a4214aSPierre Jolivet "  Publisher = {Elsevier},\n"
13c7a4214aSPierre Jolivet "  Journal = {Numerische Mathematik},\n"
14c7a4214aSPierre Jolivet "  Volume = {146},\n"
15c7a4214aSPierre Jolivet "  Pages = {597--628},\n"
16c7a4214aSPierre Jolivet "  Url = {https://github.com/htool-ddm/htool}\n"
17c7a4214aSPierre Jolivet "}\n";
18c7a4214aSPierre Jolivet static PetscBool HtoolCite = PETSC_FALSE;
19c7a4214aSPierre Jolivet 
20c7a4214aSPierre Jolivet static PetscErrorCode MatGetDiagonal_Htool(Mat A,Vec v)
21c7a4214aSPierre Jolivet {
22c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
23c7a4214aSPierre Jolivet   PetscScalar    *x;
24c7a4214aSPierre Jolivet   PetscBool      flg;
25c7a4214aSPierre Jolivet 
26c7a4214aSPierre Jolivet   PetscFunctionBegin;
279566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A,&flg));
285f80ce2aSJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only congruent layouts supported");
299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(v,&x));
30c7a4214aSPierre Jolivet   a->hmatrix->copy_local_diagonal(x);
319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(v,&x));
329566063dSJacob Faibussowitsch   PetscCall(VecScale(v,a->s));
33c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
34c7a4214aSPierre Jolivet }
35c7a4214aSPierre Jolivet 
36c7a4214aSPierre Jolivet static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A,Mat *b)
37c7a4214aSPierre Jolivet {
38c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
39c7a4214aSPierre Jolivet   Mat            B;
40c7a4214aSPierre Jolivet   PetscScalar    *ptr;
41c7a4214aSPierre Jolivet   PetscBool      flg;
42c7a4214aSPierre Jolivet 
43c7a4214aSPierre Jolivet   PetscFunctionBegin;
449566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A,&flg));
455f80ce2aSJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only congruent layouts supported");
469566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B)); /* same logic as in MatGetDiagonalBlock_MPIDense() */
47c7a4214aSPierre Jolivet   if (!B) {
489566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PETSC_COMM_SELF,A->rmap->n,A->rmap->n,A->rmap->n,A->rmap->n,NULL,&B));
499566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite(B,&ptr));
50c7a4214aSPierre Jolivet     a->hmatrix->copy_local_diagonal_block(ptr);
519566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite(B,&ptr));
529566063dSJacob Faibussowitsch     PetscCall(MatPropagateSymmetryOptions(A,B));
539566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
549566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
559566063dSJacob Faibussowitsch     PetscCall(MatScale(B,a->s));
569566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B));
57c7a4214aSPierre Jolivet     *b   = B;
589566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
59c7a4214aSPierre Jolivet   } else *b = B;
60c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
61c7a4214aSPierre Jolivet }
62c7a4214aSPierre Jolivet 
63c7a4214aSPierre Jolivet static PetscErrorCode MatMult_Htool(Mat A,Vec x,Vec y)
64c7a4214aSPierre Jolivet {
65c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
66c7a4214aSPierre Jolivet   const PetscScalar *in;
67c7a4214aSPierre Jolivet   PetscScalar       *out;
68c7a4214aSPierre Jolivet 
69c7a4214aSPierre Jolivet   PetscFunctionBegin;
709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&in));
719566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y,&out));
72c7a4214aSPierre Jolivet   a->hmatrix->mvprod_local_to_local(in,out);
739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&in));
749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y,&out));
759566063dSJacob Faibussowitsch   PetscCall(VecScale(y,a->s));
76c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
77c7a4214aSPierre Jolivet }
78c7a4214aSPierre Jolivet 
79c7a4214aSPierre Jolivet /* naive implementation of MatMultAdd() needed for FEM-BEM coupling via MATNEST */
80c7a4214aSPierre Jolivet static PetscErrorCode MatMultAdd_Htool(Mat A,Vec v1,Vec v2,Vec v3)
81c7a4214aSPierre Jolivet {
82c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
83c7a4214aSPierre Jolivet   Vec               tmp;
84c7a4214aSPierre Jolivet   const PetscScalar scale = a->s;
85c7a4214aSPierre Jolivet 
86c7a4214aSPierre Jolivet   PetscFunctionBegin;
879566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(v2,&tmp));
889566063dSJacob 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 */
89c7a4214aSPierre Jolivet   a->s = 1.0; /* set s to 1.0 since VecAXPY() may be used to scale the MatMult() output Vec */
909566063dSJacob Faibussowitsch   PetscCall(MatMult_Htool(A,v1,tmp));
919566063dSJacob Faibussowitsch   PetscCall(VecAXPY(v3,scale,tmp));
929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp));
93c7a4214aSPierre Jolivet   a->s = scale; /* set s back to its original value */
94c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
95c7a4214aSPierre Jolivet }
96c7a4214aSPierre Jolivet 
978b8ddd11SPierre Jolivet static PetscErrorCode MatMultTranspose_Htool(Mat A,Vec x,Vec y)
988b8ddd11SPierre Jolivet {
998b8ddd11SPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
1008b8ddd11SPierre Jolivet   const PetscScalar *in;
1018b8ddd11SPierre Jolivet   PetscScalar       *out;
1028b8ddd11SPierre Jolivet 
1038b8ddd11SPierre Jolivet   PetscFunctionBegin;
1049566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&in));
1059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y,&out));
1068b8ddd11SPierre Jolivet   a->hmatrix->mvprod_transp_local_to_local(in,out);
1079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&in));
1089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y,&out));
1099566063dSJacob Faibussowitsch   PetscCall(VecScale(y,a->s));
1108b8ddd11SPierre Jolivet   PetscFunctionReturn(0);
1118b8ddd11SPierre Jolivet }
1128b8ddd11SPierre Jolivet 
113c7a4214aSPierre Jolivet static PetscErrorCode MatIncreaseOverlap_Htool(Mat A,PetscInt is_max,IS is[],PetscInt ov)
114c7a4214aSPierre Jolivet {
115c7a4214aSPierre Jolivet   std::set<PetscInt> set;
116c7a4214aSPierre Jolivet   const PetscInt     *idx;
1178f308287SPierre Jolivet   PetscInt           *oidx,size,bs[2];
118c7a4214aSPierre Jolivet   PetscMPIInt        csize;
119c7a4214aSPierre Jolivet 
120c7a4214aSPierre Jolivet   PetscFunctionBegin;
1219566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(A,bs,bs+1));
1228f308287SPierre Jolivet   if (bs[0] != bs[1]) bs[0] = 1;
123c7a4214aSPierre Jolivet   for (PetscInt i=0; i<is_max; ++i) {
124c7a4214aSPierre Jolivet     /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */
125c7a4214aSPierre Jolivet     /* needed to avoid subdomain matrices to replicate A since it is dense                           */
1269566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]),&csize));
1275f80ce2aSJacob Faibussowitsch     PetscCheck(csize == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported parallel IS");
1289566063dSJacob Faibussowitsch     PetscCall(ISGetSize(is[i],&size));
1299566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is[i],&idx));
130c7a4214aSPierre Jolivet     for (PetscInt j=0; j<size; ++j) {
131c7a4214aSPierre Jolivet       set.insert(idx[j]);
132c7a4214aSPierre Jolivet       for (PetscInt k=1; k<=ov; ++k) {               /* for each layer of overlap      */
133c7a4214aSPierre Jolivet         if (idx[j] - k >= 0) set.insert(idx[j] - k); /* do not insert negative indices */
134c7a4214aSPierre 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 */
135c7a4214aSPierre Jolivet       }
136c7a4214aSPierre Jolivet     }
1379566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is[i],&idx));
1389566063dSJacob Faibussowitsch     PetscCall(ISDestroy(is+i));
1398f308287SPierre Jolivet     if (bs[0] > 1) {
1408f308287SPierre Jolivet       for (std::set<PetscInt>::iterator it=set.cbegin(); it!=set.cend(); it++) {
1418f308287SPierre Jolivet         std::vector<PetscInt> block(bs[0]);
1428f308287SPierre Jolivet         std::iota(block.begin(),block.end(),(*it/bs[0])*bs[0]);
1438f308287SPierre Jolivet         set.insert(block.cbegin(),block.cend());
1448f308287SPierre Jolivet       }
1458f308287SPierre Jolivet     }
146c7a4214aSPierre Jolivet     size = set.size(); /* size with overlap */
1479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size,&oidx));
148c7a4214aSPierre Jolivet     for (const PetscInt j : set) *oidx++ = j;
149c7a4214aSPierre Jolivet     oidx -= size;
1509566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,size,oidx,PETSC_OWN_POINTER,is+i));
151c7a4214aSPierre Jolivet   }
152c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
153c7a4214aSPierre Jolivet }
154c7a4214aSPierre Jolivet 
155c7a4214aSPierre Jolivet static PetscErrorCode MatCreateSubMatrices_Htool(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
156c7a4214aSPierre Jolivet {
157c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
158c7a4214aSPierre Jolivet   Mat               D,B,BT;
159c7a4214aSPierre Jolivet   const PetscScalar *copy;
160c7a4214aSPierre Jolivet   PetscScalar       *ptr;
161c7a4214aSPierre Jolivet   const PetscInt    *idxr,*idxc,*it;
162c7a4214aSPierre Jolivet   PetscInt          nrow,m,i;
163c7a4214aSPierre Jolivet   PetscBool         flg;
164c7a4214aSPierre Jolivet 
165c7a4214aSPierre Jolivet   PetscFunctionBegin;
166c7a4214aSPierre Jolivet   if (scall != MAT_REUSE_MATRIX) {
1679566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(n,submat));
168c7a4214aSPierre Jolivet   }
169c7a4214aSPierre Jolivet   for (i=0; i<n; ++i) {
1709566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(irow[i],&nrow));
1719566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(icol[i],&m));
1729566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(irow[i],&idxr));
1739566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(icol[i],&idxc));
174c7a4214aSPierre Jolivet     if (scall != MAT_REUSE_MATRIX) {
1759566063dSJacob Faibussowitsch       PetscCall(MatCreateDense(PETSC_COMM_SELF,nrow,m,nrow,m,NULL,(*submat)+i));
176c7a4214aSPierre Jolivet     }
1779566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite((*submat)[i],&ptr));
178c7a4214aSPierre Jolivet     if (irow[i] == icol[i]) { /* same row and column IS? */
1799566063dSJacob Faibussowitsch       PetscCall(MatHasCongruentLayouts(A,&flg));
180c7a4214aSPierre Jolivet       if (flg) {
1819566063dSJacob Faibussowitsch         PetscCall(ISSorted(irow[i],&flg));
182c7a4214aSPierre Jolivet         if (flg) { /* sorted IS? */
183c7a4214aSPierre Jolivet           it = std::lower_bound(idxr,idxr+nrow,A->rmap->rstart);
184c7a4214aSPierre Jolivet           if (it != idxr+nrow && *it == A->rmap->rstart) { /* rmap->rstart in IS? */
185c7a4214aSPierre Jolivet             if (std::distance(idxr,it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */
186c7a4214aSPierre Jolivet               for (PetscInt j=0; j<A->rmap->n && flg; ++j) if (PetscUnlikely(it[j] != A->rmap->rstart+j)) flg = PETSC_FALSE;
187c7a4214aSPierre Jolivet               if (flg) { /* complete local diagonal block in IS? */
188c7a4214aSPierre Jolivet                 /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM
189c7a4214aSPierre Jolivet                  *      [   B   C   E   ]
190c7a4214aSPierre Jolivet                  *  A = [   B   D   E   ]
191c7a4214aSPierre Jolivet                  *      [   B   F   E   ]
192c7a4214aSPierre Jolivet                  */
193c7a4214aSPierre Jolivet                 m = std::distance(idxr,it); /* shift of the coefficient (0,0) of block D from above */
1949566063dSJacob Faibussowitsch                 PetscCall(MatGetDiagonalBlock_Htool(A,&D));
1959566063dSJacob Faibussowitsch                 PetscCall(MatDenseGetArrayRead(D,&copy));
196c7a4214aSPierre Jolivet                 for (PetscInt k=0; k<A->rmap->n; ++k) {
1979566063dSJacob Faibussowitsch                   PetscCall(PetscArraycpy(ptr+(m+k)*nrow+m,copy+k*A->rmap->n,A->rmap->n)); /* block D from above */
198c7a4214aSPierre Jolivet                 }
1999566063dSJacob Faibussowitsch                 PetscCall(MatDenseRestoreArrayRead(D,&copy));
200c7a4214aSPierre Jolivet                 if (m) {
201c7a4214aSPierre Jolivet                   a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* vertical block B from above */
202c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
203c7a4214aSPierre Jolivet                   if (A->symmetric || A->hermitian) {
2049566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF,A->rmap->n,m,A->rmap->n,m,ptr+m,&B));
2059566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B,nrow));
2069566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF,m,A->rmap->n,m,A->rmap->n,ptr+m*nrow,&BT));
2079566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT,nrow));
208c7a4214aSPierre Jolivet                     if (A->hermitian && PetscDefined(USE_COMPLEX)) {
2099566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT));
210c7a4214aSPierre Jolivet                     } else {
2119566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B,MAT_REUSE_MATRIX,&BT));
212c7a4214aSPierre Jolivet                     }
2139566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2149566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
215c7a4214aSPierre Jolivet                   } else {
216c7a4214aSPierre Jolivet                     for (PetscInt k=0; k<A->rmap->n; ++k) { /* block C from above */
217c7a4214aSPierre Jolivet                       a->wrapper->copy_submatrix(m,1,idxr,idxc+m+k,ptr+(m+k)*nrow);
218c7a4214aSPierre Jolivet                     }
219c7a4214aSPierre Jolivet                   }
220c7a4214aSPierre Jolivet                 }
221c7a4214aSPierre Jolivet                 if (m+A->rmap->n != nrow) {
222c7a4214aSPierre 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 */
223c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
224c7a4214aSPierre Jolivet                   if (A->symmetric || A->hermitian) {
2259566063dSJacob 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));
2269566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B,nrow));
2279566063dSJacob 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));
2289566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT,nrow));
229c7a4214aSPierre Jolivet                     if (A->hermitian && PetscDefined(USE_COMPLEX)) {
2309566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT));
231c7a4214aSPierre Jolivet                     } else {
2329566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B,MAT_REUSE_MATRIX,&BT));
233c7a4214aSPierre Jolivet                     }
2349566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2359566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
236c7a4214aSPierre Jolivet                   } else {
237c7a4214aSPierre Jolivet                     for (PetscInt k=0; k<A->rmap->n; ++k) { /* block F from above */
238c7a4214aSPierre 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);
239c7a4214aSPierre Jolivet                     }
240c7a4214aSPierre Jolivet                   }
241c7a4214aSPierre Jolivet                 }
242c7a4214aSPierre Jolivet               } /* complete local diagonal block not in IS */
243c7a4214aSPierre Jolivet             } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */
244c7a4214aSPierre Jolivet           } else flg = PETSC_FALSE; /* rmap->rstart not in IS */
245c7a4214aSPierre Jolivet         } /* unsorted IS */
246c7a4214aSPierre Jolivet       }
247c7a4214aSPierre Jolivet     } else flg = PETSC_FALSE; /* different row and column IS */
248c7a4214aSPierre Jolivet     if (!flg) a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* reassemble everything */
2499566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(irow[i],&idxr));
2509566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(icol[i],&idxc));
2519566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite((*submat)[i],&ptr));
2529566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin((*submat)[i],MAT_FINAL_ASSEMBLY));
2539566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd((*submat)[i],MAT_FINAL_ASSEMBLY));
2549566063dSJacob Faibussowitsch     PetscCall(MatScale((*submat)[i],a->s));
255c7a4214aSPierre Jolivet   }
256c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
257c7a4214aSPierre Jolivet }
258c7a4214aSPierre Jolivet 
259c7a4214aSPierre Jolivet static PetscErrorCode MatDestroy_Htool(Mat A)
260c7a4214aSPierre Jolivet {
261c7a4214aSPierre Jolivet   Mat_Htool               *a = (Mat_Htool*)A->data;
262c7a4214aSPierre Jolivet   PetscContainer          container;
263c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
264c7a4214aSPierre Jolivet 
265c7a4214aSPierre Jolivet   PetscFunctionBegin;
2669566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A,NULL));
2679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",NULL));
2689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",NULL));
2699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",NULL));
2709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",NULL));
2719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",NULL));
2729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",NULL));
2739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",NULL));
2749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",NULL));
2759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",NULL));
2769566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A,"KernelTranspose",(PetscObject*)&container));
277c7a4214aSPierre Jolivet   if (container) { /* created in MatTranspose_Htool() */
2789566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container,(void**)&kernelt));
2799566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
2809566063dSJacob Faibussowitsch     PetscCall(PetscFree(kernelt));
2819566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&container));
2829566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)A,"KernelTranspose",NULL));
283c7a4214aSPierre Jolivet   }
284c7a4214aSPierre Jolivet   if (a->gcoords_source != a->gcoords_target) {
2859566063dSJacob Faibussowitsch     PetscCall(PetscFree(a->gcoords_source));
286c7a4214aSPierre Jolivet   }
2879566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->gcoords_target));
2889566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a->work_source,a->work_target));
289c7a4214aSPierre Jolivet   delete a->wrapper;
290c7a4214aSPierre Jolivet   delete a->hmatrix;
2919566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
292c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
293c7a4214aSPierre Jolivet }
294c7a4214aSPierre Jolivet 
295c7a4214aSPierre Jolivet static PetscErrorCode MatView_Htool(Mat A,PetscViewer pv)
296c7a4214aSPierre Jolivet {
297c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
298c7a4214aSPierre Jolivet   PetscBool      flg;
299c7a4214aSPierre Jolivet 
300c7a4214aSPierre Jolivet   PetscFunctionBegin;
3019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pv,PETSCVIEWERASCII,&flg));
302c7a4214aSPierre Jolivet   if (flg) {
3039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"symmetry: %c\n",a->hmatrix->get_symmetry_type()));
304c7a4214aSPierre Jolivet     if (PetscAbsScalar(a->s-1.0) > PETSC_MACHINE_EPSILON) {
305c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX)
3069566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(pv,"scaling: %g+%gi\n",(double)PetscRealPart(a->s),(double)PetscImaginaryPart(a->s)));
307c7a4214aSPierre Jolivet #else
3089566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(pv,"scaling: %g\n",(double)a->s));
309c7a4214aSPierre Jolivet #endif
310c7a4214aSPierre Jolivet     }
3119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"minimum cluster size: %" PetscInt_FMT "\n",a->bs[0]));
3129566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"maximum block size: %" PetscInt_FMT "\n",a->bs[1]));
3139566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"epsilon: %g\n",(double)a->epsilon));
3149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"eta: %g\n",(double)a->eta));
3159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"minimum target depth: %" PetscInt_FMT "\n",a->depth[0]));
3169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"minimum source depth: %" PetscInt_FMT "\n",a->depth[1]));
3179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"compressor: %s\n",MatHtoolCompressorTypes[a->compressor]));
3189566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"clustering: %s\n",MatHtoolClusteringTypes[a->clustering]));
3199566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"compression ratio: %s\n",a->hmatrix->get_infos("Compression_ratio").c_str()));
3209566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"space saving: %s\n",a->hmatrix->get_infos("Space_saving").c_str()));
3219566063dSJacob 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()));
3229566063dSJacob 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()));
3239566063dSJacob 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()));
3249566063dSJacob 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()));
325c7a4214aSPierre Jolivet   }
326c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
327c7a4214aSPierre Jolivet }
328c7a4214aSPierre Jolivet 
329c7a4214aSPierre Jolivet static PetscErrorCode MatScale_Htool(Mat A,PetscScalar s)
330c7a4214aSPierre Jolivet {
331c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
332c7a4214aSPierre Jolivet 
333c7a4214aSPierre Jolivet   PetscFunctionBegin;
334c7a4214aSPierre Jolivet   a->s *= s;
335c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
336c7a4214aSPierre Jolivet }
337c7a4214aSPierre Jolivet 
338c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */
339c7a4214aSPierre Jolivet static PetscErrorCode MatGetRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
340c7a4214aSPierre Jolivet {
341c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
342c7a4214aSPierre Jolivet   PetscInt       *idxc;
343c7a4214aSPierre Jolivet   PetscBLASInt   one = 1,bn;
344c7a4214aSPierre Jolivet 
345c7a4214aSPierre Jolivet   PetscFunctionBegin;
346c7a4214aSPierre Jolivet   if (nz) *nz = A->cmap->N;
347c7a4214aSPierre Jolivet   if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */
3489566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N,&idxc));
349c7a4214aSPierre Jolivet     for (PetscInt i=0; i<A->cmap->N; ++i) idxc[i] = i;
350c7a4214aSPierre Jolivet   }
351c7a4214aSPierre Jolivet   if (idx) *idx = idxc;
352c7a4214aSPierre Jolivet   if (v) {
3539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N,v));
354c7a4214aSPierre Jolivet     if (a->wrapper) a->wrapper->copy_submatrix(1,A->cmap->N,&row,idxc,*v);
355cab00e0dSPierre Jolivet     else reinterpret_cast<htool::VirtualGenerator<PetscScalar>*>(a->kernelctx)->copy_submatrix(1,A->cmap->N,&row,idxc,*v);
3569566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(A->cmap->N,&bn));
357c7a4214aSPierre Jolivet     PetscStackCallBLAS("BLASscal",BLASscal_(&bn,&a->s,*v,&one));
358c7a4214aSPierre Jolivet   }
359c7a4214aSPierre Jolivet   if (!idx) {
3609566063dSJacob Faibussowitsch     PetscCall(PetscFree(idxc));
361c7a4214aSPierre Jolivet   }
362c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
363c7a4214aSPierre Jolivet }
364c7a4214aSPierre Jolivet 
365c7a4214aSPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
366c7a4214aSPierre Jolivet {
367c7a4214aSPierre Jolivet   PetscFunctionBegin;
368c7a4214aSPierre Jolivet   if (nz) *nz = 0;
369c7a4214aSPierre Jolivet   if (idx) {
3709566063dSJacob Faibussowitsch     PetscCall(PetscFree(*idx));
371c7a4214aSPierre Jolivet   }
372c7a4214aSPierre Jolivet   if (v) {
3739566063dSJacob Faibussowitsch     PetscCall(PetscFree(*v));
374c7a4214aSPierre Jolivet   }
375c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
376c7a4214aSPierre Jolivet }
377c7a4214aSPierre Jolivet 
378c7a4214aSPierre Jolivet static PetscErrorCode MatSetFromOptions_Htool(PetscOptionItems *PetscOptionsObject,Mat A)
379c7a4214aSPierre Jolivet {
380c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
381c7a4214aSPierre Jolivet   PetscInt       n;
382c7a4214aSPierre Jolivet   PetscBool      flg;
383c7a4214aSPierre Jolivet 
384c7a4214aSPierre Jolivet   PetscFunctionBegin;
3859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"Htool options"));
3869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_min_cluster_size","Minimal leaf size in cluster tree",NULL,a->bs[0],a->bs,NULL));
3879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_max_block_size","Maximal number of coefficients in a dense block",NULL,a->bs[1],a->bs + 1,NULL));
3889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_htool_epsilon","Relative error in Frobenius norm when approximating a block",NULL,a->epsilon,&a->epsilon,NULL));
3899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_htool_eta","Admissibility condition tolerance",NULL,a->eta,&a->eta,NULL));
3909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_min_target_depth","Minimal cluster tree depth associated with the rows",NULL,a->depth[0],a->depth,NULL));
3919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_min_source_depth","Minimal cluster tree depth associated with the columns",NULL,a->depth[1],a->depth + 1,NULL));
392c7a4214aSPierre Jolivet   n = 0;
3939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-mat_htool_compressor","Type of compression","MatHtoolCompressorType",MatHtoolCompressorTypes,ALEN(MatHtoolCompressorTypes),MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA],&n,&flg));
394c7a4214aSPierre Jolivet   if (flg) a->compressor = MatHtoolCompressorType(n);
39598e73e17SPierre Jolivet   n = 0;
3969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-mat_htool_clustering","Type of clustering","MatHtoolClusteringType",MatHtoolClusteringTypes,ALEN(MatHtoolClusteringTypes),MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR],&n,&flg));
39798e73e17SPierre Jolivet   if (flg) a->clustering = MatHtoolClusteringType(n);
3989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
399c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
400c7a4214aSPierre Jolivet }
401c7a4214aSPierre Jolivet 
402c7a4214aSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A,MatAssemblyType type)
403c7a4214aSPierre Jolivet {
404c7a4214aSPierre Jolivet   Mat_Htool                                                    *a = (Mat_Htool*)A->data;
405c7a4214aSPierre Jolivet   const PetscInt                                               *ranges;
406c7a4214aSPierre Jolivet   PetscInt                                                     *offset;
407c7a4214aSPierre Jolivet   PetscMPIInt                                                  size;
408c7a4214aSPierre Jolivet   char                                                         S = PetscDefined(USE_COMPLEX) && A->hermitian ? 'H' : (A->symmetric ? 'S' : 'N'),uplo = S == 'N' ? 'N' : 'U';
409cab00e0dSPierre Jolivet   htool::VirtualGenerator<PetscScalar>                         *generator = nullptr;
41098e73e17SPierre Jolivet   std::shared_ptr<htool::VirtualCluster>                       t,s = nullptr;
4113b9338faSPierre Jolivet   std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor = nullptr;
412c7a4214aSPierre Jolivet 
413c7a4214aSPierre Jolivet   PetscFunctionBegin;
4149566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(HtoolCitation,&HtoolCite));
415c7a4214aSPierre Jolivet   delete a->wrapper;
416c7a4214aSPierre Jolivet   delete a->hmatrix;
4179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
4189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*size,&offset));
4199566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(A,&ranges));
420c7a4214aSPierre Jolivet   for (PetscInt i=0; i<size; ++i) {
421c7a4214aSPierre Jolivet     offset[2*i] = ranges[i];
422c7a4214aSPierre Jolivet     offset[2*i+1] = ranges[i+1] - ranges[i];
423c7a4214aSPierre Jolivet   }
42498e73e17SPierre Jolivet   switch (a->clustering) {
42598e73e17SPierre Jolivet   case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
42698e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
42798e73e17SPierre Jolivet     break;
42898e73e17SPierre Jolivet   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
42998e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
43098e73e17SPierre Jolivet     break;
43198e73e17SPierre Jolivet   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
43298e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim);
43398e73e17SPierre Jolivet     break;
43498e73e17SPierre Jolivet   default:
43598e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim);
43698e73e17SPierre Jolivet   }
437c7a4214aSPierre Jolivet   t->set_minclustersize(a->bs[0]);
43898e73e17SPierre Jolivet   t->build(A->rmap->N,a->gcoords_target,offset);
439c7a4214aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx);
440c7a4214aSPierre Jolivet   else {
441c7a4214aSPierre Jolivet     a->wrapper = NULL;
442cab00e0dSPierre Jolivet     generator = reinterpret_cast<htool::VirtualGenerator<PetscScalar>*>(a->kernelctx);
443c7a4214aSPierre Jolivet   }
444c7a4214aSPierre Jolivet   if (a->gcoords_target != a->gcoords_source) {
4459566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangesColumn(A,&ranges));
446c7a4214aSPierre Jolivet     for (PetscInt i=0; i<size; ++i) {
447c7a4214aSPierre Jolivet       offset[2*i] = ranges[i];
448c7a4214aSPierre Jolivet       offset[2*i+1] = ranges[i+1] - ranges[i];
449c7a4214aSPierre Jolivet     }
45098e73e17SPierre Jolivet     switch (a->clustering) {
45198e73e17SPierre Jolivet     case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
45298e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
45398e73e17SPierre Jolivet       break;
45498e73e17SPierre Jolivet     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
45598e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
45698e73e17SPierre Jolivet       break;
45798e73e17SPierre Jolivet     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
45898e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim);
45998e73e17SPierre Jolivet       break;
46098e73e17SPierre Jolivet     default:
46198e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim);
46298e73e17SPierre Jolivet     }
463c7a4214aSPierre Jolivet     s->set_minclustersize(a->bs[0]);
46498e73e17SPierre Jolivet     s->build(A->cmap->N,a->gcoords_source,offset);
465c7a4214aSPierre Jolivet     S = uplo = 'N';
466c7a4214aSPierre Jolivet   }
4679566063dSJacob Faibussowitsch   PetscCall(PetscFree(offset));
468c7a4214aSPierre Jolivet   switch (a->compressor) {
469c7a4214aSPierre Jolivet   case MAT_HTOOL_COMPRESSOR_FULL_ACA:
4703b9338faSPierre Jolivet     compressor = std::make_shared<htool::fullACA<PetscScalar>>();
471c7a4214aSPierre Jolivet     break;
472c7a4214aSPierre Jolivet   case MAT_HTOOL_COMPRESSOR_SVD:
4733b9338faSPierre Jolivet     compressor = std::make_shared<htool::SVD<PetscScalar>>();
474c7a4214aSPierre Jolivet     break;
475c7a4214aSPierre Jolivet   default:
4763b9338faSPierre Jolivet     compressor = std::make_shared<htool::sympartialACA<PetscScalar>>();
477c7a4214aSPierre Jolivet   }
4783b9338faSPierre Jolivet   a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar>*>(new htool::HMatrix<PetscScalar>(t,s ? s : t,a->epsilon,a->eta,S,uplo));
4793b9338faSPierre Jolivet   a->hmatrix->set_compression(compressor);
480c7a4214aSPierre Jolivet   a->hmatrix->set_maxblocksize(a->bs[1]);
481c7a4214aSPierre Jolivet   a->hmatrix->set_mintargetdepth(a->depth[0]);
482c7a4214aSPierre Jolivet   a->hmatrix->set_minsourcedepth(a->depth[1]);
4833b9338faSPierre Jolivet   if (s) a->hmatrix->build(a->wrapper ? *a->wrapper : *generator,a->gcoords_target,a->gcoords_source);
4843b9338faSPierre Jolivet   else   a->hmatrix->build(a->wrapper ? *a->wrapper : *generator,a->gcoords_target);
485c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
486c7a4214aSPierre Jolivet }
487c7a4214aSPierre Jolivet 
488c7a4214aSPierre Jolivet static PetscErrorCode MatProductNumeric_Htool(Mat C)
489c7a4214aSPierre Jolivet {
490c7a4214aSPierre Jolivet   Mat_Product       *product = C->product;
491c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)product->A->data;
492c7a4214aSPierre Jolivet   const PetscScalar *in;
493c7a4214aSPierre Jolivet   PetscScalar       *out;
4948b8ddd11SPierre Jolivet   PetscInt          N,lda;
495c7a4214aSPierre Jolivet 
496c7a4214aSPierre Jolivet   PetscFunctionBegin;
497c7a4214aSPierre Jolivet   MatCheckProduct(C,1);
4989566063dSJacob Faibussowitsch   PetscCall(MatGetSize(C,NULL,&N));
4999566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(C,&lda));
5005f80ce2aSJacob Faibussowitsch   PetscCheck(lda == C->rmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")",lda,C->rmap->n);
5019566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(product->B,&in));
5029566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C,&out));
5038b8ddd11SPierre Jolivet   switch (product->type) {
5048b8ddd11SPierre Jolivet   case MATPRODUCT_AB:
505c7a4214aSPierre Jolivet     a->hmatrix->mvprod_local_to_local(in,out,N);
5068b8ddd11SPierre Jolivet     break;
5078b8ddd11SPierre Jolivet   case MATPRODUCT_AtB:
5088b8ddd11SPierre Jolivet     a->hmatrix->mvprod_transp_local_to_local(in,out,N);
509c7a4214aSPierre Jolivet     break;
510c7a4214aSPierre Jolivet   default:
51198921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType %s is not supported",MatProductTypes[product->type]);
512c7a4214aSPierre Jolivet   }
5139566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C,&out));
5149566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(product->B,&in));
5159566063dSJacob Faibussowitsch   PetscCall(MatScale(C,a->s));
516c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
517c7a4214aSPierre Jolivet }
518c7a4214aSPierre Jolivet 
519c7a4214aSPierre Jolivet static PetscErrorCode MatProductSymbolic_Htool(Mat C)
520c7a4214aSPierre Jolivet {
521c7a4214aSPierre Jolivet   Mat_Product    *product = C->product;
522c7a4214aSPierre Jolivet   Mat            A,B;
523c7a4214aSPierre Jolivet   PetscBool      flg;
524c7a4214aSPierre Jolivet 
525c7a4214aSPierre Jolivet   PetscFunctionBegin;
526c7a4214aSPierre Jolivet   MatCheckProduct(C,1);
527c7a4214aSPierre Jolivet   A = product->A;
528c7a4214aSPierre Jolivet   B = product->B;
5299566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,""));
5305f80ce2aSJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MatProduct_AB not supported for %s",((PetscObject)product->B)->type_name);
531c7a4214aSPierre Jolivet   switch (product->type) {
532c7a4214aSPierre Jolivet   case MATPRODUCT_AB:
533c7a4214aSPierre Jolivet     if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) {
5349566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(C,A->rmap->n,B->cmap->n,A->rmap->N,B->cmap->N));
535c7a4214aSPierre Jolivet     }
5368b8ddd11SPierre Jolivet     break;
5378b8ddd11SPierre Jolivet   case MATPRODUCT_AtB:
5388b8ddd11SPierre Jolivet     if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) {
5399566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(C,A->cmap->n,B->cmap->n,A->cmap->N,B->cmap->N));
5408b8ddd11SPierre Jolivet     }
5418b8ddd11SPierre Jolivet     break;
5428b8ddd11SPierre Jolivet   default:
54398921bdaSJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"ProductType %s is not supported",MatProductTypes[product->type]);
5448b8ddd11SPierre Jolivet   }
5459566063dSJacob Faibussowitsch   PetscCall(MatSetType(C,MATDENSE));
5469566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
5479566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE));
5489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
5499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
550c7a4214aSPierre Jolivet   C->ops->productsymbolic = NULL;
551c7a4214aSPierre Jolivet   C->ops->productnumeric = MatProductNumeric_Htool;
552c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
553c7a4214aSPierre Jolivet }
554c7a4214aSPierre Jolivet 
555c7a4214aSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Htool(Mat C)
556c7a4214aSPierre Jolivet {
557c7a4214aSPierre Jolivet   PetscFunctionBegin;
558c7a4214aSPierre Jolivet   MatCheckProduct(C,1);
5598b8ddd11SPierre Jolivet   if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool;
560c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
561c7a4214aSPierre Jolivet }
562c7a4214aSPierre Jolivet 
56398e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix)
564c7a4214aSPierre Jolivet {
565c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
566c7a4214aSPierre Jolivet 
567c7a4214aSPierre Jolivet   PetscFunctionBegin;
568c7a4214aSPierre Jolivet   *hmatrix = a->hmatrix;
569c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
570c7a4214aSPierre Jolivet }
571c7a4214aSPierre Jolivet 
572c7a4214aSPierre Jolivet /*@C
573c7a4214aSPierre Jolivet      MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a MATHTOOL.
574c7a4214aSPierre Jolivet 
575c7a4214aSPierre Jolivet    Input Parameter:
576c7a4214aSPierre Jolivet .     A - hierarchical matrix
577c7a4214aSPierre Jolivet 
578c7a4214aSPierre Jolivet    Output Parameter:
579c7a4214aSPierre Jolivet .     hmatrix - opaque pointer to a Htool virtual matrix
580c7a4214aSPierre Jolivet 
581c7a4214aSPierre Jolivet    Level: advanced
582c7a4214aSPierre Jolivet 
583c7a4214aSPierre Jolivet .seealso:  MATHTOOL
584c7a4214aSPierre Jolivet @*/
58598e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix)
586c7a4214aSPierre Jolivet {
587c7a4214aSPierre Jolivet   PetscFunctionBegin;
588c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
589c7a4214aSPierre Jolivet   PetscValidPointer(hmatrix,2);
590*cac4c232SBarry Smith   PetscTryMethod(A,"MatHtoolGetHierarchicalMat_C",(Mat,const htool::VirtualHMatrix<PetscScalar>**),(A,hmatrix));
591c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
592c7a4214aSPierre Jolivet }
593c7a4214aSPierre Jolivet 
594c7a4214aSPierre Jolivet static PetscErrorCode MatHtoolSetKernel_Htool(Mat A,MatHtoolKernel kernel,void *kernelctx)
595c7a4214aSPierre Jolivet {
596c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
597c7a4214aSPierre Jolivet 
598c7a4214aSPierre Jolivet   PetscFunctionBegin;
599c7a4214aSPierre Jolivet   a->kernel    = kernel;
600c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
601c7a4214aSPierre Jolivet   delete a->wrapper;
602c7a4214aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx);
603c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
604c7a4214aSPierre Jolivet }
605c7a4214aSPierre Jolivet 
606c7a4214aSPierre Jolivet /*@C
60798e73e17SPierre Jolivet      MatHtoolSetKernel - Sets the kernel and context used for the assembly of a MATHTOOL.
608c7a4214aSPierre Jolivet 
609c7a4214aSPierre Jolivet    Input Parameters:
610c7a4214aSPierre Jolivet +     A - hierarchical matrix
611c7a4214aSPierre Jolivet .     kernel - computational kernel (or NULL)
612cab00e0dSPierre Jolivet -     kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
613c7a4214aSPierre Jolivet 
614c7a4214aSPierre Jolivet    Level: advanced
615c7a4214aSPierre Jolivet 
616c7a4214aSPierre Jolivet .seealso:  MATHTOOL, MatCreateHtoolFromKernel()
617c7a4214aSPierre Jolivet @*/
618c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat A,MatHtoolKernel kernel,void *kernelctx)
619c7a4214aSPierre Jolivet {
620c7a4214aSPierre Jolivet   PetscFunctionBegin;
621c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
622c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel,2);
623c7a4214aSPierre Jolivet   if (!kernel)    PetscValidPointer(kernelctx,3);
624*cac4c232SBarry Smith   PetscTryMethod(A,"MatHtoolSetKernel_C",(Mat,MatHtoolKernel,void*),(A,kernel,kernelctx));
625c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
626c7a4214aSPierre Jolivet }
627c7a4214aSPierre Jolivet 
62898e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A,IS* is)
62998e73e17SPierre Jolivet {
63098e73e17SPierre Jolivet   Mat_Htool             *a = (Mat_Htool*)A->data;
63198e73e17SPierre Jolivet   std::vector<PetscInt> source;
63298e73e17SPierre Jolivet 
63398e73e17SPierre Jolivet   PetscFunctionBegin;
634a9918087SPierre Jolivet   source = a->hmatrix->get_source_cluster()->get_local_perm();
6359566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A),source.size(),source.data(),PETSC_COPY_VALUES,is));
6369566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
63798e73e17SPierre Jolivet   PetscFunctionReturn(0);
63898e73e17SPierre Jolivet }
63998e73e17SPierre Jolivet 
64098e73e17SPierre Jolivet /*@C
64198e73e17SPierre Jolivet      MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster.
64298e73e17SPierre Jolivet 
64398e73e17SPierre Jolivet    Input Parameter:
64498e73e17SPierre Jolivet .     A - hierarchical matrix
64598e73e17SPierre Jolivet 
64698e73e17SPierre Jolivet    Output Parameter:
64798e73e17SPierre Jolivet .     is - permutation
64898e73e17SPierre Jolivet 
64998e73e17SPierre Jolivet    Level: advanced
65098e73e17SPierre Jolivet 
65198e73e17SPierre Jolivet .seealso:  MATHTOOL, MatHtoolGetPermutationTarget(), MatHtoolUsePermutation()
65298e73e17SPierre Jolivet @*/
65398e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat A,IS* is)
65498e73e17SPierre Jolivet {
65598e73e17SPierre Jolivet   PetscFunctionBegin;
65698e73e17SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
65798e73e17SPierre Jolivet   if (!is) PetscValidPointer(is,2);
658*cac4c232SBarry Smith   PetscTryMethod(A,"MatHtoolGetPermutationSource_C",(Mat,IS*),(A,is));
65998e73e17SPierre Jolivet   PetscFunctionReturn(0);
66098e73e17SPierre Jolivet }
66198e73e17SPierre Jolivet 
66298e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A,IS* is)
66398e73e17SPierre Jolivet {
66498e73e17SPierre Jolivet   Mat_Htool             *a = (Mat_Htool*)A->data;
66598e73e17SPierre Jolivet   std::vector<PetscInt> target;
66698e73e17SPierre Jolivet 
66798e73e17SPierre Jolivet   PetscFunctionBegin;
668a9918087SPierre Jolivet   target = a->hmatrix->get_target_cluster()->get_local_perm();
6699566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A),target.size(),target.data(),PETSC_COPY_VALUES,is));
6709566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
67198e73e17SPierre Jolivet   PetscFunctionReturn(0);
67298e73e17SPierre Jolivet }
67398e73e17SPierre Jolivet 
67498e73e17SPierre Jolivet /*@C
67598e73e17SPierre Jolivet      MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster.
67698e73e17SPierre Jolivet 
67798e73e17SPierre Jolivet    Input Parameter:
67898e73e17SPierre Jolivet .     A - hierarchical matrix
67998e73e17SPierre Jolivet 
68098e73e17SPierre Jolivet    Output Parameter:
68198e73e17SPierre Jolivet .     is - permutation
68298e73e17SPierre Jolivet 
68398e73e17SPierre Jolivet    Level: advanced
68498e73e17SPierre Jolivet 
68598e73e17SPierre Jolivet .seealso:  MATHTOOL, MatHtoolGetPermutationSource(), MatHtoolUsePermutation()
68698e73e17SPierre Jolivet @*/
68798e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat A,IS* is)
68898e73e17SPierre Jolivet {
68998e73e17SPierre Jolivet   PetscFunctionBegin;
69098e73e17SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
69198e73e17SPierre Jolivet   if (!is) PetscValidPointer(is,2);
692*cac4c232SBarry Smith   PetscTryMethod(A,"MatHtoolGetPermutationTarget_C",(Mat,IS*),(A,is));
69398e73e17SPierre Jolivet   PetscFunctionReturn(0);
69498e73e17SPierre Jolivet }
69598e73e17SPierre Jolivet 
69698e73e17SPierre Jolivet static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A,PetscBool use)
69798e73e17SPierre Jolivet {
69898e73e17SPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
69998e73e17SPierre Jolivet 
70098e73e17SPierre Jolivet   PetscFunctionBegin;
70198e73e17SPierre Jolivet   a->hmatrix->set_use_permutation(use);
70298e73e17SPierre Jolivet   PetscFunctionReturn(0);
70398e73e17SPierre Jolivet }
70498e73e17SPierre Jolivet 
70598e73e17SPierre Jolivet /*@C
70698e73e17SPierre Jolivet      MatHtoolUsePermutation - Sets whether MATHTOOL should permute input (resp. output) vectors following its internal source (resp. target) permutation.
70798e73e17SPierre Jolivet 
70898e73e17SPierre Jolivet    Input Parameters:
70998e73e17SPierre Jolivet +     A - hierarchical matrix
71098e73e17SPierre Jolivet -     use - Boolean value
71198e73e17SPierre Jolivet 
71298e73e17SPierre Jolivet    Level: advanced
71398e73e17SPierre Jolivet 
71498e73e17SPierre Jolivet .seealso:  MATHTOOL, MatHtoolGetPermutationSource(), MatHtoolGetPermutationTarget()
71598e73e17SPierre Jolivet @*/
71698e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat A,PetscBool use)
71798e73e17SPierre Jolivet {
71898e73e17SPierre Jolivet   PetscFunctionBegin;
71998e73e17SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
72098e73e17SPierre Jolivet   PetscValidLogicalCollectiveBool(A,use,2);
721*cac4c232SBarry Smith   PetscTryMethod(A,"MatHtoolUsePermutation_C",(Mat,PetscBool),(A,use));
72298e73e17SPierre Jolivet   PetscFunctionReturn(0);
72398e73e17SPierre Jolivet }
72498e73e17SPierre Jolivet 
725c7a4214aSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B)
726c7a4214aSPierre Jolivet {
727c7a4214aSPierre Jolivet   Mat            C;
728c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
729c7a4214aSPierre Jolivet   PetscInt       lda;
730c7a4214aSPierre Jolivet   PetscScalar    *array;
731c7a4214aSPierre Jolivet 
732c7a4214aSPierre Jolivet   PetscFunctionBegin;
733c7a4214aSPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
734c7a4214aSPierre Jolivet     C = *B;
7355f80ce2aSJacob Faibussowitsch     PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible dimensions");
7369566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(C,&lda));
7375f80ce2aSJacob Faibussowitsch     PetscCheck(lda == C->rmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")",lda,C->rmap->n);
738c7a4214aSPierre Jolivet   } else {
7399566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C));
7409566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N));
7419566063dSJacob Faibussowitsch     PetscCall(MatSetType(C,MATDENSE));
7429566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
7439566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE));
744c7a4214aSPierre Jolivet   }
7459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
7469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
7479566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C,&array));
748c7a4214aSPierre Jolivet   a->hmatrix->copy_local_dense_perm(array);
7499566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C,&array));
7509566063dSJacob Faibussowitsch   PetscCall(MatScale(C,a->s));
751c7a4214aSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
7529566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&C));
753c7a4214aSPierre Jolivet   } else *B = C;
754c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
755c7a4214aSPierre Jolivet }
756c7a4214aSPierre Jolivet 
75798e73e17SPierre Jolivet static PetscErrorCode GenEntriesTranspose(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *rows,const PetscInt *cols,PetscScalar *ptr,void *ctx)
758c7a4214aSPierre Jolivet {
759c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose*)ctx;
76098e73e17SPierre Jolivet   PetscScalar             *tmp;
76198e73e17SPierre Jolivet 
76298e73e17SPierre Jolivet   PetscFunctionBegin;
76398e73e17SPierre Jolivet   generator->kernel(sdim,N,M,cols,rows,ptr,generator->kernelctx);
7649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M*N,&tmp));
7659566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tmp,ptr,M*N));
76698e73e17SPierre Jolivet   for (PetscInt i=0; i<M; ++i) {
76798e73e17SPierre Jolivet     for (PetscInt j=0; j<N; ++j) ptr[i+j*M] = tmp[j+i*N];
76898e73e17SPierre Jolivet   }
7699566063dSJacob Faibussowitsch   PetscCall(PetscFree(tmp));
77098e73e17SPierre Jolivet   PetscFunctionReturn(0);
771c7a4214aSPierre Jolivet }
772c7a4214aSPierre Jolivet 
773c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */
774c7a4214aSPierre Jolivet static PetscErrorCode MatTranspose_Htool(Mat A,MatReuse reuse,Mat *B)
775c7a4214aSPierre Jolivet {
776c7a4214aSPierre Jolivet   Mat                     C;
777c7a4214aSPierre Jolivet   Mat_Htool               *a = (Mat_Htool*)A->data,*c;
778c7a4214aSPierre Jolivet   PetscInt                M = A->rmap->N,N = A->cmap->N,m = A->rmap->n,n = A->cmap->n;
779c7a4214aSPierre Jolivet   PetscContainer          container;
780c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
781c7a4214aSPierre Jolivet 
782c7a4214aSPierre Jolivet   PetscFunctionBegin;
7835f80ce2aSJacob Faibussowitsch   PetscCheck(reuse != MAT_INPLACE_MATRIX,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatTranspose() with MAT_INPLACE_MATRIX not supported");
784c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
7859566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C));
7869566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C,n,m,N,M));
7879566063dSJacob Faibussowitsch     PetscCall(MatSetType(C,((PetscObject)A)->type_name));
7889566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
7899566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)C),&container));
7909566063dSJacob Faibussowitsch     PetscCall(PetscNew(&kernelt));
7919566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(container,kernelt));
7929566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)C,"KernelTranspose",(PetscObject)container));
793c7a4214aSPierre Jolivet   } else {
794c7a4214aSPierre Jolivet     C = *B;
7959566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)C,"KernelTranspose",(PetscObject*)&container));
7965f80ce2aSJacob Faibussowitsch     PetscCheck(container,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatTranspose() with MAT_INITIAL_MATRIX first");
7979566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container,(void**)&kernelt));
798c7a4214aSPierre Jolivet   }
799c7a4214aSPierre Jolivet   c                  = (Mat_Htool*)C->data;
800c7a4214aSPierre Jolivet   c->dim             = a->dim;
801c7a4214aSPierre Jolivet   c->s               = a->s;
80298e73e17SPierre Jolivet   c->kernel          = GenEntriesTranspose;
803c7a4214aSPierre Jolivet   if (kernelt->A != A) {
8049566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
805c7a4214aSPierre Jolivet     kernelt->A       = A;
8069566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
807c7a4214aSPierre Jolivet   }
808c7a4214aSPierre Jolivet   kernelt->kernel    = a->kernel;
809c7a4214aSPierre Jolivet   kernelt->kernelctx = a->kernelctx;
810c7a4214aSPierre Jolivet   c->kernelctx       = kernelt;
811c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
8129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(N*c->dim,&c->gcoords_target));
8139566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(c->gcoords_target,a->gcoords_source,N*c->dim));
814c7a4214aSPierre Jolivet     if (a->gcoords_target != a->gcoords_source) {
8159566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(M*c->dim,&c->gcoords_source));
8169566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->gcoords_source,a->gcoords_target,M*c->dim));
817c7a4214aSPierre Jolivet     } else c->gcoords_source = c->gcoords_target;
8189566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(M,&c->work_source,N,&c->work_target));
819c7a4214aSPierre Jolivet   }
8209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
8219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
822c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) *B = C;
823c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
824c7a4214aSPierre Jolivet }
825c7a4214aSPierre Jolivet 
826c7a4214aSPierre Jolivet /*@C
827c7a4214aSPierre Jolivet      MatCreateHtoolFromKernel - Creates a MATHTOOL from a user-supplied kernel.
828c7a4214aSPierre Jolivet 
829c7a4214aSPierre Jolivet    Input Parameters:
830c7a4214aSPierre Jolivet +     comm - MPI communicator
831c7a4214aSPierre Jolivet .     m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
832c7a4214aSPierre Jolivet .     n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
833c7a4214aSPierre Jolivet .     M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
834c7a4214aSPierre Jolivet .     N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
835c7a4214aSPierre Jolivet .     spacedim - dimension of the space coordinates
836c7a4214aSPierre Jolivet .     coords_target - coordinates of the target
837c7a4214aSPierre Jolivet .     coords_source - coordinates of the source
838c7a4214aSPierre Jolivet .     kernel - computational kernel (or NULL)
839cab00e0dSPierre Jolivet -     kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
840c7a4214aSPierre Jolivet 
841c7a4214aSPierre Jolivet    Output Parameter:
842c7a4214aSPierre Jolivet .     B - matrix
843c7a4214aSPierre Jolivet 
844c7a4214aSPierre Jolivet    Options Database Keys:
845c7a4214aSPierre Jolivet +     -mat_htool_min_cluster_size <PetscInt> - minimal leaf size in cluster tree
846c7a4214aSPierre Jolivet .     -mat_htool_max_block_size <PetscInt> - maximal number of coefficients in a dense block
847c7a4214aSPierre Jolivet .     -mat_htool_epsilon <PetscReal> - relative error in Frobenius norm when approximating a block
848c7a4214aSPierre Jolivet .     -mat_htool_eta <PetscReal> - admissibility condition tolerance
849c7a4214aSPierre Jolivet .     -mat_htool_min_target_depth <PetscInt> - minimal cluster tree depth associated with the rows
850c7a4214aSPierre Jolivet .     -mat_htool_min_source_depth <PetscInt> - minimal cluster tree depth associated with the columns
85198e73e17SPierre Jolivet .     -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression
85298e73e17SPierre Jolivet -     -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering
853c7a4214aSPierre Jolivet 
854c7a4214aSPierre Jolivet    Level: intermediate
855c7a4214aSPierre Jolivet 
85653022affSStefano Zampini .seealso:  MatCreate(), MATHTOOL, PCSetCoordinates(), MatHtoolSetKernel(), MatHtoolCompressorType, MATH2OPUS, MatCreateH2OpusFromKernel()
857c7a4214aSPierre Jolivet @*/
858c7a4214aSPierre 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)
859c7a4214aSPierre Jolivet {
860c7a4214aSPierre Jolivet   Mat            A;
861c7a4214aSPierre Jolivet   Mat_Htool      *a;
862c7a4214aSPierre Jolivet 
863c7a4214aSPierre Jolivet   PetscFunctionBegin;
8649566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,&A));
865c7a4214aSPierre Jolivet   PetscValidLogicalCollectiveInt(A,spacedim,6);
866c7a4214aSPierre Jolivet   PetscValidRealPointer(coords_target,7);
867c7a4214aSPierre Jolivet   PetscValidRealPointer(coords_source,8);
868c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel,9);
869c7a4214aSPierre Jolivet   if (!kernel)    PetscValidPointer(kernelctx,10);
8709566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,m,n,M,N));
8719566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATHTOOL));
8729566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
873c7a4214aSPierre Jolivet   a            = (Mat_Htool*)A->data;
874c7a4214aSPierre Jolivet   a->dim       = spacedim;
875c7a4214aSPierre Jolivet   a->s         = 1.0;
876c7a4214aSPierre Jolivet   a->kernel    = kernel;
877c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
8789566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(A->rmap->N*spacedim,&a->gcoords_target));
8799566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->gcoords_target+A->rmap->rstart*spacedim,coords_target,A->rmap->n*spacedim));
8809566063dSJacob Faibussowitsch   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE,a->gcoords_target,A->rmap->N*spacedim,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)A))); /* global target coordinates */
881c7a4214aSPierre Jolivet   if (coords_target != coords_source) {
8829566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(A->cmap->N*spacedim,&a->gcoords_source));
8839566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->gcoords_source+A->cmap->rstart*spacedim,coords_source,A->cmap->n*spacedim));
8849566063dSJacob Faibussowitsch     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE,a->gcoords_source,A->cmap->N*spacedim,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)A))); /* global source coordinates */
885c7a4214aSPierre Jolivet   } else a->gcoords_source = a->gcoords_target;
8869566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(A->cmap->N,&a->work_source,A->rmap->N,&a->work_target));
887c7a4214aSPierre Jolivet   *B = A;
888c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
889c7a4214aSPierre Jolivet }
890c7a4214aSPierre Jolivet 
891c7a4214aSPierre Jolivet /*MC
892c7a4214aSPierre Jolivet      MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package.
893c7a4214aSPierre Jolivet 
894c7a4214aSPierre Jolivet   Use ./configure --download-htool to install PETSc to use Htool.
895c7a4214aSPierre Jolivet 
896c7a4214aSPierre Jolivet    Options Database Keys:
897e85a1b63SPatrick Sanan .     -mat_type htool - matrix type to "htool" during a call to MatSetFromOptions()
898c7a4214aSPierre Jolivet 
899c7a4214aSPierre Jolivet    Level: beginner
900c7a4214aSPierre Jolivet 
90153022affSStefano Zampini .seealso: MATH2OPUS, MATDENSE, MatCreateHtoolFromKernel(), MatHtoolSetKernel()
902c7a4214aSPierre Jolivet M*/
903c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A)
904c7a4214aSPierre Jolivet {
905c7a4214aSPierre Jolivet   Mat_Htool *a;
906c7a4214aSPierre Jolivet 
907c7a4214aSPierre Jolivet   PetscFunctionBegin;
9089566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(A,&a));
909c7a4214aSPierre Jolivet   A->data = (void*)a;
9109566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATHTOOL));
9119566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(A->ops,sizeof(struct _MatOps)));
912c7a4214aSPierre Jolivet   A->ops->getdiagonal       = MatGetDiagonal_Htool;
913c7a4214aSPierre Jolivet   A->ops->getdiagonalblock  = MatGetDiagonalBlock_Htool;
914c7a4214aSPierre Jolivet   A->ops->mult              = MatMult_Htool;
915c7a4214aSPierre Jolivet   A->ops->multadd           = MatMultAdd_Htool;
9168b8ddd11SPierre Jolivet   A->ops->multtranspose     = MatMultTranspose_Htool;
9178b8ddd11SPierre Jolivet   if (!PetscDefined(USE_COMPLEX)) A->ops->multhermitiantranspose = MatMultTranspose_Htool;
918c7a4214aSPierre Jolivet   A->ops->increaseoverlap   = MatIncreaseOverlap_Htool;
919c7a4214aSPierre Jolivet   A->ops->createsubmatrices = MatCreateSubMatrices_Htool;
920c7a4214aSPierre Jolivet   A->ops->transpose         = MatTranspose_Htool;
921c7a4214aSPierre Jolivet   A->ops->destroy           = MatDestroy_Htool;
922c7a4214aSPierre Jolivet   A->ops->view              = MatView_Htool;
923c7a4214aSPierre Jolivet   A->ops->setfromoptions    = MatSetFromOptions_Htool;
924c7a4214aSPierre Jolivet   A->ops->scale             = MatScale_Htool;
925c7a4214aSPierre Jolivet   A->ops->getrow            = MatGetRow_Htool;
926c7a4214aSPierre Jolivet   A->ops->restorerow        = MatRestoreRow_Htool;
927c7a4214aSPierre Jolivet   A->ops->assemblyend       = MatAssemblyEnd_Htool;
928c7a4214aSPierre Jolivet   a->dim                    = 0;
929c7a4214aSPierre Jolivet   a->gcoords_target         = NULL;
930c7a4214aSPierre Jolivet   a->gcoords_source         = NULL;
931c7a4214aSPierre Jolivet   a->s                      = 1.0;
932c7a4214aSPierre Jolivet   a->bs[0]                  = 10;
933c7a4214aSPierre Jolivet   a->bs[1]                  = 1000000;
934c7a4214aSPierre Jolivet   a->epsilon                = PetscSqrtReal(PETSC_SMALL);
935c7a4214aSPierre Jolivet   a->eta                    = 10.0;
936c7a4214aSPierre Jolivet   a->depth[0]               = 0;
937c7a4214aSPierre Jolivet   a->depth[1]               = 0;
938c7a4214aSPierre Jolivet   a->compressor             = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA;
9399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",MatProductSetFromOptions_Htool));
9409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",MatProductSetFromOptions_Htool));
9419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",MatConvert_Htool_Dense));
9429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",MatConvert_Htool_Dense));
9439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",MatHtoolGetHierarchicalMat_Htool));
9449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",MatHtoolSetKernel_Htool));
9459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",MatHtoolGetPermutationSource_Htool));
9469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",MatHtoolGetPermutationTarget_Htool));
9479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",MatHtoolUsePermutation_Htool));
948c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
949c7a4214aSPierre Jolivet }
950