xref: /petsc/src/mat/impls/htool/htool.cxx (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
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   PetscErrorCode ierr;
26c7a4214aSPierre Jolivet 
27c7a4214aSPierre Jolivet   PetscFunctionBegin;
28c7a4214aSPierre Jolivet   ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr);
29c7a4214aSPierre Jolivet   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only congruent layouts supported");
30c7a4214aSPierre Jolivet   ierr = VecGetArrayWrite(v,&x);CHKERRQ(ierr);
31c7a4214aSPierre Jolivet   a->hmatrix->copy_local_diagonal(x);
32c7a4214aSPierre Jolivet   ierr = VecRestoreArrayWrite(v,&x);CHKERRQ(ierr);
33c7a4214aSPierre Jolivet   ierr = VecScale(v,a->s);CHKERRQ(ierr);
34c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
35c7a4214aSPierre Jolivet }
36c7a4214aSPierre Jolivet 
37c7a4214aSPierre Jolivet static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A,Mat *b)
38c7a4214aSPierre Jolivet {
39c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
40c7a4214aSPierre Jolivet   Mat            B;
41c7a4214aSPierre Jolivet   PetscScalar    *ptr;
42c7a4214aSPierre Jolivet   PetscBool      flg;
43c7a4214aSPierre Jolivet   PetscErrorCode ierr;
44c7a4214aSPierre Jolivet 
45c7a4214aSPierre Jolivet   PetscFunctionBegin;
46c7a4214aSPierre Jolivet   ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr);
47c7a4214aSPierre Jolivet   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only congruent layouts supported");
48c7a4214aSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B);CHKERRQ(ierr); /* same logic as in MatGetDiagonalBlock_MPIDense() */
49c7a4214aSPierre Jolivet   if (!B) {
50c7a4214aSPierre Jolivet     ierr = MatCreateDense(PETSC_COMM_SELF,A->rmap->n,A->rmap->n,A->rmap->n,A->rmap->n,NULL,&B);CHKERRQ(ierr);
51c7a4214aSPierre Jolivet     ierr = MatDenseGetArrayWrite(B,&ptr);CHKERRQ(ierr);
52c7a4214aSPierre Jolivet     a->hmatrix->copy_local_diagonal_block(ptr);
53c7a4214aSPierre Jolivet     ierr = MatDenseRestoreArrayWrite(B,&ptr);CHKERRQ(ierr);
54c7a4214aSPierre Jolivet     ierr = MatPropagateSymmetryOptions(A,B);CHKERRQ(ierr);
55c7a4214aSPierre Jolivet     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56c7a4214aSPierre Jolivet     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
57c7a4214aSPierre Jolivet     ierr = MatScale(B,a->s);CHKERRQ(ierr);
58c7a4214aSPierre Jolivet     ierr = PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B);CHKERRQ(ierr);
59c7a4214aSPierre Jolivet     *b   = B;
60c7a4214aSPierre Jolivet     ierr = MatDestroy(&B);CHKERRQ(ierr);
61c7a4214aSPierre Jolivet   } else *b = B;
62c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
63c7a4214aSPierre Jolivet }
64c7a4214aSPierre Jolivet 
65c7a4214aSPierre Jolivet static PetscErrorCode MatMult_Htool(Mat A,Vec x,Vec y)
66c7a4214aSPierre Jolivet {
67c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
68c7a4214aSPierre Jolivet   const PetscScalar *in;
69c7a4214aSPierre Jolivet   PetscScalar       *out;
70c7a4214aSPierre Jolivet   PetscErrorCode    ierr;
71c7a4214aSPierre Jolivet 
72c7a4214aSPierre Jolivet   PetscFunctionBegin;
73c7a4214aSPierre Jolivet   ierr = VecGetArrayRead(x,&in);CHKERRQ(ierr);
74c7a4214aSPierre Jolivet   ierr = VecGetArrayWrite(y,&out);CHKERRQ(ierr);
75c7a4214aSPierre Jolivet   a->hmatrix->mvprod_local_to_local(in,out);
76c7a4214aSPierre Jolivet   ierr = VecRestoreArrayRead(x,&in);CHKERRQ(ierr);
77c7a4214aSPierre Jolivet   ierr = VecRestoreArrayWrite(y,&out);CHKERRQ(ierr);
78c7a4214aSPierre Jolivet   ierr = VecScale(y,a->s);CHKERRQ(ierr);
79c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
80c7a4214aSPierre Jolivet }
81c7a4214aSPierre Jolivet 
82c7a4214aSPierre Jolivet /* naive implementation of MatMultAdd() needed for FEM-BEM coupling via MATNEST */
83c7a4214aSPierre Jolivet static PetscErrorCode MatMultAdd_Htool(Mat A,Vec v1,Vec v2,Vec v3)
84c7a4214aSPierre Jolivet {
85c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
86c7a4214aSPierre Jolivet   Vec               tmp;
87c7a4214aSPierre Jolivet   const PetscScalar scale = a->s;
88c7a4214aSPierre Jolivet   PetscErrorCode    ierr;
89c7a4214aSPierre Jolivet 
90c7a4214aSPierre Jolivet   PetscFunctionBegin;
91c7a4214aSPierre Jolivet   ierr = VecDuplicate(v2,&tmp);CHKERRQ(ierr);
92c7a4214aSPierre Jolivet   ierr = VecCopy(v2,v3);CHKERRQ(ierr); /* no-op in MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]) since VecCopy() checks for x == y */
93c7a4214aSPierre Jolivet   a->s = 1.0; /* set s to 1.0 since VecAXPY() may be used to scale the MatMult() output Vec */
94c7a4214aSPierre Jolivet   ierr = MatMult_Htool(A,v1,tmp);CHKERRQ(ierr);
95c7a4214aSPierre Jolivet   ierr = VecAXPY(v3,scale,tmp);CHKERRQ(ierr);
96c7a4214aSPierre Jolivet   ierr = VecDestroy(&tmp);CHKERRQ(ierr);
97c7a4214aSPierre Jolivet   a->s = scale; /* set s back to its original value */
98c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
99c7a4214aSPierre Jolivet }
100c7a4214aSPierre Jolivet 
1018b8ddd11SPierre Jolivet static PetscErrorCode MatMultTranspose_Htool(Mat A,Vec x,Vec y)
1028b8ddd11SPierre Jolivet {
1038b8ddd11SPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
1048b8ddd11SPierre Jolivet   const PetscScalar *in;
1058b8ddd11SPierre Jolivet   PetscScalar       *out;
1068b8ddd11SPierre Jolivet   PetscErrorCode    ierr;
1078b8ddd11SPierre Jolivet 
1088b8ddd11SPierre Jolivet   PetscFunctionBegin;
1098b8ddd11SPierre Jolivet   ierr = VecGetArrayRead(x,&in);CHKERRQ(ierr);
1108b8ddd11SPierre Jolivet   ierr = VecGetArrayWrite(y,&out);CHKERRQ(ierr);
1118b8ddd11SPierre Jolivet   a->hmatrix->mvprod_transp_local_to_local(in,out);
1128b8ddd11SPierre Jolivet   ierr = VecRestoreArrayRead(x,&in);CHKERRQ(ierr);
1138b8ddd11SPierre Jolivet   ierr = VecRestoreArrayWrite(y,&out);CHKERRQ(ierr);
1148b8ddd11SPierre Jolivet   ierr = VecScale(y,a->s);CHKERRQ(ierr);
1158b8ddd11SPierre Jolivet   PetscFunctionReturn(0);
1168b8ddd11SPierre Jolivet }
1178b8ddd11SPierre Jolivet 
118c7a4214aSPierre Jolivet static PetscErrorCode MatIncreaseOverlap_Htool(Mat A,PetscInt is_max,IS is[],PetscInt ov)
119c7a4214aSPierre Jolivet {
120c7a4214aSPierre Jolivet   std::set<PetscInt> set;
121c7a4214aSPierre Jolivet   const PetscInt     *idx;
1228f308287SPierre Jolivet   PetscInt           *oidx,size,bs[2];
123c7a4214aSPierre Jolivet   PetscMPIInt        csize;
124c7a4214aSPierre Jolivet   PetscErrorCode     ierr;
125c7a4214aSPierre Jolivet 
126c7a4214aSPierre Jolivet   PetscFunctionBegin;
1278f308287SPierre Jolivet   ierr = MatGetBlockSizes(A,bs,bs+1);CHKERRQ(ierr);
1288f308287SPierre Jolivet   if (bs[0] != bs[1]) bs[0] = 1;
129c7a4214aSPierre Jolivet   for (PetscInt i=0; i<is_max; ++i) {
130c7a4214aSPierre Jolivet     /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */
131c7a4214aSPierre Jolivet     /* needed to avoid subdomain matrices to replicate A since it is dense                           */
132c7a4214aSPierre Jolivet     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)is[i]),&csize);CHKERRMPI(ierr);
133c7a4214aSPierre Jolivet     if (csize != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported parallel IS");
134c7a4214aSPierre Jolivet     ierr = ISGetSize(is[i],&size);CHKERRQ(ierr);
135c7a4214aSPierre Jolivet     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
136c7a4214aSPierre Jolivet     for (PetscInt j=0; j<size; ++j) {
137c7a4214aSPierre Jolivet       set.insert(idx[j]);
138c7a4214aSPierre Jolivet       for (PetscInt k=1; k<=ov; ++k) {               /* for each layer of overlap      */
139c7a4214aSPierre Jolivet         if (idx[j] - k >= 0) set.insert(idx[j] - k); /* do not insert negative indices */
140c7a4214aSPierre 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 */
141c7a4214aSPierre Jolivet       }
142c7a4214aSPierre Jolivet     }
143c7a4214aSPierre Jolivet     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
144c7a4214aSPierre Jolivet     ierr = ISDestroy(is+i);CHKERRQ(ierr);
1458f308287SPierre Jolivet     if (bs[0] > 1) {
1468f308287SPierre Jolivet       for (std::set<PetscInt>::iterator it=set.cbegin(); it!=set.cend(); it++) {
1478f308287SPierre Jolivet         std::vector<PetscInt> block(bs[0]);
1488f308287SPierre Jolivet         std::iota(block.begin(),block.end(),(*it/bs[0])*bs[0]);
1498f308287SPierre Jolivet         set.insert(block.cbegin(),block.cend());
1508f308287SPierre Jolivet       }
1518f308287SPierre Jolivet     }
152c7a4214aSPierre Jolivet     size = set.size(); /* size with overlap */
153c7a4214aSPierre Jolivet     ierr = PetscMalloc1(size,&oidx);CHKERRQ(ierr);
154c7a4214aSPierre Jolivet     for (const PetscInt j : set) *oidx++ = j;
155c7a4214aSPierre Jolivet     oidx -= size;
156c7a4214aSPierre Jolivet     ierr = ISCreateGeneral(PETSC_COMM_SELF,size,oidx,PETSC_OWN_POINTER,is+i);CHKERRQ(ierr);
157c7a4214aSPierre Jolivet   }
158c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
159c7a4214aSPierre Jolivet }
160c7a4214aSPierre Jolivet 
161c7a4214aSPierre Jolivet static PetscErrorCode MatCreateSubMatrices_Htool(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
162c7a4214aSPierre Jolivet {
163c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)A->data;
164c7a4214aSPierre Jolivet   Mat               D,B,BT;
165c7a4214aSPierre Jolivet   const PetscScalar *copy;
166c7a4214aSPierre Jolivet   PetscScalar       *ptr;
167c7a4214aSPierre Jolivet   const PetscInt    *idxr,*idxc,*it;
168c7a4214aSPierre Jolivet   PetscInt          nrow,m,i;
169c7a4214aSPierre Jolivet   PetscBool         flg;
170c7a4214aSPierre Jolivet   PetscErrorCode    ierr;
171c7a4214aSPierre Jolivet 
172c7a4214aSPierre Jolivet   PetscFunctionBegin;
173c7a4214aSPierre Jolivet   if (scall != MAT_REUSE_MATRIX) {
174c7a4214aSPierre Jolivet     ierr = PetscCalloc1(n,submat);CHKERRQ(ierr);
175c7a4214aSPierre Jolivet   }
176c7a4214aSPierre Jolivet   for (i=0; i<n; ++i) {
177c7a4214aSPierre Jolivet     ierr = ISGetLocalSize(irow[i],&nrow);CHKERRQ(ierr);
178c7a4214aSPierre Jolivet     ierr = ISGetLocalSize(icol[i],&m);CHKERRQ(ierr);
179c7a4214aSPierre Jolivet     ierr = ISGetIndices(irow[i],&idxr);CHKERRQ(ierr);
180c7a4214aSPierre Jolivet     ierr = ISGetIndices(icol[i],&idxc);CHKERRQ(ierr);
181c7a4214aSPierre Jolivet     if (scall != MAT_REUSE_MATRIX) {
182c7a4214aSPierre Jolivet       ierr = MatCreateDense(PETSC_COMM_SELF,nrow,m,nrow,m,NULL,(*submat)+i);CHKERRQ(ierr);
183c7a4214aSPierre Jolivet     }
184c7a4214aSPierre Jolivet     ierr = MatDenseGetArrayWrite((*submat)[i],&ptr);CHKERRQ(ierr);
185c7a4214aSPierre Jolivet     if (irow[i] == icol[i]) { /* same row and column IS? */
186c7a4214aSPierre Jolivet       ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr);
187c7a4214aSPierre Jolivet       if (flg) {
188c7a4214aSPierre Jolivet         ierr = ISSorted(irow[i],&flg);CHKERRQ(ierr);
189c7a4214aSPierre Jolivet         if (flg) { /* sorted IS? */
190c7a4214aSPierre Jolivet           it = std::lower_bound(idxr,idxr+nrow,A->rmap->rstart);
191c7a4214aSPierre Jolivet           if (it != idxr+nrow && *it == A->rmap->rstart) { /* rmap->rstart in IS? */
192c7a4214aSPierre Jolivet             if (std::distance(idxr,it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */
193c7a4214aSPierre Jolivet               for (PetscInt j=0; j<A->rmap->n && flg; ++j) if (PetscUnlikely(it[j] != A->rmap->rstart+j)) flg = PETSC_FALSE;
194c7a4214aSPierre Jolivet               if (flg) { /* complete local diagonal block in IS? */
195c7a4214aSPierre Jolivet                 /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM
196c7a4214aSPierre Jolivet                  *      [   B   C   E   ]
197c7a4214aSPierre Jolivet                  *  A = [   B   D   E   ]
198c7a4214aSPierre Jolivet                  *      [   B   F   E   ]
199c7a4214aSPierre Jolivet                  */
200c7a4214aSPierre Jolivet                 m = std::distance(idxr,it); /* shift of the coefficient (0,0) of block D from above */
201c7a4214aSPierre Jolivet                 ierr = MatGetDiagonalBlock_Htool(A,&D);CHKERRQ(ierr);
202c7a4214aSPierre Jolivet                 ierr = MatDenseGetArrayRead(D,&copy);CHKERRQ(ierr);
203c7a4214aSPierre Jolivet                 for (PetscInt k=0; k<A->rmap->n; ++k) {
204c7a4214aSPierre Jolivet                   ierr = PetscArraycpy(ptr+(m+k)*nrow+m,copy+k*A->rmap->n,A->rmap->n);CHKERRQ(ierr); /* block D from above */
205c7a4214aSPierre Jolivet                 }
206c7a4214aSPierre Jolivet                 ierr = MatDenseRestoreArrayRead(D,&copy);CHKERRQ(ierr);
207c7a4214aSPierre Jolivet                 if (m) {
208c7a4214aSPierre Jolivet                   a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* vertical block B from above */
209c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
210c7a4214aSPierre Jolivet                   if (A->symmetric || A->hermitian) {
211c7a4214aSPierre Jolivet                     ierr = MatCreateDense(PETSC_COMM_SELF,A->rmap->n,m,A->rmap->n,m,ptr+m,&B);CHKERRQ(ierr);
212c7a4214aSPierre Jolivet                     ierr = MatDenseSetLDA(B,nrow);CHKERRQ(ierr);
213c7a4214aSPierre Jolivet                     ierr = MatCreateDense(PETSC_COMM_SELF,m,A->rmap->n,m,A->rmap->n,ptr+m*nrow,&BT);CHKERRQ(ierr);
214c7a4214aSPierre Jolivet                     ierr = MatDenseSetLDA(BT,nrow);CHKERRQ(ierr);
215c7a4214aSPierre Jolivet                     if (A->hermitian && PetscDefined(USE_COMPLEX)) {
216c7a4214aSPierre Jolivet                       ierr = MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr);
217c7a4214aSPierre Jolivet                     } else {
218c7a4214aSPierre Jolivet                       ierr = MatTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr);
219c7a4214aSPierre Jolivet                     }
220c7a4214aSPierre Jolivet                     ierr = MatDestroy(&B);CHKERRQ(ierr);
221c7a4214aSPierre Jolivet                     ierr = MatDestroy(&BT);CHKERRQ(ierr);
222c7a4214aSPierre Jolivet                   } else {
223c7a4214aSPierre Jolivet                     for (PetscInt k=0; k<A->rmap->n; ++k) { /* block C from above */
224c7a4214aSPierre Jolivet                       a->wrapper->copy_submatrix(m,1,idxr,idxc+m+k,ptr+(m+k)*nrow);
225c7a4214aSPierre Jolivet                     }
226c7a4214aSPierre Jolivet                   }
227c7a4214aSPierre Jolivet                 }
228c7a4214aSPierre Jolivet                 if (m+A->rmap->n != nrow) {
229c7a4214aSPierre 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 */
230c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
231c7a4214aSPierre Jolivet                   if (A->symmetric || A->hermitian) {
232c7a4214aSPierre Jolivet                     ierr = 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);CHKERRQ(ierr);
233c7a4214aSPierre Jolivet                     ierr = MatDenseSetLDA(B,nrow);CHKERRQ(ierr);
234c7a4214aSPierre Jolivet                     ierr = 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);CHKERRQ(ierr);
235c7a4214aSPierre Jolivet                     ierr = MatDenseSetLDA(BT,nrow);CHKERRQ(ierr);
236c7a4214aSPierre Jolivet                     if (A->hermitian && PetscDefined(USE_COMPLEX)) {
237c7a4214aSPierre Jolivet                       ierr = MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr);
238c7a4214aSPierre Jolivet                     } else {
239c7a4214aSPierre Jolivet                       ierr = MatTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr);
240c7a4214aSPierre Jolivet                     }
241c7a4214aSPierre Jolivet                     ierr = MatDestroy(&B);CHKERRQ(ierr);
242c7a4214aSPierre Jolivet                     ierr = MatDestroy(&BT);CHKERRQ(ierr);
243c7a4214aSPierre Jolivet                   } else {
244c7a4214aSPierre Jolivet                     for (PetscInt k=0; k<A->rmap->n; ++k) { /* block F from above */
245c7a4214aSPierre 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);
246c7a4214aSPierre Jolivet                     }
247c7a4214aSPierre Jolivet                   }
248c7a4214aSPierre Jolivet                 }
249c7a4214aSPierre Jolivet               } /* complete local diagonal block not in IS */
250c7a4214aSPierre Jolivet             } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */
251c7a4214aSPierre Jolivet           } else flg = PETSC_FALSE; /* rmap->rstart not in IS */
252c7a4214aSPierre Jolivet         } /* unsorted IS */
253c7a4214aSPierre Jolivet       }
254c7a4214aSPierre Jolivet     } else flg = PETSC_FALSE; /* different row and column IS */
255c7a4214aSPierre Jolivet     if (!flg) a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* reassemble everything */
256c7a4214aSPierre Jolivet     ierr = ISRestoreIndices(irow[i],&idxr);CHKERRQ(ierr);
257c7a4214aSPierre Jolivet     ierr = ISRestoreIndices(icol[i],&idxc);CHKERRQ(ierr);
258c7a4214aSPierre Jolivet     ierr = MatDenseRestoreArrayWrite((*submat)[i],&ptr);CHKERRQ(ierr);
259c7a4214aSPierre Jolivet     ierr = MatAssemblyBegin((*submat)[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
260c7a4214aSPierre Jolivet     ierr = MatAssemblyEnd((*submat)[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
261c7a4214aSPierre Jolivet     ierr = MatScale((*submat)[i],a->s);CHKERRQ(ierr);
262c7a4214aSPierre Jolivet   }
263c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
264c7a4214aSPierre Jolivet }
265c7a4214aSPierre Jolivet 
266c7a4214aSPierre Jolivet static PetscErrorCode MatDestroy_Htool(Mat A)
267c7a4214aSPierre Jolivet {
268c7a4214aSPierre Jolivet   Mat_Htool               *a = (Mat_Htool*)A->data;
269c7a4214aSPierre Jolivet   PetscContainer          container;
270c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
271c7a4214aSPierre Jolivet   PetscErrorCode          ierr;
272c7a4214aSPierre Jolivet 
273c7a4214aSPierre Jolivet   PetscFunctionBegin;
274c7a4214aSPierre Jolivet   ierr = PetscObjectChangeTypeName((PetscObject)A,NULL);CHKERRQ(ierr);
275c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",NULL);CHKERRQ(ierr);
276c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",NULL);CHKERRQ(ierr);
277c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",NULL);CHKERRQ(ierr);
278c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",NULL);CHKERRQ(ierr);
279c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",NULL);CHKERRQ(ierr);
280c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",NULL);CHKERRQ(ierr);
28198e73e17SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",NULL);CHKERRQ(ierr);
28298e73e17SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",NULL);CHKERRQ(ierr);
28398e73e17SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",NULL);CHKERRQ(ierr);
284c7a4214aSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)A,"KernelTranspose",(PetscObject*)&container);CHKERRQ(ierr);
285c7a4214aSPierre Jolivet   if (container) { /* created in MatTranspose_Htool() */
286c7a4214aSPierre Jolivet     ierr = PetscContainerGetPointer(container,(void**)&kernelt);CHKERRQ(ierr);
287c7a4214aSPierre Jolivet     ierr = MatDestroy(&kernelt->A);CHKERRQ(ierr);
288c7a4214aSPierre Jolivet     ierr = PetscFree(kernelt);CHKERRQ(ierr);
289c7a4214aSPierre Jolivet     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
290c7a4214aSPierre Jolivet     ierr = PetscObjectCompose((PetscObject)A,"KernelTranspose",NULL);CHKERRQ(ierr);
291c7a4214aSPierre Jolivet   }
292c7a4214aSPierre Jolivet   if (a->gcoords_source != a->gcoords_target) {
293c7a4214aSPierre Jolivet     ierr = PetscFree(a->gcoords_source);CHKERRQ(ierr);
294c7a4214aSPierre Jolivet   }
295c7a4214aSPierre Jolivet   ierr = PetscFree(a->gcoords_target);CHKERRQ(ierr);
296c7a4214aSPierre Jolivet   ierr = PetscFree2(a->work_source,a->work_target);CHKERRQ(ierr);
297c7a4214aSPierre Jolivet   delete a->wrapper;
298c7a4214aSPierre Jolivet   delete a->hmatrix;
299c7a4214aSPierre Jolivet   ierr = PetscFree(A->data);CHKERRQ(ierr);
300c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
301c7a4214aSPierre Jolivet }
302c7a4214aSPierre Jolivet 
303c7a4214aSPierre Jolivet static PetscErrorCode MatView_Htool(Mat A,PetscViewer pv)
304c7a4214aSPierre Jolivet {
305c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
306c7a4214aSPierre Jolivet   PetscBool      flg;
307c7a4214aSPierre Jolivet   PetscErrorCode ierr;
308c7a4214aSPierre Jolivet 
309c7a4214aSPierre Jolivet   PetscFunctionBegin;
310c7a4214aSPierre Jolivet   ierr = PetscObjectTypeCompare((PetscObject)pv,PETSCVIEWERASCII,&flg);CHKERRQ(ierr);
311c7a4214aSPierre Jolivet   if (flg) {
312c7a4214aSPierre Jolivet     ierr = PetscViewerASCIIPrintf(pv,"symmetry: %c\n",a->hmatrix->get_symmetry_type());CHKERRQ(ierr);
313c7a4214aSPierre Jolivet     if (PetscAbsScalar(a->s-1.0) > PETSC_MACHINE_EPSILON) {
314c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX)
315c7a4214aSPierre Jolivet       ierr = PetscViewerASCIIPrintf(pv,"scaling: %g+%gi\n",(double)PetscRealPart(a->s),(double)PetscImaginaryPart(a->s));CHKERRQ(ierr);
316c7a4214aSPierre Jolivet #else
317c7a4214aSPierre Jolivet       ierr = PetscViewerASCIIPrintf(pv,"scaling: %g\n",(double)a->s);CHKERRQ(ierr);
318c7a4214aSPierre Jolivet #endif
319c7a4214aSPierre Jolivet     }
320c0aa6a63SJacob Faibussowitsch     ierr = PetscViewerASCIIPrintf(pv,"minimum cluster size: %" PetscInt_FMT "\n",a->bs[0]);CHKERRQ(ierr);
321c0aa6a63SJacob Faibussowitsch     ierr = PetscViewerASCIIPrintf(pv,"maximum block size: %" PetscInt_FMT "\n",a->bs[1]);CHKERRQ(ierr);
322c7a4214aSPierre Jolivet     ierr = PetscViewerASCIIPrintf(pv,"epsilon: %g\n",(double)a->epsilon);CHKERRQ(ierr);
323c7a4214aSPierre Jolivet     ierr = PetscViewerASCIIPrintf(pv,"eta: %g\n",(double)a->eta);CHKERRQ(ierr);
324c0aa6a63SJacob Faibussowitsch     ierr = PetscViewerASCIIPrintf(pv,"minimum target depth: %" PetscInt_FMT "\n",a->depth[0]);CHKERRQ(ierr);
325c0aa6a63SJacob Faibussowitsch     ierr = PetscViewerASCIIPrintf(pv,"minimum source depth: %" PetscInt_FMT "\n",a->depth[1]);CHKERRQ(ierr);
326c7a4214aSPierre Jolivet     ierr = PetscViewerASCIIPrintf(pv,"compressor: %s\n",MatHtoolCompressorTypes[a->compressor]);CHKERRQ(ierr);
32798e73e17SPierre Jolivet     ierr = PetscViewerASCIIPrintf(pv,"clustering: %s\n",MatHtoolClusteringTypes[a->clustering]);CHKERRQ(ierr);
328cab00e0dSPierre Jolivet     ierr = PetscViewerASCIIPrintf(pv,"compression ratio: %s\n",a->hmatrix->get_infos("Compression_ratio").c_str());CHKERRQ(ierr);
329cab00e0dSPierre Jolivet     ierr = PetscViewerASCIIPrintf(pv,"space saving: %s\n",a->hmatrix->get_infos("Space_saving").c_str());CHKERRQ(ierr);
330c7a4214aSPierre Jolivet     ierr = 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());CHKERRQ(ierr);
331c7a4214aSPierre Jolivet     ierr = 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());CHKERRQ(ierr);
332c7a4214aSPierre Jolivet     ierr = 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());CHKERRQ(ierr);
333c7a4214aSPierre Jolivet     ierr = 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());CHKERRQ(ierr);
334c7a4214aSPierre Jolivet   }
335c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
336c7a4214aSPierre Jolivet }
337c7a4214aSPierre Jolivet 
338c7a4214aSPierre Jolivet static PetscErrorCode MatScale_Htool(Mat A,PetscScalar s)
339c7a4214aSPierre Jolivet {
340c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
341c7a4214aSPierre Jolivet 
342c7a4214aSPierre Jolivet   PetscFunctionBegin;
343c7a4214aSPierre Jolivet   a->s *= s;
344c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
345c7a4214aSPierre Jolivet }
346c7a4214aSPierre Jolivet 
347c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */
348c7a4214aSPierre Jolivet static PetscErrorCode MatGetRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
349c7a4214aSPierre Jolivet {
350c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
351c7a4214aSPierre Jolivet   PetscInt       *idxc;
352c7a4214aSPierre Jolivet   PetscBLASInt   one = 1,bn;
353c7a4214aSPierre Jolivet   PetscErrorCode ierr;
354c7a4214aSPierre Jolivet 
355c7a4214aSPierre Jolivet   PetscFunctionBegin;
356c7a4214aSPierre Jolivet   if (nz) *nz = A->cmap->N;
357c7a4214aSPierre Jolivet   if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */
358c7a4214aSPierre Jolivet     ierr = PetscMalloc1(A->cmap->N,&idxc);CHKERRQ(ierr);
359c7a4214aSPierre Jolivet     for (PetscInt i=0; i<A->cmap->N; ++i) idxc[i] = i;
360c7a4214aSPierre Jolivet   }
361c7a4214aSPierre Jolivet   if (idx) *idx = idxc;
362c7a4214aSPierre Jolivet   if (v) {
363c7a4214aSPierre Jolivet     ierr = PetscMalloc1(A->cmap->N,v);CHKERRQ(ierr);
364c7a4214aSPierre Jolivet     if (a->wrapper) a->wrapper->copy_submatrix(1,A->cmap->N,&row,idxc,*v);
365cab00e0dSPierre Jolivet     else reinterpret_cast<htool::VirtualGenerator<PetscScalar>*>(a->kernelctx)->copy_submatrix(1,A->cmap->N,&row,idxc,*v);
366c7a4214aSPierre Jolivet     ierr = PetscBLASIntCast(A->cmap->N,&bn);CHKERRQ(ierr);
367c7a4214aSPierre Jolivet     PetscStackCallBLAS("BLASscal",BLASscal_(&bn,&a->s,*v,&one));
368c7a4214aSPierre Jolivet   }
369c7a4214aSPierre Jolivet   if (!idx) {
370c7a4214aSPierre Jolivet     ierr = PetscFree(idxc);CHKERRQ(ierr);
371c7a4214aSPierre Jolivet   }
372c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
373c7a4214aSPierre Jolivet }
374c7a4214aSPierre Jolivet 
375c7a4214aSPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
376c7a4214aSPierre Jolivet {
377c7a4214aSPierre Jolivet   PetscErrorCode ierr;
378c7a4214aSPierre Jolivet 
379c7a4214aSPierre Jolivet   PetscFunctionBegin;
380c7a4214aSPierre Jolivet   if (nz) *nz = 0;
381c7a4214aSPierre Jolivet   if (idx) {
382c7a4214aSPierre Jolivet     ierr = PetscFree(*idx);CHKERRQ(ierr);
383c7a4214aSPierre Jolivet   }
384c7a4214aSPierre Jolivet   if (v) {
385c7a4214aSPierre Jolivet     ierr = PetscFree(*v);CHKERRQ(ierr);
386c7a4214aSPierre Jolivet   }
387c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
388c7a4214aSPierre Jolivet }
389c7a4214aSPierre Jolivet 
390c7a4214aSPierre Jolivet static PetscErrorCode MatSetFromOptions_Htool(PetscOptionItems *PetscOptionsObject,Mat A)
391c7a4214aSPierre Jolivet {
392c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
393c7a4214aSPierre Jolivet   PetscInt       n;
394c7a4214aSPierre Jolivet   PetscBool      flg;
395c7a4214aSPierre Jolivet   PetscErrorCode ierr;
396c7a4214aSPierre Jolivet 
397c7a4214aSPierre Jolivet   PetscFunctionBegin;
398c7a4214aSPierre Jolivet   ierr = PetscOptionsHead(PetscOptionsObject,"Htool options");CHKERRQ(ierr);
399c7a4214aSPierre Jolivet   ierr = PetscOptionsInt("-mat_htool_min_cluster_size","Minimal leaf size in cluster tree",NULL,a->bs[0],a->bs,NULL);CHKERRQ(ierr);
400c7a4214aSPierre Jolivet   ierr = PetscOptionsInt("-mat_htool_max_block_size","Maximal number of coefficients in a dense block",NULL,a->bs[1],a->bs + 1,NULL);CHKERRQ(ierr);
401c7a4214aSPierre Jolivet   ierr = PetscOptionsReal("-mat_htool_epsilon","Relative error in Frobenius norm when approximating a block",NULL,a->epsilon,&a->epsilon,NULL);CHKERRQ(ierr);
402c7a4214aSPierre Jolivet   ierr = PetscOptionsReal("-mat_htool_eta","Admissibility condition tolerance",NULL,a->eta,&a->eta,NULL);CHKERRQ(ierr);
403c7a4214aSPierre Jolivet   ierr = PetscOptionsInt("-mat_htool_min_target_depth","Minimal cluster tree depth associated with the rows",NULL,a->depth[0],a->depth,NULL);CHKERRQ(ierr);
404c7a4214aSPierre Jolivet   ierr = PetscOptionsInt("-mat_htool_min_source_depth","Minimal cluster tree depth associated with the columns",NULL,a->depth[1],a->depth + 1,NULL);CHKERRQ(ierr);
405c7a4214aSPierre Jolivet   n = 0;
40698e73e17SPierre Jolivet   ierr = PetscOptionsEList("-mat_htool_compressor","Type of compression","MatHtoolCompressorType",MatHtoolCompressorTypes,ALEN(MatHtoolCompressorTypes),MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA],&n,&flg);CHKERRQ(ierr);
407c7a4214aSPierre Jolivet   if (flg) a->compressor = MatHtoolCompressorType(n);
40898e73e17SPierre Jolivet   n = 0;
40998e73e17SPierre Jolivet   ierr = PetscOptionsEList("-mat_htool_clustering","Type of clustering","MatHtoolClusteringType",MatHtoolClusteringTypes,ALEN(MatHtoolClusteringTypes),MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR],&n,&flg);CHKERRQ(ierr);
41098e73e17SPierre Jolivet   if (flg) a->clustering = MatHtoolClusteringType(n);
411c7a4214aSPierre Jolivet   ierr = PetscOptionsTail();CHKERRQ(ierr);
412c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
413c7a4214aSPierre Jolivet }
414c7a4214aSPierre Jolivet 
415c7a4214aSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A,MatAssemblyType type)
416c7a4214aSPierre Jolivet {
417c7a4214aSPierre Jolivet   Mat_Htool                                                    *a = (Mat_Htool*)A->data;
418c7a4214aSPierre Jolivet   const PetscInt                                               *ranges;
419c7a4214aSPierre Jolivet   PetscInt                                                     *offset;
420c7a4214aSPierre Jolivet   PetscMPIInt                                                  size;
421c7a4214aSPierre Jolivet   char                                                         S = PetscDefined(USE_COMPLEX) && A->hermitian ? 'H' : (A->symmetric ? 'S' : 'N'),uplo = S == 'N' ? 'N' : 'U';
422cab00e0dSPierre Jolivet   htool::VirtualGenerator<PetscScalar>                         *generator = nullptr;
42398e73e17SPierre Jolivet   std::shared_ptr<htool::VirtualCluster>                       t,s = nullptr;
4243b9338faSPierre Jolivet   std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor = nullptr;
425c7a4214aSPierre Jolivet   PetscErrorCode                                               ierr;
426c7a4214aSPierre Jolivet 
427c7a4214aSPierre Jolivet   PetscFunctionBegin;
428c7a4214aSPierre Jolivet   ierr = PetscCitationsRegister(HtoolCitation,&HtoolCite);CHKERRQ(ierr);
429c7a4214aSPierre Jolivet   delete a->wrapper;
430c7a4214aSPierre Jolivet   delete a->hmatrix;
431c7a4214aSPierre Jolivet   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr);
432c7a4214aSPierre Jolivet   ierr = PetscMalloc1(2*size,&offset);CHKERRQ(ierr);
433c7a4214aSPierre Jolivet   ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr);
434c7a4214aSPierre Jolivet   for (PetscInt i=0; i<size; ++i) {
435c7a4214aSPierre Jolivet     offset[2*i] = ranges[i];
436c7a4214aSPierre Jolivet     offset[2*i+1] = ranges[i+1] - ranges[i];
437c7a4214aSPierre Jolivet   }
43898e73e17SPierre Jolivet   switch (a->clustering) {
43998e73e17SPierre Jolivet   case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
44098e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
44198e73e17SPierre Jolivet     break;
44298e73e17SPierre Jolivet   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
44398e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
44498e73e17SPierre Jolivet     break;
44598e73e17SPierre Jolivet   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
44698e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim);
44798e73e17SPierre Jolivet     break;
44898e73e17SPierre Jolivet   default:
44998e73e17SPierre Jolivet     t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim);
45098e73e17SPierre Jolivet   }
451c7a4214aSPierre Jolivet   t->set_minclustersize(a->bs[0]);
45298e73e17SPierre Jolivet   t->build(A->rmap->N,a->gcoords_target,offset);
453c7a4214aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx);
454c7a4214aSPierre Jolivet   else {
455c7a4214aSPierre Jolivet     a->wrapper = NULL;
456cab00e0dSPierre Jolivet     generator = reinterpret_cast<htool::VirtualGenerator<PetscScalar>*>(a->kernelctx);
457c7a4214aSPierre Jolivet   }
458c7a4214aSPierre Jolivet   if (a->gcoords_target != a->gcoords_source) {
459c7a4214aSPierre Jolivet     ierr = MatGetOwnershipRangesColumn(A,&ranges);CHKERRQ(ierr);
460c7a4214aSPierre Jolivet     for (PetscInt i=0; i<size; ++i) {
461c7a4214aSPierre Jolivet       offset[2*i] = ranges[i];
462c7a4214aSPierre Jolivet       offset[2*i+1] = ranges[i+1] - ranges[i];
463c7a4214aSPierre Jolivet     }
46498e73e17SPierre Jolivet     switch (a->clustering) {
46598e73e17SPierre Jolivet     case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
46698e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
46798e73e17SPierre Jolivet       break;
46898e73e17SPierre Jolivet     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
46998e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
47098e73e17SPierre Jolivet       break;
47198e73e17SPierre Jolivet     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
47298e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim);
47398e73e17SPierre Jolivet       break;
47498e73e17SPierre Jolivet     default:
47598e73e17SPierre Jolivet       s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim);
47698e73e17SPierre Jolivet     }
477c7a4214aSPierre Jolivet     s->set_minclustersize(a->bs[0]);
47898e73e17SPierre Jolivet     s->build(A->cmap->N,a->gcoords_source,offset);
479c7a4214aSPierre Jolivet     S = uplo = 'N';
480c7a4214aSPierre Jolivet   }
481c7a4214aSPierre Jolivet   ierr = PetscFree(offset);CHKERRQ(ierr);
482c7a4214aSPierre Jolivet   switch (a->compressor) {
483c7a4214aSPierre Jolivet   case MAT_HTOOL_COMPRESSOR_FULL_ACA:
4843b9338faSPierre Jolivet     compressor = std::make_shared<htool::fullACA<PetscScalar>>();
485c7a4214aSPierre Jolivet     break;
486c7a4214aSPierre Jolivet   case MAT_HTOOL_COMPRESSOR_SVD:
4873b9338faSPierre Jolivet     compressor = std::make_shared<htool::SVD<PetscScalar>>();
488c7a4214aSPierre Jolivet     break;
489c7a4214aSPierre Jolivet   default:
4903b9338faSPierre Jolivet     compressor = std::make_shared<htool::sympartialACA<PetscScalar>>();
491c7a4214aSPierre Jolivet   }
4923b9338faSPierre Jolivet   a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar>*>(new htool::HMatrix<PetscScalar>(t,s ? s : t,a->epsilon,a->eta,S,uplo));
4933b9338faSPierre Jolivet   a->hmatrix->set_compression(compressor);
494c7a4214aSPierre Jolivet   a->hmatrix->set_maxblocksize(a->bs[1]);
495c7a4214aSPierre Jolivet   a->hmatrix->set_mintargetdepth(a->depth[0]);
496c7a4214aSPierre Jolivet   a->hmatrix->set_minsourcedepth(a->depth[1]);
4973b9338faSPierre Jolivet   if (s) a->hmatrix->build(a->wrapper ? *a->wrapper : *generator,a->gcoords_target,a->gcoords_source);
4983b9338faSPierre Jolivet   else   a->hmatrix->build(a->wrapper ? *a->wrapper : *generator,a->gcoords_target);
499c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
500c7a4214aSPierre Jolivet }
501c7a4214aSPierre Jolivet 
502c7a4214aSPierre Jolivet static PetscErrorCode MatProductNumeric_Htool(Mat C)
503c7a4214aSPierre Jolivet {
504c7a4214aSPierre Jolivet   Mat_Product       *product = C->product;
505c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool*)product->A->data;
506c7a4214aSPierre Jolivet   const PetscScalar *in;
507c7a4214aSPierre Jolivet   PetscScalar       *out;
5088b8ddd11SPierre Jolivet   PetscInt          N,lda;
509c7a4214aSPierre Jolivet   PetscErrorCode    ierr;
510c7a4214aSPierre Jolivet 
511c7a4214aSPierre Jolivet   PetscFunctionBegin;
512c7a4214aSPierre Jolivet   MatCheckProduct(C,1);
513c7a4214aSPierre Jolivet   ierr = MatGetSize(C,NULL,&N);CHKERRQ(ierr);
514c7a4214aSPierre Jolivet   ierr = MatDenseGetLDA(C,&lda);CHKERRQ(ierr);
515*98921bdaSJacob Faibussowitsch   if (lda != C->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")",lda,C->rmap->n);
516c7a4214aSPierre Jolivet   ierr = MatDenseGetArrayRead(product->B,&in);CHKERRQ(ierr);
517c7a4214aSPierre Jolivet   ierr = MatDenseGetArrayWrite(C,&out);CHKERRQ(ierr);
5188b8ddd11SPierre Jolivet   switch (product->type) {
5198b8ddd11SPierre Jolivet   case MATPRODUCT_AB:
520c7a4214aSPierre Jolivet     a->hmatrix->mvprod_local_to_local(in,out,N);
5218b8ddd11SPierre Jolivet     break;
5228b8ddd11SPierre Jolivet   case MATPRODUCT_AtB:
5238b8ddd11SPierre Jolivet     a->hmatrix->mvprod_transp_local_to_local(in,out,N);
524c7a4214aSPierre Jolivet     break;
525c7a4214aSPierre Jolivet   default:
526*98921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType %s is not supported",MatProductTypes[product->type]);
527c7a4214aSPierre Jolivet   }
5288b8ddd11SPierre Jolivet   ierr = MatDenseRestoreArrayWrite(C,&out);CHKERRQ(ierr);
5298b8ddd11SPierre Jolivet   ierr = MatDenseRestoreArrayRead(product->B,&in);CHKERRQ(ierr);
5308b8ddd11SPierre Jolivet   ierr = MatScale(C,a->s);CHKERRQ(ierr);
531c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
532c7a4214aSPierre Jolivet }
533c7a4214aSPierre Jolivet 
534c7a4214aSPierre Jolivet static PetscErrorCode MatProductSymbolic_Htool(Mat C)
535c7a4214aSPierre Jolivet {
536c7a4214aSPierre Jolivet   Mat_Product    *product = C->product;
537c7a4214aSPierre Jolivet   Mat            A,B;
538c7a4214aSPierre Jolivet   PetscBool      flg;
539c7a4214aSPierre Jolivet   PetscErrorCode ierr;
540c7a4214aSPierre Jolivet 
541c7a4214aSPierre Jolivet   PetscFunctionBegin;
542c7a4214aSPierre Jolivet   MatCheckProduct(C,1);
543c7a4214aSPierre Jolivet   A = product->A;
544c7a4214aSPierre Jolivet   B = product->B;
545c7a4214aSPierre Jolivet   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
546*98921bdaSJacob Faibussowitsch   if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MatProduct_AB not supported for %s",((PetscObject)product->B)->type_name);
547c7a4214aSPierre Jolivet   switch (product->type) {
548c7a4214aSPierre Jolivet   case MATPRODUCT_AB:
549c7a4214aSPierre Jolivet     if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) {
550c7a4214aSPierre Jolivet       ierr = MatSetSizes(C,A->rmap->n,B->cmap->n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
551c7a4214aSPierre Jolivet     }
5528b8ddd11SPierre Jolivet     break;
5538b8ddd11SPierre Jolivet   case MATPRODUCT_AtB:
5548b8ddd11SPierre Jolivet     if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) {
5558b8ddd11SPierre Jolivet       ierr = MatSetSizes(C,A->cmap->n,B->cmap->n,A->cmap->N,B->cmap->N);CHKERRQ(ierr);
5568b8ddd11SPierre Jolivet     }
5578b8ddd11SPierre Jolivet     break;
5588b8ddd11SPierre Jolivet   default:
559*98921bdaSJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"ProductType %s is not supported",MatProductTypes[product->type]);
5608b8ddd11SPierre Jolivet   }
561c7a4214aSPierre Jolivet   ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr);
562c7a4214aSPierre Jolivet   ierr = MatSetUp(C);CHKERRQ(ierr);
563c7a4214aSPierre Jolivet   ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
564c7a4214aSPierre Jolivet   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
565c7a4214aSPierre Jolivet   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
566c7a4214aSPierre Jolivet   C->ops->productsymbolic = NULL;
567c7a4214aSPierre Jolivet   C->ops->productnumeric = MatProductNumeric_Htool;
568c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
569c7a4214aSPierre Jolivet }
570c7a4214aSPierre Jolivet 
571c7a4214aSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Htool(Mat C)
572c7a4214aSPierre Jolivet {
573c7a4214aSPierre Jolivet   PetscFunctionBegin;
574c7a4214aSPierre Jolivet   MatCheckProduct(C,1);
5758b8ddd11SPierre Jolivet   if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool;
576c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
577c7a4214aSPierre Jolivet }
578c7a4214aSPierre Jolivet 
57998e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix)
580c7a4214aSPierre Jolivet {
581c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
582c7a4214aSPierre Jolivet 
583c7a4214aSPierre Jolivet   PetscFunctionBegin;
584c7a4214aSPierre Jolivet   *hmatrix = a->hmatrix;
585c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
586c7a4214aSPierre Jolivet }
587c7a4214aSPierre Jolivet 
588c7a4214aSPierre Jolivet /*@C
589c7a4214aSPierre Jolivet      MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a MATHTOOL.
590c7a4214aSPierre Jolivet 
591c7a4214aSPierre Jolivet    Input Parameter:
592c7a4214aSPierre Jolivet .     A - hierarchical matrix
593c7a4214aSPierre Jolivet 
594c7a4214aSPierre Jolivet    Output Parameter:
595c7a4214aSPierre Jolivet .     hmatrix - opaque pointer to a Htool virtual matrix
596c7a4214aSPierre Jolivet 
597c7a4214aSPierre Jolivet    Level: advanced
598c7a4214aSPierre Jolivet 
599c7a4214aSPierre Jolivet .seealso:  MATHTOOL
600c7a4214aSPierre Jolivet @*/
60198e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix)
602c7a4214aSPierre Jolivet {
603c7a4214aSPierre Jolivet   PetscErrorCode ierr;
604c7a4214aSPierre Jolivet 
605c7a4214aSPierre Jolivet   PetscFunctionBegin;
606c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
607c7a4214aSPierre Jolivet   PetscValidPointer(hmatrix,2);
60898e73e17SPierre Jolivet   ierr = PetscTryMethod(A,"MatHtoolGetHierarchicalMat_C",(Mat,const htool::VirtualHMatrix<PetscScalar>**),(A,hmatrix));CHKERRQ(ierr);
609c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
610c7a4214aSPierre Jolivet }
611c7a4214aSPierre Jolivet 
612c7a4214aSPierre Jolivet static PetscErrorCode MatHtoolSetKernel_Htool(Mat A,MatHtoolKernel kernel,void *kernelctx)
613c7a4214aSPierre Jolivet {
614c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
615c7a4214aSPierre Jolivet 
616c7a4214aSPierre Jolivet   PetscFunctionBegin;
617c7a4214aSPierre Jolivet   a->kernel    = kernel;
618c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
619c7a4214aSPierre Jolivet   delete a->wrapper;
620c7a4214aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx);
621c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
622c7a4214aSPierre Jolivet }
623c7a4214aSPierre Jolivet 
624c7a4214aSPierre Jolivet /*@C
62598e73e17SPierre Jolivet      MatHtoolSetKernel - Sets the kernel and context used for the assembly of a MATHTOOL.
626c7a4214aSPierre Jolivet 
627c7a4214aSPierre Jolivet    Input Parameters:
628c7a4214aSPierre Jolivet +     A - hierarchical matrix
629c7a4214aSPierre Jolivet .     kernel - computational kernel (or NULL)
630cab00e0dSPierre Jolivet -     kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
631c7a4214aSPierre Jolivet 
632c7a4214aSPierre Jolivet    Level: advanced
633c7a4214aSPierre Jolivet 
634c7a4214aSPierre Jolivet .seealso:  MATHTOOL, MatCreateHtoolFromKernel()
635c7a4214aSPierre Jolivet @*/
636c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat A,MatHtoolKernel kernel,void *kernelctx)
637c7a4214aSPierre Jolivet {
638c7a4214aSPierre Jolivet   PetscErrorCode ierr;
639c7a4214aSPierre Jolivet 
640c7a4214aSPierre Jolivet   PetscFunctionBegin;
641c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
642c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel,2);
643c7a4214aSPierre Jolivet   if (!kernel)    PetscValidPointer(kernelctx,3);
644c7a4214aSPierre Jolivet   ierr = PetscTryMethod(A,"MatHtoolSetKernel_C",(Mat,MatHtoolKernel,void*),(A,kernel,kernelctx));CHKERRQ(ierr);
645c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
646c7a4214aSPierre Jolivet }
647c7a4214aSPierre Jolivet 
64898e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A,IS* is)
64998e73e17SPierre Jolivet {
65098e73e17SPierre Jolivet   Mat_Htool             *a = (Mat_Htool*)A->data;
65198e73e17SPierre Jolivet   std::vector<PetscInt> source;
65298e73e17SPierre Jolivet   PetscErrorCode        ierr;
65398e73e17SPierre Jolivet 
65498e73e17SPierre Jolivet   PetscFunctionBegin;
655a9918087SPierre Jolivet   source = a->hmatrix->get_source_cluster()->get_local_perm();
65698e73e17SPierre Jolivet   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),source.size(),source.data(),PETSC_COPY_VALUES,is);CHKERRQ(ierr);
65798e73e17SPierre Jolivet   ierr = ISSetPermutation(*is);CHKERRQ(ierr);
65898e73e17SPierre Jolivet   PetscFunctionReturn(0);
65998e73e17SPierre Jolivet }
66098e73e17SPierre Jolivet 
66198e73e17SPierre Jolivet /*@C
66298e73e17SPierre Jolivet      MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster.
66398e73e17SPierre Jolivet 
66498e73e17SPierre Jolivet    Input Parameter:
66598e73e17SPierre Jolivet .     A - hierarchical matrix
66698e73e17SPierre Jolivet 
66798e73e17SPierre Jolivet    Output Parameter:
66898e73e17SPierre Jolivet .     is - permutation
66998e73e17SPierre Jolivet 
67098e73e17SPierre Jolivet    Level: advanced
67198e73e17SPierre Jolivet 
67298e73e17SPierre Jolivet .seealso:  MATHTOOL, MatHtoolGetPermutationTarget(), MatHtoolUsePermutation()
67398e73e17SPierre Jolivet @*/
67498e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat A,IS* is)
67598e73e17SPierre Jolivet {
67698e73e17SPierre Jolivet   PetscErrorCode ierr;
67798e73e17SPierre Jolivet 
67898e73e17SPierre Jolivet   PetscFunctionBegin;
67998e73e17SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
68098e73e17SPierre Jolivet   if (!is) PetscValidPointer(is,2);
68198e73e17SPierre Jolivet   ierr = PetscTryMethod(A,"MatHtoolGetPermutationSource_C",(Mat,IS*),(A,is));CHKERRQ(ierr);
68298e73e17SPierre Jolivet   PetscFunctionReturn(0);
68398e73e17SPierre Jolivet }
68498e73e17SPierre Jolivet 
68598e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A,IS* is)
68698e73e17SPierre Jolivet {
68798e73e17SPierre Jolivet   Mat_Htool             *a = (Mat_Htool*)A->data;
68898e73e17SPierre Jolivet   std::vector<PetscInt> target;
68998e73e17SPierre Jolivet   PetscErrorCode        ierr;
69098e73e17SPierre Jolivet 
69198e73e17SPierre Jolivet   PetscFunctionBegin;
692a9918087SPierre Jolivet   target = a->hmatrix->get_target_cluster()->get_local_perm();
69398e73e17SPierre Jolivet   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),target.size(),target.data(),PETSC_COPY_VALUES,is);CHKERRQ(ierr);
69498e73e17SPierre Jolivet   ierr = ISSetPermutation(*is);CHKERRQ(ierr);
69598e73e17SPierre Jolivet   PetscFunctionReturn(0);
69698e73e17SPierre Jolivet }
69798e73e17SPierre Jolivet 
69898e73e17SPierre Jolivet /*@C
69998e73e17SPierre Jolivet      MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster.
70098e73e17SPierre Jolivet 
70198e73e17SPierre Jolivet    Input Parameter:
70298e73e17SPierre Jolivet .     A - hierarchical matrix
70398e73e17SPierre Jolivet 
70498e73e17SPierre Jolivet    Output Parameter:
70598e73e17SPierre Jolivet .     is - permutation
70698e73e17SPierre Jolivet 
70798e73e17SPierre Jolivet    Level: advanced
70898e73e17SPierre Jolivet 
70998e73e17SPierre Jolivet .seealso:  MATHTOOL, MatHtoolGetPermutationSource(), MatHtoolUsePermutation()
71098e73e17SPierre Jolivet @*/
71198e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat A,IS* is)
71298e73e17SPierre Jolivet {
71398e73e17SPierre Jolivet   PetscErrorCode ierr;
71498e73e17SPierre Jolivet 
71598e73e17SPierre Jolivet   PetscFunctionBegin;
71698e73e17SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
71798e73e17SPierre Jolivet   if (!is) PetscValidPointer(is,2);
71898e73e17SPierre Jolivet   ierr = PetscTryMethod(A,"MatHtoolGetPermutationTarget_C",(Mat,IS*),(A,is));CHKERRQ(ierr);
71998e73e17SPierre Jolivet   PetscFunctionReturn(0);
72098e73e17SPierre Jolivet }
72198e73e17SPierre Jolivet 
72298e73e17SPierre Jolivet static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A,PetscBool use)
72398e73e17SPierre Jolivet {
72498e73e17SPierre Jolivet   Mat_Htool *a = (Mat_Htool*)A->data;
72598e73e17SPierre Jolivet 
72698e73e17SPierre Jolivet   PetscFunctionBegin;
72798e73e17SPierre Jolivet   a->hmatrix->set_use_permutation(use);
72898e73e17SPierre Jolivet   PetscFunctionReturn(0);
72998e73e17SPierre Jolivet }
73098e73e17SPierre Jolivet 
73198e73e17SPierre Jolivet /*@C
73298e73e17SPierre Jolivet      MatHtoolUsePermutation - Sets whether MATHTOOL should permute input (resp. output) vectors following its internal source (resp. target) permutation.
73398e73e17SPierre Jolivet 
73498e73e17SPierre Jolivet    Input Parameters:
73598e73e17SPierre Jolivet +     A - hierarchical matrix
73698e73e17SPierre Jolivet -     use - Boolean value
73798e73e17SPierre Jolivet 
73898e73e17SPierre Jolivet    Level: advanced
73998e73e17SPierre Jolivet 
74098e73e17SPierre Jolivet .seealso:  MATHTOOL, MatHtoolGetPermutationSource(), MatHtoolGetPermutationTarget()
74198e73e17SPierre Jolivet @*/
74298e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat A,PetscBool use)
74398e73e17SPierre Jolivet {
74498e73e17SPierre Jolivet   PetscErrorCode ierr;
74598e73e17SPierre Jolivet 
74698e73e17SPierre Jolivet   PetscFunctionBegin;
74798e73e17SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
74898e73e17SPierre Jolivet   PetscValidLogicalCollectiveBool(A,use,2);
74998e73e17SPierre Jolivet   ierr = PetscTryMethod(A,"MatHtoolUsePermutation_C",(Mat,PetscBool),(A,use));CHKERRQ(ierr);
75098e73e17SPierre Jolivet   PetscFunctionReturn(0);
75198e73e17SPierre Jolivet }
75298e73e17SPierre Jolivet 
753c7a4214aSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B)
754c7a4214aSPierre Jolivet {
755c7a4214aSPierre Jolivet   Mat            C;
756c7a4214aSPierre Jolivet   Mat_Htool      *a = (Mat_Htool*)A->data;
757c7a4214aSPierre Jolivet   PetscInt       lda;
758c7a4214aSPierre Jolivet   PetscScalar    *array;
759c7a4214aSPierre Jolivet   PetscErrorCode ierr;
760c7a4214aSPierre Jolivet 
761c7a4214aSPierre Jolivet   PetscFunctionBegin;
762c7a4214aSPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
763c7a4214aSPierre Jolivet     C = *B;
764c7a4214aSPierre Jolivet     if (C->rmap->n != A->rmap->n || C->cmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible dimensions");
765c7a4214aSPierre Jolivet     ierr = MatDenseGetLDA(C,&lda);CHKERRQ(ierr);
766*98921bdaSJacob Faibussowitsch     if (lda != C->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")",lda,C->rmap->n);
767c7a4214aSPierre Jolivet   } else {
768c7a4214aSPierre Jolivet     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
769c7a4214aSPierre Jolivet     ierr = MatSetSizes(C,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
770c7a4214aSPierre Jolivet     ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr);
771c7a4214aSPierre Jolivet     ierr = MatSetUp(C);CHKERRQ(ierr);
772c7a4214aSPierre Jolivet     ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
773c7a4214aSPierre Jolivet   }
774c7a4214aSPierre Jolivet   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
775c7a4214aSPierre Jolivet   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
776c7a4214aSPierre Jolivet   ierr = MatDenseGetArrayWrite(C,&array);CHKERRQ(ierr);
777c7a4214aSPierre Jolivet   a->hmatrix->copy_local_dense_perm(array);
778c7a4214aSPierre Jolivet   ierr = MatDenseRestoreArrayWrite(C,&array);CHKERRQ(ierr);
779c7a4214aSPierre Jolivet   ierr = MatScale(C,a->s);CHKERRQ(ierr);
780c7a4214aSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
781c7a4214aSPierre Jolivet     ierr = MatHeaderReplace(A,&C);CHKERRQ(ierr);
782c7a4214aSPierre Jolivet   } else *B = C;
783c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
784c7a4214aSPierre Jolivet }
785c7a4214aSPierre Jolivet 
78698e73e17SPierre Jolivet static PetscErrorCode GenEntriesTranspose(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *rows,const PetscInt *cols,PetscScalar *ptr,void *ctx)
787c7a4214aSPierre Jolivet {
788c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose*)ctx;
78998e73e17SPierre Jolivet   PetscScalar             *tmp;
79098e73e17SPierre Jolivet   PetscErrorCode          ierr;
79198e73e17SPierre Jolivet 
79298e73e17SPierre Jolivet   PetscFunctionBegin;
79398e73e17SPierre Jolivet   generator->kernel(sdim,N,M,cols,rows,ptr,generator->kernelctx);
79498e73e17SPierre Jolivet   ierr = PetscMalloc1(M*N,&tmp);CHKERRQ(ierr);
79598e73e17SPierre Jolivet   ierr = PetscArraycpy(tmp,ptr,M*N);CHKERRQ(ierr);
79698e73e17SPierre Jolivet   for (PetscInt i=0; i<M; ++i) {
79798e73e17SPierre Jolivet     for (PetscInt j=0; j<N; ++j) ptr[i+j*M] = tmp[j+i*N];
79898e73e17SPierre Jolivet   }
79998e73e17SPierre Jolivet   ierr = PetscFree(tmp);CHKERRQ(ierr);
80098e73e17SPierre Jolivet   PetscFunctionReturn(0);
801c7a4214aSPierre Jolivet }
802c7a4214aSPierre Jolivet 
803c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */
804c7a4214aSPierre Jolivet static PetscErrorCode MatTranspose_Htool(Mat A,MatReuse reuse,Mat *B)
805c7a4214aSPierre Jolivet {
806c7a4214aSPierre Jolivet   Mat                     C;
807c7a4214aSPierre Jolivet   Mat_Htool               *a = (Mat_Htool*)A->data,*c;
808c7a4214aSPierre Jolivet   PetscInt                M = A->rmap->N,N = A->cmap->N,m = A->rmap->n,n = A->cmap->n;
809c7a4214aSPierre Jolivet   PetscContainer          container;
810c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
811c7a4214aSPierre Jolivet   PetscErrorCode          ierr;
812c7a4214aSPierre Jolivet 
813c7a4214aSPierre Jolivet   PetscFunctionBegin;
814c7a4214aSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatTranspose() with MAT_INPLACE_MATRIX not supported");
815c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
816c7a4214aSPierre Jolivet     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
817c7a4214aSPierre Jolivet     ierr = MatSetSizes(C,n,m,N,M);CHKERRQ(ierr);
818c7a4214aSPierre Jolivet     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
819c7a4214aSPierre Jolivet     ierr = MatSetUp(C);CHKERRQ(ierr);
820c7a4214aSPierre Jolivet     ierr = PetscContainerCreate(PetscObjectComm((PetscObject)C),&container);CHKERRQ(ierr);
821c7a4214aSPierre Jolivet     ierr = PetscNew(&kernelt);CHKERRQ(ierr);
822c7a4214aSPierre Jolivet     ierr = PetscContainerSetPointer(container,kernelt);CHKERRQ(ierr);
823c7a4214aSPierre Jolivet     ierr = PetscObjectCompose((PetscObject)C,"KernelTranspose",(PetscObject)container);CHKERRQ(ierr);
824c7a4214aSPierre Jolivet   } else {
825c7a4214aSPierre Jolivet     C = *B;
826c7a4214aSPierre Jolivet     ierr = PetscObjectQuery((PetscObject)C,"KernelTranspose",(PetscObject*)&container);CHKERRQ(ierr);
827c7a4214aSPierre Jolivet     if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatTranspose() with MAT_INITIAL_MATRIX first");
828c7a4214aSPierre Jolivet     ierr = PetscContainerGetPointer(container,(void**)&kernelt);CHKERRQ(ierr);
829c7a4214aSPierre Jolivet   }
830c7a4214aSPierre Jolivet   c                  = (Mat_Htool*)C->data;
831c7a4214aSPierre Jolivet   c->dim             = a->dim;
832c7a4214aSPierre Jolivet   c->s               = a->s;
83398e73e17SPierre Jolivet   c->kernel          = GenEntriesTranspose;
834c7a4214aSPierre Jolivet   if (kernelt->A != A) {
835c7a4214aSPierre Jolivet     ierr = MatDestroy(&kernelt->A);CHKERRQ(ierr);
836c7a4214aSPierre Jolivet     kernelt->A       = A;
837c7a4214aSPierre Jolivet     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
838c7a4214aSPierre Jolivet   }
839c7a4214aSPierre Jolivet   kernelt->kernel    = a->kernel;
840c7a4214aSPierre Jolivet   kernelt->kernelctx = a->kernelctx;
841c7a4214aSPierre Jolivet   c->kernelctx       = kernelt;
842c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
843c7a4214aSPierre Jolivet     ierr = PetscMalloc1(N*c->dim,&c->gcoords_target);CHKERRQ(ierr);
844c7a4214aSPierre Jolivet     ierr = PetscArraycpy(c->gcoords_target,a->gcoords_source,N*c->dim);CHKERRQ(ierr);
845c7a4214aSPierre Jolivet     if (a->gcoords_target != a->gcoords_source) {
846c7a4214aSPierre Jolivet       ierr = PetscMalloc1(M*c->dim,&c->gcoords_source);CHKERRQ(ierr);
847c7a4214aSPierre Jolivet       ierr = PetscArraycpy(c->gcoords_source,a->gcoords_target,M*c->dim);CHKERRQ(ierr);
848c7a4214aSPierre Jolivet     } else c->gcoords_source = c->gcoords_target;
849c7a4214aSPierre Jolivet     ierr = PetscCalloc2(M,&c->work_source,N,&c->work_target);CHKERRQ(ierr);
850c7a4214aSPierre Jolivet   }
851c7a4214aSPierre Jolivet   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
852c7a4214aSPierre Jolivet   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
853c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) *B = C;
854c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
855c7a4214aSPierre Jolivet }
856c7a4214aSPierre Jolivet 
857c7a4214aSPierre Jolivet /*@C
858c7a4214aSPierre Jolivet      MatCreateHtoolFromKernel - Creates a MATHTOOL from a user-supplied kernel.
859c7a4214aSPierre Jolivet 
860c7a4214aSPierre Jolivet    Input Parameters:
861c7a4214aSPierre Jolivet +     comm - MPI communicator
862c7a4214aSPierre Jolivet .     m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
863c7a4214aSPierre Jolivet .     n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
864c7a4214aSPierre Jolivet .     M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
865c7a4214aSPierre Jolivet .     N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
866c7a4214aSPierre Jolivet .     spacedim - dimension of the space coordinates
867c7a4214aSPierre Jolivet .     coords_target - coordinates of the target
868c7a4214aSPierre Jolivet .     coords_source - coordinates of the source
869c7a4214aSPierre Jolivet .     kernel - computational kernel (or NULL)
870cab00e0dSPierre Jolivet -     kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
871c7a4214aSPierre Jolivet 
872c7a4214aSPierre Jolivet    Output Parameter:
873c7a4214aSPierre Jolivet .     B - matrix
874c7a4214aSPierre Jolivet 
875c7a4214aSPierre Jolivet    Options Database Keys:
876c7a4214aSPierre Jolivet +     -mat_htool_min_cluster_size <PetscInt> - minimal leaf size in cluster tree
877c7a4214aSPierre Jolivet .     -mat_htool_max_block_size <PetscInt> - maximal number of coefficients in a dense block
878c7a4214aSPierre Jolivet .     -mat_htool_epsilon <PetscReal> - relative error in Frobenius norm when approximating a block
879c7a4214aSPierre Jolivet .     -mat_htool_eta <PetscReal> - admissibility condition tolerance
880c7a4214aSPierre Jolivet .     -mat_htool_min_target_depth <PetscInt> - minimal cluster tree depth associated with the rows
881c7a4214aSPierre Jolivet .     -mat_htool_min_source_depth <PetscInt> - minimal cluster tree depth associated with the columns
88298e73e17SPierre Jolivet .     -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression
88398e73e17SPierre Jolivet -     -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering
884c7a4214aSPierre Jolivet 
885c7a4214aSPierre Jolivet    Level: intermediate
886c7a4214aSPierre Jolivet 
88753022affSStefano Zampini .seealso:  MatCreate(), MATHTOOL, PCSetCoordinates(), MatHtoolSetKernel(), MatHtoolCompressorType, MATH2OPUS, MatCreateH2OpusFromKernel()
888c7a4214aSPierre Jolivet @*/
889c7a4214aSPierre 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)
890c7a4214aSPierre Jolivet {
891c7a4214aSPierre Jolivet   Mat            A;
892c7a4214aSPierre Jolivet   Mat_Htool      *a;
893c7a4214aSPierre Jolivet   PetscErrorCode ierr;
894c7a4214aSPierre Jolivet 
895c7a4214aSPierre Jolivet   PetscFunctionBegin;
896c7a4214aSPierre Jolivet   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
897c7a4214aSPierre Jolivet   PetscValidLogicalCollectiveInt(A,spacedim,6);
898c7a4214aSPierre Jolivet   PetscValidRealPointer(coords_target,7);
899c7a4214aSPierre Jolivet   PetscValidRealPointer(coords_source,8);
900c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel,9);
901c7a4214aSPierre Jolivet   if (!kernel)    PetscValidPointer(kernelctx,10);
902c7a4214aSPierre Jolivet   ierr = MatSetSizes(A,m,n,M,N);CHKERRQ(ierr);
903c7a4214aSPierre Jolivet   ierr = MatSetType(A,MATHTOOL);CHKERRQ(ierr);
904c7a4214aSPierre Jolivet   ierr = MatSetUp(A);CHKERRQ(ierr);
905c7a4214aSPierre Jolivet   a            = (Mat_Htool*)A->data;
906c7a4214aSPierre Jolivet   a->dim       = spacedim;
907c7a4214aSPierre Jolivet   a->s         = 1.0;
908c7a4214aSPierre Jolivet   a->kernel    = kernel;
909c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
910c7a4214aSPierre Jolivet   ierr = PetscCalloc1(A->rmap->N*spacedim,&a->gcoords_target);CHKERRQ(ierr);
911c7a4214aSPierre Jolivet   ierr = PetscArraycpy(a->gcoords_target+A->rmap->rstart*spacedim,coords_target,A->rmap->n*spacedim);CHKERRQ(ierr);
912c7a4214aSPierre Jolivet   ierr = MPIU_Allreduce(MPI_IN_PLACE,a->gcoords_target,A->rmap->N*spacedim,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); /* global target coordinates */
913c7a4214aSPierre Jolivet   if (coords_target != coords_source) {
914c7a4214aSPierre Jolivet     ierr = PetscCalloc1(A->cmap->N*spacedim,&a->gcoords_source);CHKERRQ(ierr);
915c7a4214aSPierre Jolivet     ierr = PetscArraycpy(a->gcoords_source+A->cmap->rstart*spacedim,coords_source,A->cmap->n*spacedim);CHKERRQ(ierr);
916c7a4214aSPierre Jolivet     ierr = MPIU_Allreduce(MPI_IN_PLACE,a->gcoords_source,A->cmap->N*spacedim,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); /* global source coordinates */
917c7a4214aSPierre Jolivet   } else a->gcoords_source = a->gcoords_target;
918c7a4214aSPierre Jolivet   ierr = PetscCalloc2(A->cmap->N,&a->work_source,A->rmap->N,&a->work_target);CHKERRQ(ierr);
919c7a4214aSPierre Jolivet   *B = A;
920c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
921c7a4214aSPierre Jolivet }
922c7a4214aSPierre Jolivet 
923c7a4214aSPierre Jolivet /*MC
924c7a4214aSPierre Jolivet      MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package.
925c7a4214aSPierre Jolivet 
926c7a4214aSPierre Jolivet   Use ./configure --download-htool to install PETSc to use Htool.
927c7a4214aSPierre Jolivet 
928c7a4214aSPierre Jolivet    Options Database Keys:
929e85a1b63SPatrick Sanan .     -mat_type htool - matrix type to "htool" during a call to MatSetFromOptions()
930c7a4214aSPierre Jolivet 
931c7a4214aSPierre Jolivet    Level: beginner
932c7a4214aSPierre Jolivet 
93353022affSStefano Zampini .seealso: MATH2OPUS, MATDENSE, MatCreateHtoolFromKernel(), MatHtoolSetKernel()
934c7a4214aSPierre Jolivet M*/
935c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A)
936c7a4214aSPierre Jolivet {
937c7a4214aSPierre Jolivet   Mat_Htool      *a;
938c7a4214aSPierre Jolivet   PetscErrorCode ierr;
939c7a4214aSPierre Jolivet 
940c7a4214aSPierre Jolivet   PetscFunctionBegin;
941c7a4214aSPierre Jolivet   ierr = PetscNewLog(A,&a);CHKERRQ(ierr);
942c7a4214aSPierre Jolivet   A->data = (void*)a;
943c7a4214aSPierre Jolivet   ierr = PetscObjectChangeTypeName((PetscObject)A,MATHTOOL);CHKERRQ(ierr);
944c7a4214aSPierre Jolivet   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
945c7a4214aSPierre Jolivet   A->ops->getdiagonal       = MatGetDiagonal_Htool;
946c7a4214aSPierre Jolivet   A->ops->getdiagonalblock  = MatGetDiagonalBlock_Htool;
947c7a4214aSPierre Jolivet   A->ops->mult              = MatMult_Htool;
948c7a4214aSPierre Jolivet   A->ops->multadd           = MatMultAdd_Htool;
9498b8ddd11SPierre Jolivet   A->ops->multtranspose     = MatMultTranspose_Htool;
9508b8ddd11SPierre Jolivet   if (!PetscDefined(USE_COMPLEX)) A->ops->multhermitiantranspose = MatMultTranspose_Htool;
951c7a4214aSPierre Jolivet   A->ops->increaseoverlap   = MatIncreaseOverlap_Htool;
952c7a4214aSPierre Jolivet   A->ops->createsubmatrices = MatCreateSubMatrices_Htool;
953c7a4214aSPierre Jolivet   A->ops->transpose         = MatTranspose_Htool;
954c7a4214aSPierre Jolivet   A->ops->destroy           = MatDestroy_Htool;
955c7a4214aSPierre Jolivet   A->ops->view              = MatView_Htool;
956c7a4214aSPierre Jolivet   A->ops->setfromoptions    = MatSetFromOptions_Htool;
957c7a4214aSPierre Jolivet   A->ops->scale             = MatScale_Htool;
958c7a4214aSPierre Jolivet   A->ops->getrow            = MatGetRow_Htool;
959c7a4214aSPierre Jolivet   A->ops->restorerow        = MatRestoreRow_Htool;
960c7a4214aSPierre Jolivet   A->ops->assemblyend       = MatAssemblyEnd_Htool;
961c7a4214aSPierre Jolivet   a->dim                    = 0;
962c7a4214aSPierre Jolivet   a->gcoords_target         = NULL;
963c7a4214aSPierre Jolivet   a->gcoords_source         = NULL;
964c7a4214aSPierre Jolivet   a->s                      = 1.0;
965c7a4214aSPierre Jolivet   a->bs[0]                  = 10;
966c7a4214aSPierre Jolivet   a->bs[1]                  = 1000000;
967c7a4214aSPierre Jolivet   a->epsilon                = PetscSqrtReal(PETSC_SMALL);
968c7a4214aSPierre Jolivet   a->eta                    = 10.0;
969c7a4214aSPierre Jolivet   a->depth[0]               = 0;
970c7a4214aSPierre Jolivet   a->depth[1]               = 0;
971c7a4214aSPierre Jolivet   a->compressor             = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA;
972c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",MatProductSetFromOptions_Htool);CHKERRQ(ierr);
973c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",MatProductSetFromOptions_Htool);CHKERRQ(ierr);
974c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",MatConvert_Htool_Dense);CHKERRQ(ierr);
975c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",MatConvert_Htool_Dense);CHKERRQ(ierr);
976c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",MatHtoolGetHierarchicalMat_Htool);CHKERRQ(ierr);
977c7a4214aSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",MatHtoolSetKernel_Htool);CHKERRQ(ierr);
97898e73e17SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",MatHtoolGetPermutationSource_Htool);CHKERRQ(ierr);
97998e73e17SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",MatHtoolGetPermutationTarget_Htool);CHKERRQ(ierr);
98098e73e17SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",MatHtoolUsePermutation_Htool);CHKERRQ(ierr);
981c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
982c7a4214aSPierre Jolivet }
983