1c7a4214aSPierre Jolivet #include <../src/mat/impls/htool/htool.hpp> /*I "petscmat.h" I*/ 2c7a4214aSPierre Jolivet #include <petscblaslapack.h> 3c7a4214aSPierre Jolivet #include <set> 4c7a4214aSPierre Jolivet 5*98e73e17SPierre Jolivet #define ALEN(a) (sizeof(a)/sizeof((a)[0])) 6c7a4214aSPierre Jolivet const char *const MatHtoolCompressorTypes[] = { "sympartialACA", "fullACA", "SVD" }; 7*98e73e17SPierre 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 101c7a4214aSPierre Jolivet static PetscErrorCode MatIncreaseOverlap_Htool(Mat A,PetscInt is_max,IS is[],PetscInt ov) 102c7a4214aSPierre Jolivet { 103c7a4214aSPierre Jolivet std::set<PetscInt> set; 104c7a4214aSPierre Jolivet const PetscInt *idx; 105c7a4214aSPierre Jolivet PetscInt *oidx,size; 106c7a4214aSPierre Jolivet PetscMPIInt csize; 107c7a4214aSPierre Jolivet PetscErrorCode ierr; 108c7a4214aSPierre Jolivet 109c7a4214aSPierre Jolivet PetscFunctionBegin; 110c7a4214aSPierre Jolivet for (PetscInt i=0; i<is_max; ++i) { 111c7a4214aSPierre Jolivet /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */ 112c7a4214aSPierre Jolivet /* needed to avoid subdomain matrices to replicate A since it is dense */ 113c7a4214aSPierre Jolivet ierr = MPI_Comm_size(PetscObjectComm((PetscObject)is[i]),&csize);CHKERRMPI(ierr); 114c7a4214aSPierre Jolivet if (csize != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported parallel IS"); 115c7a4214aSPierre Jolivet ierr = ISGetSize(is[i],&size);CHKERRQ(ierr); 116c7a4214aSPierre Jolivet ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr); 117c7a4214aSPierre Jolivet for (PetscInt j=0; j<size; ++j) { 118c7a4214aSPierre Jolivet set.insert(idx[j]); 119c7a4214aSPierre Jolivet for (PetscInt k=1; k<=ov; ++k) { /* for each layer of overlap */ 120c7a4214aSPierre Jolivet if (idx[j] - k >= 0) set.insert(idx[j] - k); /* do not insert negative indices */ 121c7a4214aSPierre 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 */ 122c7a4214aSPierre Jolivet } 123c7a4214aSPierre Jolivet } 124c7a4214aSPierre Jolivet ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr); 125c7a4214aSPierre Jolivet ierr = ISDestroy(is+i);CHKERRQ(ierr); 126c7a4214aSPierre Jolivet size = set.size(); /* size with overlap */ 127c7a4214aSPierre Jolivet ierr = PetscMalloc1(size,&oidx);CHKERRQ(ierr); 128c7a4214aSPierre Jolivet for (const PetscInt j : set) *oidx++ = j; 129c7a4214aSPierre Jolivet oidx -= size; 130c7a4214aSPierre Jolivet ierr = ISCreateGeneral(PETSC_COMM_SELF,size,oidx,PETSC_OWN_POINTER,is+i);CHKERRQ(ierr); 131c7a4214aSPierre Jolivet } 132c7a4214aSPierre Jolivet PetscFunctionReturn(0); 133c7a4214aSPierre Jolivet } 134c7a4214aSPierre Jolivet 135c7a4214aSPierre Jolivet static PetscErrorCode MatCreateSubMatrices_Htool(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 136c7a4214aSPierre Jolivet { 137c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 138c7a4214aSPierre Jolivet Mat D,B,BT; 139c7a4214aSPierre Jolivet const PetscScalar *copy; 140c7a4214aSPierre Jolivet PetscScalar *ptr; 141c7a4214aSPierre Jolivet const PetscInt *idxr,*idxc,*it; 142c7a4214aSPierre Jolivet PetscInt nrow,m,i; 143c7a4214aSPierre Jolivet PetscBool flg; 144c7a4214aSPierre Jolivet PetscErrorCode ierr; 145c7a4214aSPierre Jolivet 146c7a4214aSPierre Jolivet PetscFunctionBegin; 147c7a4214aSPierre Jolivet if (scall != MAT_REUSE_MATRIX) { 148c7a4214aSPierre Jolivet ierr = PetscCalloc1(n,submat);CHKERRQ(ierr); 149c7a4214aSPierre Jolivet } 150c7a4214aSPierre Jolivet for (i=0; i<n; ++i) { 151c7a4214aSPierre Jolivet ierr = ISGetLocalSize(irow[i],&nrow);CHKERRQ(ierr); 152c7a4214aSPierre Jolivet ierr = ISGetLocalSize(icol[i],&m);CHKERRQ(ierr); 153c7a4214aSPierre Jolivet ierr = ISGetIndices(irow[i],&idxr);CHKERRQ(ierr); 154c7a4214aSPierre Jolivet ierr = ISGetIndices(icol[i],&idxc);CHKERRQ(ierr); 155c7a4214aSPierre Jolivet if (scall != MAT_REUSE_MATRIX) { 156c7a4214aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_SELF,nrow,m,nrow,m,NULL,(*submat)+i);CHKERRQ(ierr); 157c7a4214aSPierre Jolivet } 158c7a4214aSPierre Jolivet ierr = MatDenseGetArrayWrite((*submat)[i],&ptr);CHKERRQ(ierr); 159c7a4214aSPierre Jolivet if (irow[i] == icol[i]) { /* same row and column IS? */ 160c7a4214aSPierre Jolivet ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr); 161c7a4214aSPierre Jolivet if (flg) { 162c7a4214aSPierre Jolivet ierr = ISSorted(irow[i],&flg);CHKERRQ(ierr); 163c7a4214aSPierre Jolivet if (flg) { /* sorted IS? */ 164c7a4214aSPierre Jolivet it = std::lower_bound(idxr,idxr+nrow,A->rmap->rstart); 165c7a4214aSPierre Jolivet if (it != idxr+nrow && *it == A->rmap->rstart) { /* rmap->rstart in IS? */ 166c7a4214aSPierre Jolivet if (std::distance(idxr,it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */ 167c7a4214aSPierre Jolivet for (PetscInt j=0; j<A->rmap->n && flg; ++j) if (PetscUnlikely(it[j] != A->rmap->rstart+j)) flg = PETSC_FALSE; 168c7a4214aSPierre Jolivet if (flg) { /* complete local diagonal block in IS? */ 169c7a4214aSPierre Jolivet /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM 170c7a4214aSPierre Jolivet * [ B C E ] 171c7a4214aSPierre Jolivet * A = [ B D E ] 172c7a4214aSPierre Jolivet * [ B F E ] 173c7a4214aSPierre Jolivet */ 174c7a4214aSPierre Jolivet m = std::distance(idxr,it); /* shift of the coefficient (0,0) of block D from above */ 175c7a4214aSPierre Jolivet ierr = MatGetDiagonalBlock_Htool(A,&D);CHKERRQ(ierr); 176c7a4214aSPierre Jolivet ierr = MatDenseGetArrayRead(D,©);CHKERRQ(ierr); 177c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { 178c7a4214aSPierre Jolivet ierr = PetscArraycpy(ptr+(m+k)*nrow+m,copy+k*A->rmap->n,A->rmap->n);CHKERRQ(ierr); /* block D from above */ 179c7a4214aSPierre Jolivet } 180c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayRead(D,©);CHKERRQ(ierr); 181c7a4214aSPierre Jolivet if (m) { 182c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* vertical block B from above */ 183c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 184c7a4214aSPierre Jolivet if (A->symmetric || A->hermitian) { 185c7a4214aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_SELF,A->rmap->n,m,A->rmap->n,m,ptr+m,&B);CHKERRQ(ierr); 186c7a4214aSPierre Jolivet ierr = MatDenseSetLDA(B,nrow);CHKERRQ(ierr); 187c7a4214aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_SELF,m,A->rmap->n,m,A->rmap->n,ptr+m*nrow,&BT);CHKERRQ(ierr); 188c7a4214aSPierre Jolivet ierr = MatDenseSetLDA(BT,nrow);CHKERRQ(ierr); 189c7a4214aSPierre Jolivet if (A->hermitian && PetscDefined(USE_COMPLEX)) { 190c7a4214aSPierre Jolivet ierr = MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr); 191c7a4214aSPierre Jolivet } else { 192c7a4214aSPierre Jolivet ierr = MatTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr); 193c7a4214aSPierre Jolivet } 194c7a4214aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 195c7a4214aSPierre Jolivet ierr = MatDestroy(&BT);CHKERRQ(ierr); 196c7a4214aSPierre Jolivet } else { 197c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { /* block C from above */ 198c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(m,1,idxr,idxc+m+k,ptr+(m+k)*nrow); 199c7a4214aSPierre Jolivet } 200c7a4214aSPierre Jolivet } 201c7a4214aSPierre Jolivet } 202c7a4214aSPierre Jolivet if (m+A->rmap->n != nrow) { 203c7a4214aSPierre 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 */ 204c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 205c7a4214aSPierre Jolivet if (A->symmetric || A->hermitian) { 206c7a4214aSPierre 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); 207c7a4214aSPierre Jolivet ierr = MatDenseSetLDA(B,nrow);CHKERRQ(ierr); 208c7a4214aSPierre 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); 209c7a4214aSPierre Jolivet ierr = MatDenseSetLDA(BT,nrow);CHKERRQ(ierr); 210c7a4214aSPierre Jolivet if (A->hermitian && PetscDefined(USE_COMPLEX)) { 211c7a4214aSPierre Jolivet ierr = MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr); 212c7a4214aSPierre Jolivet } else { 213c7a4214aSPierre Jolivet ierr = MatTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr); 214c7a4214aSPierre Jolivet } 215c7a4214aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 216c7a4214aSPierre Jolivet ierr = MatDestroy(&BT);CHKERRQ(ierr); 217c7a4214aSPierre Jolivet } else { 218c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { /* block F from above */ 219c7a4214aSPierre 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); 220c7a4214aSPierre Jolivet } 221c7a4214aSPierre Jolivet } 222c7a4214aSPierre Jolivet } 223c7a4214aSPierre Jolivet } /* complete local diagonal block not in IS */ 224c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */ 225c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* rmap->rstart not in IS */ 226c7a4214aSPierre Jolivet } /* unsorted IS */ 227c7a4214aSPierre Jolivet } 228c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* different row and column IS */ 229c7a4214aSPierre Jolivet if (!flg) a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* reassemble everything */ 230c7a4214aSPierre Jolivet ierr = ISRestoreIndices(irow[i],&idxr);CHKERRQ(ierr); 231c7a4214aSPierre Jolivet ierr = ISRestoreIndices(icol[i],&idxc);CHKERRQ(ierr); 232c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayWrite((*submat)[i],&ptr);CHKERRQ(ierr); 233c7a4214aSPierre Jolivet ierr = MatAssemblyBegin((*submat)[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 234c7a4214aSPierre Jolivet ierr = MatAssemblyEnd((*submat)[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 235c7a4214aSPierre Jolivet ierr = MatScale((*submat)[i],a->s);CHKERRQ(ierr); 236c7a4214aSPierre Jolivet } 237c7a4214aSPierre Jolivet PetscFunctionReturn(0); 238c7a4214aSPierre Jolivet } 239c7a4214aSPierre Jolivet 240c7a4214aSPierre Jolivet static PetscErrorCode MatDestroy_Htool(Mat A) 241c7a4214aSPierre Jolivet { 242c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 243c7a4214aSPierre Jolivet PetscContainer container; 244c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 245c7a4214aSPierre Jolivet PetscErrorCode ierr; 246c7a4214aSPierre Jolivet 247c7a4214aSPierre Jolivet PetscFunctionBegin; 248c7a4214aSPierre Jolivet ierr = PetscObjectChangeTypeName((PetscObject)A,NULL);CHKERRQ(ierr); 249c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",NULL);CHKERRQ(ierr); 250c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",NULL);CHKERRQ(ierr); 251c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",NULL);CHKERRQ(ierr); 252c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",NULL);CHKERRQ(ierr); 253c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",NULL);CHKERRQ(ierr); 254c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",NULL);CHKERRQ(ierr); 255*98e73e17SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",NULL);CHKERRQ(ierr); 256*98e73e17SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",NULL);CHKERRQ(ierr); 257*98e73e17SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",NULL);CHKERRQ(ierr); 258c7a4214aSPierre Jolivet ierr = PetscObjectQuery((PetscObject)A,"KernelTranspose",(PetscObject*)&container);CHKERRQ(ierr); 259c7a4214aSPierre Jolivet if (container) { /* created in MatTranspose_Htool() */ 260c7a4214aSPierre Jolivet ierr = PetscContainerGetPointer(container,(void**)&kernelt);CHKERRQ(ierr); 261c7a4214aSPierre Jolivet ierr = MatDestroy(&kernelt->A);CHKERRQ(ierr); 262c7a4214aSPierre Jolivet ierr = PetscFree(kernelt);CHKERRQ(ierr); 263c7a4214aSPierre Jolivet ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 264c7a4214aSPierre Jolivet ierr = PetscObjectCompose((PetscObject)A,"KernelTranspose",NULL);CHKERRQ(ierr); 265c7a4214aSPierre Jolivet } 266c7a4214aSPierre Jolivet if (a->gcoords_source != a->gcoords_target) { 267c7a4214aSPierre Jolivet ierr = PetscFree(a->gcoords_source);CHKERRQ(ierr); 268c7a4214aSPierre Jolivet } 269c7a4214aSPierre Jolivet ierr = PetscFree(a->gcoords_target);CHKERRQ(ierr); 270c7a4214aSPierre Jolivet ierr = PetscFree2(a->work_source,a->work_target);CHKERRQ(ierr); 271c7a4214aSPierre Jolivet delete a->wrapper; 272c7a4214aSPierre Jolivet delete a->hmatrix; 273c7a4214aSPierre Jolivet ierr = PetscFree(A->data);CHKERRQ(ierr); 274c7a4214aSPierre Jolivet PetscFunctionReturn(0); 275c7a4214aSPierre Jolivet } 276c7a4214aSPierre Jolivet 277c7a4214aSPierre Jolivet static PetscErrorCode MatView_Htool(Mat A,PetscViewer pv) 278c7a4214aSPierre Jolivet { 279c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 280c7a4214aSPierre Jolivet PetscBool flg; 281c7a4214aSPierre Jolivet PetscErrorCode ierr; 282c7a4214aSPierre Jolivet 283c7a4214aSPierre Jolivet PetscFunctionBegin; 284c7a4214aSPierre Jolivet ierr = PetscObjectTypeCompare((PetscObject)pv,PETSCVIEWERASCII,&flg);CHKERRQ(ierr); 285c7a4214aSPierre Jolivet if (flg) { 286c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"symmetry: %c\n",a->hmatrix->get_symmetry_type());CHKERRQ(ierr); 287c7a4214aSPierre Jolivet if (PetscAbsScalar(a->s-1.0) > PETSC_MACHINE_EPSILON) { 288c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX) 289c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"scaling: %g+%gi\n",(double)PetscRealPart(a->s),(double)PetscImaginaryPart(a->s));CHKERRQ(ierr); 290c7a4214aSPierre Jolivet #else 291c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"scaling: %g\n",(double)a->s);CHKERRQ(ierr); 292c7a4214aSPierre Jolivet #endif 293c7a4214aSPierre Jolivet } 294c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"minimum cluster size: %D\n",a->bs[0]);CHKERRQ(ierr); 295c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"maximum block size: %D\n",a->bs[1]);CHKERRQ(ierr); 296c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"epsilon: %g\n",(double)a->epsilon);CHKERRQ(ierr); 297c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"eta: %g\n",(double)a->eta);CHKERRQ(ierr); 298c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"minimum target depth: %D\n",a->depth[0]);CHKERRQ(ierr); 299c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"minimum source depth: %D\n",a->depth[1]);CHKERRQ(ierr); 300c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"compressor: %s\n",MatHtoolCompressorTypes[a->compressor]);CHKERRQ(ierr); 301*98e73e17SPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"clustering: %s\n",MatHtoolClusteringTypes[a->clustering]);CHKERRQ(ierr); 302c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"compression: %s\n",a->hmatrix->get_infos("Compression").c_str());CHKERRQ(ierr); 303c7a4214aSPierre 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); 304c7a4214aSPierre 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); 305c7a4214aSPierre 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); 306c7a4214aSPierre 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); 307c7a4214aSPierre Jolivet } 308c7a4214aSPierre Jolivet PetscFunctionReturn(0); 309c7a4214aSPierre Jolivet } 310c7a4214aSPierre Jolivet 311c7a4214aSPierre Jolivet static PetscErrorCode MatScale_Htool(Mat A,PetscScalar s) 312c7a4214aSPierre Jolivet { 313c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 314c7a4214aSPierre Jolivet 315c7a4214aSPierre Jolivet PetscFunctionBegin; 316c7a4214aSPierre Jolivet a->s *= s; 317c7a4214aSPierre Jolivet PetscFunctionReturn(0); 318c7a4214aSPierre Jolivet } 319c7a4214aSPierre Jolivet 320c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */ 321c7a4214aSPierre Jolivet static PetscErrorCode MatGetRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 322c7a4214aSPierre Jolivet { 323c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 324c7a4214aSPierre Jolivet PetscInt *idxc; 325c7a4214aSPierre Jolivet PetscBLASInt one = 1,bn; 326c7a4214aSPierre Jolivet PetscErrorCode ierr; 327c7a4214aSPierre Jolivet 328c7a4214aSPierre Jolivet PetscFunctionBegin; 329c7a4214aSPierre Jolivet if (nz) *nz = A->cmap->N; 330c7a4214aSPierre Jolivet if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */ 331c7a4214aSPierre Jolivet ierr = PetscMalloc1(A->cmap->N,&idxc);CHKERRQ(ierr); 332c7a4214aSPierre Jolivet for (PetscInt i=0; i<A->cmap->N; ++i) idxc[i] = i; 333c7a4214aSPierre Jolivet } 334c7a4214aSPierre Jolivet if (idx) *idx = idxc; 335c7a4214aSPierre Jolivet if (v) { 336c7a4214aSPierre Jolivet ierr = PetscMalloc1(A->cmap->N,v);CHKERRQ(ierr); 337c7a4214aSPierre Jolivet if (a->wrapper) a->wrapper->copy_submatrix(1,A->cmap->N,&row,idxc,*v); 338c7a4214aSPierre Jolivet else reinterpret_cast<htool::IMatrix<PetscScalar>*>(a->kernelctx)->copy_submatrix(1,A->cmap->N,&row,idxc,*v); 339c7a4214aSPierre Jolivet ierr = PetscBLASIntCast(A->cmap->N,&bn);CHKERRQ(ierr); 340c7a4214aSPierre Jolivet PetscStackCallBLAS("BLASscal",BLASscal_(&bn,&a->s,*v,&one)); 341c7a4214aSPierre Jolivet } 342c7a4214aSPierre Jolivet if (!idx) { 343c7a4214aSPierre Jolivet ierr = PetscFree(idxc);CHKERRQ(ierr); 344c7a4214aSPierre Jolivet } 345c7a4214aSPierre Jolivet PetscFunctionReturn(0); 346c7a4214aSPierre Jolivet } 347c7a4214aSPierre Jolivet 348c7a4214aSPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 349c7a4214aSPierre Jolivet { 350c7a4214aSPierre Jolivet PetscErrorCode ierr; 351c7a4214aSPierre Jolivet 352c7a4214aSPierre Jolivet PetscFunctionBegin; 353c7a4214aSPierre Jolivet if (nz) *nz = 0; 354c7a4214aSPierre Jolivet if (idx) { 355c7a4214aSPierre Jolivet ierr = PetscFree(*idx);CHKERRQ(ierr); 356c7a4214aSPierre Jolivet } 357c7a4214aSPierre Jolivet if (v) { 358c7a4214aSPierre Jolivet ierr = PetscFree(*v);CHKERRQ(ierr); 359c7a4214aSPierre Jolivet } 360c7a4214aSPierre Jolivet PetscFunctionReturn(0); 361c7a4214aSPierre Jolivet } 362c7a4214aSPierre Jolivet 363c7a4214aSPierre Jolivet static PetscErrorCode MatSetFromOptions_Htool(PetscOptionItems *PetscOptionsObject,Mat A) 364c7a4214aSPierre Jolivet { 365c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 366c7a4214aSPierre Jolivet PetscInt n; 367c7a4214aSPierre Jolivet PetscBool flg; 368c7a4214aSPierre Jolivet PetscErrorCode ierr; 369c7a4214aSPierre Jolivet 370c7a4214aSPierre Jolivet PetscFunctionBegin; 371c7a4214aSPierre Jolivet ierr = PetscOptionsHead(PetscOptionsObject,"Htool options");CHKERRQ(ierr); 372c7a4214aSPierre Jolivet ierr = PetscOptionsInt("-mat_htool_min_cluster_size","Minimal leaf size in cluster tree",NULL,a->bs[0],a->bs,NULL);CHKERRQ(ierr); 373c7a4214aSPierre 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); 374c7a4214aSPierre Jolivet ierr = PetscOptionsReal("-mat_htool_epsilon","Relative error in Frobenius norm when approximating a block",NULL,a->epsilon,&a->epsilon,NULL);CHKERRQ(ierr); 375c7a4214aSPierre Jolivet ierr = PetscOptionsReal("-mat_htool_eta","Admissibility condition tolerance",NULL,a->eta,&a->eta,NULL);CHKERRQ(ierr); 376c7a4214aSPierre 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); 377c7a4214aSPierre 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); 378c7a4214aSPierre Jolivet n = 0; 379*98e73e17SPierre Jolivet ierr = PetscOptionsEList("-mat_htool_compressor","Type of compression","MatHtoolCompressorType",MatHtoolCompressorTypes,ALEN(MatHtoolCompressorTypes),MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA],&n,&flg);CHKERRQ(ierr); 380c7a4214aSPierre Jolivet if (flg) a->compressor = MatHtoolCompressorType(n); 381*98e73e17SPierre Jolivet n = 0; 382*98e73e17SPierre Jolivet ierr = PetscOptionsEList("-mat_htool_clustering","Type of clustering","MatHtoolClusteringType",MatHtoolClusteringTypes,ALEN(MatHtoolClusteringTypes),MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR],&n,&flg);CHKERRQ(ierr); 383*98e73e17SPierre Jolivet if (flg) a->clustering = MatHtoolClusteringType(n); 384c7a4214aSPierre Jolivet ierr = PetscOptionsTail();CHKERRQ(ierr); 385c7a4214aSPierre Jolivet PetscFunctionReturn(0); 386c7a4214aSPierre Jolivet } 387c7a4214aSPierre Jolivet 388c7a4214aSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A,MatAssemblyType type) 389c7a4214aSPierre Jolivet { 390c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 391c7a4214aSPierre Jolivet const PetscInt *ranges; 392c7a4214aSPierre Jolivet PetscInt *offset; 393c7a4214aSPierre Jolivet PetscMPIInt size; 394c7a4214aSPierre Jolivet char S = PetscDefined(USE_COMPLEX) && A->hermitian ? 'H' : (A->symmetric ? 'S' : 'N'),uplo = S == 'N' ? 'N' : 'U'; 395c7a4214aSPierre Jolivet htool::IMatrix<PetscScalar> *generator = nullptr; 396*98e73e17SPierre Jolivet std::shared_ptr<htool::VirtualCluster> t,s = nullptr; 397c7a4214aSPierre Jolivet PetscErrorCode ierr; 398c7a4214aSPierre Jolivet 399c7a4214aSPierre Jolivet PetscFunctionBegin; 400c7a4214aSPierre Jolivet ierr = PetscCitationsRegister(HtoolCitation,&HtoolCite);CHKERRQ(ierr); 401c7a4214aSPierre Jolivet delete a->wrapper; 402c7a4214aSPierre Jolivet delete a->hmatrix; 403c7a4214aSPierre Jolivet ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr); 404c7a4214aSPierre Jolivet ierr = PetscMalloc1(2*size,&offset);CHKERRQ(ierr); 405c7a4214aSPierre Jolivet ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 406c7a4214aSPierre Jolivet for (PetscInt i=0; i<size; ++i) { 407c7a4214aSPierre Jolivet offset[2*i] = ranges[i]; 408c7a4214aSPierre Jolivet offset[2*i+1] = ranges[i+1] - ranges[i]; 409c7a4214aSPierre Jolivet } 410*98e73e17SPierre Jolivet switch (a->clustering) { 411*98e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 412*98e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 413*98e73e17SPierre Jolivet break; 414*98e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 415*98e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 416*98e73e17SPierre Jolivet break; 417*98e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 418*98e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 419*98e73e17SPierre Jolivet break; 420*98e73e17SPierre Jolivet default: 421*98e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 422*98e73e17SPierre Jolivet } 423c7a4214aSPierre Jolivet t->set_minclustersize(a->bs[0]); 424*98e73e17SPierre Jolivet t->build(A->rmap->N,a->gcoords_target,offset); 425c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx); 426c7a4214aSPierre Jolivet else { 427c7a4214aSPierre Jolivet a->wrapper = NULL; 428c7a4214aSPierre Jolivet generator = reinterpret_cast<htool::IMatrix<PetscScalar>*>(a->kernelctx); 429c7a4214aSPierre Jolivet } 430c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 431c7a4214aSPierre Jolivet ierr = MatGetOwnershipRangesColumn(A,&ranges);CHKERRQ(ierr); 432c7a4214aSPierre Jolivet for (PetscInt i=0; i<size; ++i) { 433c7a4214aSPierre Jolivet offset[2*i] = ranges[i]; 434c7a4214aSPierre Jolivet offset[2*i+1] = ranges[i+1] - ranges[i]; 435c7a4214aSPierre Jolivet } 436*98e73e17SPierre Jolivet switch (a->clustering) { 437*98e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 438*98e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 439*98e73e17SPierre Jolivet break; 440*98e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 441*98e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 442*98e73e17SPierre Jolivet break; 443*98e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 444*98e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 445*98e73e17SPierre Jolivet break; 446*98e73e17SPierre Jolivet default: 447*98e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 448*98e73e17SPierre Jolivet } 449c7a4214aSPierre Jolivet s->set_minclustersize(a->bs[0]); 450*98e73e17SPierre Jolivet s->build(A->cmap->N,a->gcoords_source,offset); 451c7a4214aSPierre Jolivet S = uplo = 'N'; 452c7a4214aSPierre Jolivet } 453c7a4214aSPierre Jolivet ierr = PetscFree(offset);CHKERRQ(ierr); 454c7a4214aSPierre Jolivet switch (a->compressor) { 455c7a4214aSPierre Jolivet case MAT_HTOOL_COMPRESSOR_FULL_ACA: 456*98e73e17SPierre Jolivet a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar>*>(new htool::HMatrix<PetscScalar,htool::fullACA,htool::RjasanowSteinbach>(t,s?s:t,a->epsilon,a->eta,S,uplo)); 457c7a4214aSPierre Jolivet break; 458c7a4214aSPierre Jolivet case MAT_HTOOL_COMPRESSOR_SVD: 459*98e73e17SPierre Jolivet a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar>*>(new htool::HMatrix<PetscScalar,htool::SVD,htool::RjasanowSteinbach>(t,s?s:t,a->epsilon,a->eta,S,uplo)); 460c7a4214aSPierre Jolivet break; 461c7a4214aSPierre Jolivet default: 462*98e73e17SPierre Jolivet a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar>*>(new htool::HMatrix<PetscScalar,htool::sympartialACA,htool::RjasanowSteinbach>(t,s?s:t,a->epsilon,a->eta,S,uplo)); 463c7a4214aSPierre Jolivet } 464c7a4214aSPierre Jolivet a->hmatrix->set_maxblocksize(a->bs[1]); 465c7a4214aSPierre Jolivet a->hmatrix->set_mintargetdepth(a->depth[0]); 466c7a4214aSPierre Jolivet a->hmatrix->set_minsourcedepth(a->depth[1]); 467c7a4214aSPierre Jolivet if (s) a->hmatrix->build_auto(a->wrapper ? *a->wrapper : *generator,a->gcoords_target,a->gcoords_source); 468c7a4214aSPierre Jolivet else a->hmatrix->build_auto_sym(a->wrapper ? *a->wrapper : *generator,a->gcoords_target); 469c7a4214aSPierre Jolivet PetscFunctionReturn(0); 470c7a4214aSPierre Jolivet } 471c7a4214aSPierre Jolivet 472c7a4214aSPierre Jolivet static PetscErrorCode MatProductNumeric_Htool(Mat C) 473c7a4214aSPierre Jolivet { 474c7a4214aSPierre Jolivet Mat_Product *product = C->product; 475c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)product->A->data; 476c7a4214aSPierre Jolivet const PetscScalar *in; 477c7a4214aSPierre Jolivet PetscScalar *out; 478c7a4214aSPierre Jolivet PetscInt lda; 479c7a4214aSPierre Jolivet PetscErrorCode ierr; 480c7a4214aSPierre Jolivet 481c7a4214aSPierre Jolivet PetscFunctionBegin; 482c7a4214aSPierre Jolivet MatCheckProduct(C,1); 483c7a4214aSPierre Jolivet switch (product->type) { 484c7a4214aSPierre Jolivet case MATPRODUCT_AB: 485c7a4214aSPierre Jolivet PetscInt N; 486c7a4214aSPierre Jolivet ierr = MatGetSize(C,NULL,&N);CHKERRQ(ierr); 487c7a4214aSPierre Jolivet ierr = MatDenseGetLDA(C,&lda);CHKERRQ(ierr); 488c7a4214aSPierre Jolivet if (lda != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%D != %D)",lda,C->rmap->n); 489c7a4214aSPierre Jolivet ierr = MatDenseGetArrayRead(product->B,&in);CHKERRQ(ierr); 490c7a4214aSPierre Jolivet ierr = MatDenseGetArrayWrite(C,&out);CHKERRQ(ierr); 491c7a4214aSPierre Jolivet a->hmatrix->mvprod_local_to_local(in,out,N); 492c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayWrite(C,&out);CHKERRQ(ierr); 493c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayRead(product->B,&in);CHKERRQ(ierr); 494c7a4214aSPierre Jolivet ierr = MatScale(C,a->s);CHKERRQ(ierr); 495c7a4214aSPierre Jolivet break; 496c7a4214aSPierre Jolivet default: 497c7a4214aSPierre Jolivet SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType %s is not supported",MatProductTypes[product->type]); 498c7a4214aSPierre Jolivet } 499c7a4214aSPierre Jolivet PetscFunctionReturn(0); 500c7a4214aSPierre Jolivet } 501c7a4214aSPierre Jolivet 502c7a4214aSPierre Jolivet static PetscErrorCode MatProductSymbolic_Htool(Mat C) 503c7a4214aSPierre Jolivet { 504c7a4214aSPierre Jolivet Mat_Product *product = C->product; 505c7a4214aSPierre Jolivet Mat A,B; 506c7a4214aSPierre Jolivet PetscBool flg; 507c7a4214aSPierre Jolivet PetscErrorCode ierr; 508c7a4214aSPierre Jolivet 509c7a4214aSPierre Jolivet PetscFunctionBegin; 510c7a4214aSPierre Jolivet MatCheckProduct(C,1); 511c7a4214aSPierre Jolivet A = product->A; 512c7a4214aSPierre Jolivet B = product->B; 513c7a4214aSPierre Jolivet ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 514c7a4214aSPierre Jolivet if (!flg) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MatProduct_AB not supported for %s",((PetscObject)product->B)->type_name); 515c7a4214aSPierre Jolivet switch (product->type) { 516c7a4214aSPierre Jolivet case MATPRODUCT_AB: 517c7a4214aSPierre Jolivet if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) { 518c7a4214aSPierre Jolivet ierr = MatSetSizes(C,A->rmap->n,B->cmap->n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 519c7a4214aSPierre Jolivet } 520c7a4214aSPierre Jolivet ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr); 521c7a4214aSPierre Jolivet ierr = MatSetUp(C);CHKERRQ(ierr); 522c7a4214aSPierre Jolivet ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 523c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 524c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 525c7a4214aSPierre Jolivet break; 526c7a4214aSPierre Jolivet default: 527c7a4214aSPierre Jolivet SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"ProductType %s is not supported",MatProductTypes[product->type]); 528c7a4214aSPierre Jolivet } 529c7a4214aSPierre Jolivet C->ops->productsymbolic = NULL; 530c7a4214aSPierre Jolivet C->ops->productnumeric = MatProductNumeric_Htool; 531c7a4214aSPierre Jolivet PetscFunctionReturn(0); 532c7a4214aSPierre Jolivet } 533c7a4214aSPierre Jolivet 534c7a4214aSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Htool(Mat C) 535c7a4214aSPierre Jolivet { 536c7a4214aSPierre Jolivet PetscFunctionBegin; 537c7a4214aSPierre Jolivet MatCheckProduct(C,1); 538c7a4214aSPierre Jolivet if (C->product->type == MATPRODUCT_AB) C->ops->productsymbolic = MatProductSymbolic_Htool; 539c7a4214aSPierre Jolivet PetscFunctionReturn(0); 540c7a4214aSPierre Jolivet } 541c7a4214aSPierre Jolivet 542*98e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix) 543c7a4214aSPierre Jolivet { 544c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 545c7a4214aSPierre Jolivet 546c7a4214aSPierre Jolivet PetscFunctionBegin; 547c7a4214aSPierre Jolivet *hmatrix = a->hmatrix; 548c7a4214aSPierre Jolivet PetscFunctionReturn(0); 549c7a4214aSPierre Jolivet } 550c7a4214aSPierre Jolivet 551c7a4214aSPierre Jolivet /*@C 552c7a4214aSPierre Jolivet MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a MATHTOOL. 553c7a4214aSPierre Jolivet 554c7a4214aSPierre Jolivet Input Parameter: 555c7a4214aSPierre Jolivet . A - hierarchical matrix 556c7a4214aSPierre Jolivet 557c7a4214aSPierre Jolivet Output Parameter: 558c7a4214aSPierre Jolivet . hmatrix - opaque pointer to a Htool virtual matrix 559c7a4214aSPierre Jolivet 560c7a4214aSPierre Jolivet Level: advanced 561c7a4214aSPierre Jolivet 562c7a4214aSPierre Jolivet .seealso: MATHTOOL 563c7a4214aSPierre Jolivet @*/ 564*98e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix) 565c7a4214aSPierre Jolivet { 566c7a4214aSPierre Jolivet PetscErrorCode ierr; 567c7a4214aSPierre Jolivet 568c7a4214aSPierre Jolivet PetscFunctionBegin; 569c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 570c7a4214aSPierre Jolivet PetscValidPointer(hmatrix,2); 571*98e73e17SPierre Jolivet ierr = PetscTryMethod(A,"MatHtoolGetHierarchicalMat_C",(Mat,const htool::VirtualHMatrix<PetscScalar>**),(A,hmatrix));CHKERRQ(ierr); 572c7a4214aSPierre Jolivet PetscFunctionReturn(0); 573c7a4214aSPierre Jolivet } 574c7a4214aSPierre Jolivet 575c7a4214aSPierre Jolivet static PetscErrorCode MatHtoolSetKernel_Htool(Mat A,MatHtoolKernel kernel,void *kernelctx) 576c7a4214aSPierre Jolivet { 577c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 578c7a4214aSPierre Jolivet 579c7a4214aSPierre Jolivet PetscFunctionBegin; 580c7a4214aSPierre Jolivet a->kernel = kernel; 581c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 582c7a4214aSPierre Jolivet delete a->wrapper; 583c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx); 584c7a4214aSPierre Jolivet PetscFunctionReturn(0); 585c7a4214aSPierre Jolivet } 586c7a4214aSPierre Jolivet 587c7a4214aSPierre Jolivet /*@C 588*98e73e17SPierre Jolivet MatHtoolSetKernel - Sets the kernel and context used for the assembly of a MATHTOOL. 589c7a4214aSPierre Jolivet 590c7a4214aSPierre Jolivet Input Parameters: 591c7a4214aSPierre Jolivet + A - hierarchical matrix 592c7a4214aSPierre Jolivet . kernel - computational kernel (or NULL) 593c7a4214aSPierre Jolivet - kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::IMatrix<PetscScalar>*) 594c7a4214aSPierre Jolivet 595c7a4214aSPierre Jolivet Level: advanced 596c7a4214aSPierre Jolivet 597c7a4214aSPierre Jolivet .seealso: MATHTOOL, MatCreateHtoolFromKernel() 598c7a4214aSPierre Jolivet @*/ 599c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat A,MatHtoolKernel kernel,void *kernelctx) 600c7a4214aSPierre Jolivet { 601c7a4214aSPierre Jolivet PetscErrorCode ierr; 602c7a4214aSPierre Jolivet 603c7a4214aSPierre Jolivet PetscFunctionBegin; 604c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 605c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel,2); 606c7a4214aSPierre Jolivet if (!kernel) PetscValidPointer(kernelctx,3); 607c7a4214aSPierre Jolivet ierr = PetscTryMethod(A,"MatHtoolSetKernel_C",(Mat,MatHtoolKernel,void*),(A,kernel,kernelctx));CHKERRQ(ierr); 608c7a4214aSPierre Jolivet PetscFunctionReturn(0); 609c7a4214aSPierre Jolivet } 610c7a4214aSPierre Jolivet 611*98e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A,IS* is) 612*98e73e17SPierre Jolivet { 613*98e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 614*98e73e17SPierre Jolivet std::vector<PetscInt> source; 615*98e73e17SPierre Jolivet PetscErrorCode ierr; 616*98e73e17SPierre Jolivet 617*98e73e17SPierre Jolivet PetscFunctionBegin; 618*98e73e17SPierre Jolivet source = a->hmatrix->get_local_perm_source(); 619*98e73e17SPierre Jolivet ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),source.size(),source.data(),PETSC_COPY_VALUES,is);CHKERRQ(ierr); 620*98e73e17SPierre Jolivet ierr = ISSetPermutation(*is);CHKERRQ(ierr); 621*98e73e17SPierre Jolivet PetscFunctionReturn(0); 622*98e73e17SPierre Jolivet } 623*98e73e17SPierre Jolivet 624*98e73e17SPierre Jolivet /*@C 625*98e73e17SPierre Jolivet MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster. 626*98e73e17SPierre Jolivet 627*98e73e17SPierre Jolivet Input Parameter: 628*98e73e17SPierre Jolivet . A - hierarchical matrix 629*98e73e17SPierre Jolivet 630*98e73e17SPierre Jolivet Output Parameter: 631*98e73e17SPierre Jolivet . is - permutation 632*98e73e17SPierre Jolivet 633*98e73e17SPierre Jolivet Level: advanced 634*98e73e17SPierre Jolivet 635*98e73e17SPierre Jolivet .seealso: MATHTOOL, MatHtoolGetPermutationTarget(), MatHtoolUsePermutation() 636*98e73e17SPierre Jolivet @*/ 637*98e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat A,IS* is) 638*98e73e17SPierre Jolivet { 639*98e73e17SPierre Jolivet PetscErrorCode ierr; 640*98e73e17SPierre Jolivet 641*98e73e17SPierre Jolivet PetscFunctionBegin; 642*98e73e17SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 643*98e73e17SPierre Jolivet if (!is) PetscValidPointer(is,2); 644*98e73e17SPierre Jolivet ierr = PetscTryMethod(A,"MatHtoolGetPermutationSource_C",(Mat,IS*),(A,is));CHKERRQ(ierr); 645*98e73e17SPierre Jolivet PetscFunctionReturn(0); 646*98e73e17SPierre Jolivet } 647*98e73e17SPierre Jolivet 648*98e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A,IS* is) 649*98e73e17SPierre Jolivet { 650*98e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 651*98e73e17SPierre Jolivet std::vector<PetscInt> target; 652*98e73e17SPierre Jolivet PetscErrorCode ierr; 653*98e73e17SPierre Jolivet 654*98e73e17SPierre Jolivet PetscFunctionBegin; 655*98e73e17SPierre Jolivet target = a->hmatrix->get_local_perm_target(); 656*98e73e17SPierre Jolivet ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),target.size(),target.data(),PETSC_COPY_VALUES,is);CHKERRQ(ierr); 657*98e73e17SPierre Jolivet ierr = ISSetPermutation(*is);CHKERRQ(ierr); 658*98e73e17SPierre Jolivet PetscFunctionReturn(0); 659*98e73e17SPierre Jolivet } 660*98e73e17SPierre Jolivet 661*98e73e17SPierre Jolivet /*@C 662*98e73e17SPierre Jolivet MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster. 663*98e73e17SPierre Jolivet 664*98e73e17SPierre Jolivet Input Parameter: 665*98e73e17SPierre Jolivet . A - hierarchical matrix 666*98e73e17SPierre Jolivet 667*98e73e17SPierre Jolivet Output Parameter: 668*98e73e17SPierre Jolivet . is - permutation 669*98e73e17SPierre Jolivet 670*98e73e17SPierre Jolivet Level: advanced 671*98e73e17SPierre Jolivet 672*98e73e17SPierre Jolivet .seealso: MATHTOOL, MatHtoolGetPermutationSource(), MatHtoolUsePermutation() 673*98e73e17SPierre Jolivet @*/ 674*98e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat A,IS* is) 675*98e73e17SPierre Jolivet { 676*98e73e17SPierre Jolivet PetscErrorCode ierr; 677*98e73e17SPierre Jolivet 678*98e73e17SPierre Jolivet PetscFunctionBegin; 679*98e73e17SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 680*98e73e17SPierre Jolivet if (!is) PetscValidPointer(is,2); 681*98e73e17SPierre Jolivet ierr = PetscTryMethod(A,"MatHtoolGetPermutationTarget_C",(Mat,IS*),(A,is));CHKERRQ(ierr); 682*98e73e17SPierre Jolivet PetscFunctionReturn(0); 683*98e73e17SPierre Jolivet } 684*98e73e17SPierre Jolivet 685*98e73e17SPierre Jolivet static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A,PetscBool use) 686*98e73e17SPierre Jolivet { 687*98e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 688*98e73e17SPierre Jolivet 689*98e73e17SPierre Jolivet PetscFunctionBegin; 690*98e73e17SPierre Jolivet a->hmatrix->set_use_permutation(use); 691*98e73e17SPierre Jolivet PetscFunctionReturn(0); 692*98e73e17SPierre Jolivet } 693*98e73e17SPierre Jolivet 694*98e73e17SPierre Jolivet /*@C 695*98e73e17SPierre Jolivet MatHtoolUsePermutation - Sets whether MATHTOOL should permute input (resp. output) vectors following its internal source (resp. target) permutation. 696*98e73e17SPierre Jolivet 697*98e73e17SPierre Jolivet Input Parameters: 698*98e73e17SPierre Jolivet + A - hierarchical matrix 699*98e73e17SPierre Jolivet - use - Boolean value 700*98e73e17SPierre Jolivet 701*98e73e17SPierre Jolivet Level: advanced 702*98e73e17SPierre Jolivet 703*98e73e17SPierre Jolivet .seealso: MATHTOOL, MatHtoolGetPermutationSource(), MatHtoolGetPermutationTarget() 704*98e73e17SPierre Jolivet @*/ 705*98e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat A,PetscBool use) 706*98e73e17SPierre Jolivet { 707*98e73e17SPierre Jolivet PetscErrorCode ierr; 708*98e73e17SPierre Jolivet 709*98e73e17SPierre Jolivet PetscFunctionBegin; 710*98e73e17SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 711*98e73e17SPierre Jolivet PetscValidLogicalCollectiveBool(A,use,2); 712*98e73e17SPierre Jolivet ierr = PetscTryMethod(A,"MatHtoolUsePermutation_C",(Mat,PetscBool),(A,use));CHKERRQ(ierr); 713*98e73e17SPierre Jolivet PetscFunctionReturn(0); 714*98e73e17SPierre Jolivet } 715*98e73e17SPierre Jolivet 716c7a4214aSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B) 717c7a4214aSPierre Jolivet { 718c7a4214aSPierre Jolivet Mat C; 719c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 720c7a4214aSPierre Jolivet PetscInt lda; 721c7a4214aSPierre Jolivet PetscScalar *array; 722c7a4214aSPierre Jolivet PetscErrorCode ierr; 723c7a4214aSPierre Jolivet 724c7a4214aSPierre Jolivet PetscFunctionBegin; 725c7a4214aSPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 726c7a4214aSPierre Jolivet C = *B; 727c7a4214aSPierre Jolivet if (C->rmap->n != A->rmap->n || C->cmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible dimensions"); 728c7a4214aSPierre Jolivet ierr = MatDenseGetLDA(C,&lda);CHKERRQ(ierr); 729c7a4214aSPierre Jolivet if (lda != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%D != %D)",lda,C->rmap->n); 730c7a4214aSPierre Jolivet } else { 731c7a4214aSPierre Jolivet ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 732c7a4214aSPierre Jolivet ierr = MatSetSizes(C,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 733c7a4214aSPierre Jolivet ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr); 734c7a4214aSPierre Jolivet ierr = MatSetUp(C);CHKERRQ(ierr); 735c7a4214aSPierre Jolivet ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 736c7a4214aSPierre Jolivet } 737c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 738c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 739c7a4214aSPierre Jolivet ierr = MatDenseGetArrayWrite(C,&array);CHKERRQ(ierr); 740c7a4214aSPierre Jolivet a->hmatrix->copy_local_dense_perm(array); 741c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayWrite(C,&array);CHKERRQ(ierr); 742c7a4214aSPierre Jolivet ierr = MatScale(C,a->s);CHKERRQ(ierr); 743c7a4214aSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 744c7a4214aSPierre Jolivet ierr = MatHeaderReplace(A,&C);CHKERRQ(ierr); 745c7a4214aSPierre Jolivet } else *B = C; 746c7a4214aSPierre Jolivet PetscFunctionReturn(0); 747c7a4214aSPierre Jolivet } 748c7a4214aSPierre Jolivet 749*98e73e17SPierre Jolivet static PetscErrorCode GenEntriesTranspose(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *rows,const PetscInt *cols,PetscScalar *ptr,void *ctx) 750c7a4214aSPierre Jolivet { 751c7a4214aSPierre Jolivet MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose*)ctx; 752*98e73e17SPierre Jolivet PetscScalar *tmp; 753*98e73e17SPierre Jolivet PetscErrorCode ierr; 754*98e73e17SPierre Jolivet 755*98e73e17SPierre Jolivet PetscFunctionBegin; 756*98e73e17SPierre Jolivet generator->kernel(sdim,N,M,cols,rows,ptr,generator->kernelctx); 757*98e73e17SPierre Jolivet ierr = PetscMalloc1(M*N,&tmp);CHKERRQ(ierr); 758*98e73e17SPierre Jolivet ierr = PetscArraycpy(tmp,ptr,M*N);CHKERRQ(ierr); 759*98e73e17SPierre Jolivet for (PetscInt i=0; i<M; ++i) { 760*98e73e17SPierre Jolivet for (PetscInt j=0; j<N; ++j) ptr[i+j*M] = tmp[j+i*N]; 761*98e73e17SPierre Jolivet } 762*98e73e17SPierre Jolivet ierr = PetscFree(tmp);CHKERRQ(ierr); 763*98e73e17SPierre Jolivet PetscFunctionReturn(0); 764c7a4214aSPierre Jolivet } 765c7a4214aSPierre Jolivet 766c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */ 767c7a4214aSPierre Jolivet static PetscErrorCode MatTranspose_Htool(Mat A,MatReuse reuse,Mat *B) 768c7a4214aSPierre Jolivet { 769c7a4214aSPierre Jolivet Mat C; 770c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data,*c; 771c7a4214aSPierre Jolivet PetscInt M = A->rmap->N,N = A->cmap->N,m = A->rmap->n,n = A->cmap->n; 772c7a4214aSPierre Jolivet PetscContainer container; 773c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 774c7a4214aSPierre Jolivet PetscErrorCode ierr; 775c7a4214aSPierre Jolivet 776c7a4214aSPierre Jolivet PetscFunctionBegin; 777c7a4214aSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatTranspose() with MAT_INPLACE_MATRIX not supported"); 778c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 779c7a4214aSPierre Jolivet ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 780c7a4214aSPierre Jolivet ierr = MatSetSizes(C,n,m,N,M);CHKERRQ(ierr); 781c7a4214aSPierre Jolivet ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr); 782c7a4214aSPierre Jolivet ierr = MatSetUp(C);CHKERRQ(ierr); 783c7a4214aSPierre Jolivet ierr = PetscContainerCreate(PetscObjectComm((PetscObject)C),&container);CHKERRQ(ierr); 784c7a4214aSPierre Jolivet ierr = PetscNew(&kernelt);CHKERRQ(ierr); 785c7a4214aSPierre Jolivet ierr = PetscContainerSetPointer(container,kernelt);CHKERRQ(ierr); 786c7a4214aSPierre Jolivet ierr = PetscObjectCompose((PetscObject)C,"KernelTranspose",(PetscObject)container);CHKERRQ(ierr); 787c7a4214aSPierre Jolivet } else { 788c7a4214aSPierre Jolivet C = *B; 789c7a4214aSPierre Jolivet ierr = PetscObjectQuery((PetscObject)C,"KernelTranspose",(PetscObject*)&container);CHKERRQ(ierr); 790c7a4214aSPierre Jolivet if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatTranspose() with MAT_INITIAL_MATRIX first"); 791c7a4214aSPierre Jolivet ierr = PetscContainerGetPointer(container,(void**)&kernelt);CHKERRQ(ierr); 792c7a4214aSPierre Jolivet } 793c7a4214aSPierre Jolivet c = (Mat_Htool*)C->data; 794c7a4214aSPierre Jolivet c->dim = a->dim; 795c7a4214aSPierre Jolivet c->s = a->s; 796*98e73e17SPierre Jolivet c->kernel = GenEntriesTranspose; 797c7a4214aSPierre Jolivet if (kernelt->A != A) { 798c7a4214aSPierre Jolivet ierr = MatDestroy(&kernelt->A);CHKERRQ(ierr); 799c7a4214aSPierre Jolivet kernelt->A = A; 800c7a4214aSPierre Jolivet ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 801c7a4214aSPierre Jolivet } 802c7a4214aSPierre Jolivet kernelt->kernel = a->kernel; 803c7a4214aSPierre Jolivet kernelt->kernelctx = a->kernelctx; 804c7a4214aSPierre Jolivet c->kernelctx = kernelt; 805c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 806c7a4214aSPierre Jolivet ierr = PetscMalloc1(N*c->dim,&c->gcoords_target);CHKERRQ(ierr); 807c7a4214aSPierre Jolivet ierr = PetscArraycpy(c->gcoords_target,a->gcoords_source,N*c->dim);CHKERRQ(ierr); 808c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 809c7a4214aSPierre Jolivet ierr = PetscMalloc1(M*c->dim,&c->gcoords_source);CHKERRQ(ierr); 810c7a4214aSPierre Jolivet ierr = PetscArraycpy(c->gcoords_source,a->gcoords_target,M*c->dim);CHKERRQ(ierr); 811c7a4214aSPierre Jolivet } else c->gcoords_source = c->gcoords_target; 812c7a4214aSPierre Jolivet ierr = PetscCalloc2(M,&c->work_source,N,&c->work_target);CHKERRQ(ierr); 813c7a4214aSPierre Jolivet } 814c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 815c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 816c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) *B = C; 817c7a4214aSPierre Jolivet PetscFunctionReturn(0); 818c7a4214aSPierre Jolivet } 819c7a4214aSPierre Jolivet 820c7a4214aSPierre Jolivet /*@C 821c7a4214aSPierre Jolivet MatCreateHtoolFromKernel - Creates a MATHTOOL from a user-supplied kernel. 822c7a4214aSPierre Jolivet 823c7a4214aSPierre Jolivet Input Parameters: 824c7a4214aSPierre Jolivet + comm - MPI communicator 825c7a4214aSPierre Jolivet . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 826c7a4214aSPierre Jolivet . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 827c7a4214aSPierre Jolivet . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 828c7a4214aSPierre Jolivet . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 829c7a4214aSPierre Jolivet . spacedim - dimension of the space coordinates 830c7a4214aSPierre Jolivet . coords_target - coordinates of the target 831c7a4214aSPierre Jolivet . coords_source - coordinates of the source 832c7a4214aSPierre Jolivet . kernel - computational kernel (or NULL) 833c7a4214aSPierre Jolivet - kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::IMatrix<PetscScalar>*) 834c7a4214aSPierre Jolivet 835c7a4214aSPierre Jolivet Output Parameter: 836c7a4214aSPierre Jolivet . B - matrix 837c7a4214aSPierre Jolivet 838c7a4214aSPierre Jolivet Options Database Keys: 839c7a4214aSPierre Jolivet + -mat_htool_min_cluster_size <PetscInt> - minimal leaf size in cluster tree 840c7a4214aSPierre Jolivet . -mat_htool_max_block_size <PetscInt> - maximal number of coefficients in a dense block 841c7a4214aSPierre Jolivet . -mat_htool_epsilon <PetscReal> - relative error in Frobenius norm when approximating a block 842c7a4214aSPierre Jolivet . -mat_htool_eta <PetscReal> - admissibility condition tolerance 843c7a4214aSPierre Jolivet . -mat_htool_min_target_depth <PetscInt> - minimal cluster tree depth associated with the rows 844c7a4214aSPierre Jolivet . -mat_htool_min_source_depth <PetscInt> - minimal cluster tree depth associated with the columns 845*98e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression 846*98e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering 847c7a4214aSPierre Jolivet 848c7a4214aSPierre Jolivet Level: intermediate 849c7a4214aSPierre Jolivet 850c7a4214aSPierre Jolivet .seealso: MatCreate(), MATHTOOL, PCSetCoordinates(), MatHtoolSetKernel(), MatHtoolCompressorType, MATHARA, MatCreateHaraFromKernel() 851c7a4214aSPierre Jolivet @*/ 852c7a4214aSPierre 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) 853c7a4214aSPierre Jolivet { 854c7a4214aSPierre Jolivet Mat A; 855c7a4214aSPierre Jolivet Mat_Htool *a; 856c7a4214aSPierre Jolivet PetscErrorCode ierr; 857c7a4214aSPierre Jolivet 858c7a4214aSPierre Jolivet PetscFunctionBegin; 859c7a4214aSPierre Jolivet ierr = MatCreate(comm,&A);CHKERRQ(ierr); 860c7a4214aSPierre Jolivet PetscValidLogicalCollectiveInt(A,spacedim,6); 861c7a4214aSPierre Jolivet PetscValidRealPointer(coords_target,7); 862c7a4214aSPierre Jolivet PetscValidRealPointer(coords_source,8); 863c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel,9); 864c7a4214aSPierre Jolivet if (!kernel) PetscValidPointer(kernelctx,10); 865c7a4214aSPierre Jolivet ierr = MatSetSizes(A,m,n,M,N);CHKERRQ(ierr); 866c7a4214aSPierre Jolivet ierr = MatSetType(A,MATHTOOL);CHKERRQ(ierr); 867c7a4214aSPierre Jolivet ierr = MatSetUp(A);CHKERRQ(ierr); 868c7a4214aSPierre Jolivet a = (Mat_Htool*)A->data; 869c7a4214aSPierre Jolivet a->dim = spacedim; 870c7a4214aSPierre Jolivet a->s = 1.0; 871c7a4214aSPierre Jolivet a->kernel = kernel; 872c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 873c7a4214aSPierre Jolivet ierr = PetscCalloc1(A->rmap->N*spacedim,&a->gcoords_target);CHKERRQ(ierr); 874c7a4214aSPierre Jolivet ierr = PetscArraycpy(a->gcoords_target+A->rmap->rstart*spacedim,coords_target,A->rmap->n*spacedim);CHKERRQ(ierr); 875c7a4214aSPierre 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 */ 876c7a4214aSPierre Jolivet if (coords_target != coords_source) { 877c7a4214aSPierre Jolivet ierr = PetscCalloc1(A->cmap->N*spacedim,&a->gcoords_source);CHKERRQ(ierr); 878c7a4214aSPierre Jolivet ierr = PetscArraycpy(a->gcoords_source+A->cmap->rstart*spacedim,coords_source,A->cmap->n*spacedim);CHKERRQ(ierr); 879c7a4214aSPierre 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 */ 880c7a4214aSPierre Jolivet } else a->gcoords_source = a->gcoords_target; 881c7a4214aSPierre Jolivet ierr = PetscCalloc2(A->cmap->N,&a->work_source,A->rmap->N,&a->work_target);CHKERRQ(ierr); 882c7a4214aSPierre Jolivet *B = A; 883c7a4214aSPierre Jolivet PetscFunctionReturn(0); 884c7a4214aSPierre Jolivet } 885c7a4214aSPierre Jolivet 886c7a4214aSPierre Jolivet /*MC 887c7a4214aSPierre Jolivet MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package. 888c7a4214aSPierre Jolivet 889c7a4214aSPierre Jolivet Use ./configure --download-htool to install PETSc to use Htool. 890c7a4214aSPierre Jolivet 891c7a4214aSPierre Jolivet Options Database Keys: 892e85a1b63SPatrick Sanan . -mat_type htool - matrix type to "htool" during a call to MatSetFromOptions() 893c7a4214aSPierre Jolivet 894c7a4214aSPierre Jolivet Level: beginner 895c7a4214aSPierre Jolivet 896c7a4214aSPierre Jolivet .seealso: MATHARA, MATDENSE, MatCreateHtoolFromKernel(), MatHtoolSetKernel() 897c7a4214aSPierre Jolivet M*/ 898c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A) 899c7a4214aSPierre Jolivet { 900c7a4214aSPierre Jolivet Mat_Htool *a; 901c7a4214aSPierre Jolivet PetscErrorCode ierr; 902c7a4214aSPierre Jolivet 903c7a4214aSPierre Jolivet PetscFunctionBegin; 904c7a4214aSPierre Jolivet ierr = PetscNewLog(A,&a);CHKERRQ(ierr); 905c7a4214aSPierre Jolivet A->data = (void*)a; 906c7a4214aSPierre Jolivet ierr = PetscObjectChangeTypeName((PetscObject)A,MATHTOOL);CHKERRQ(ierr); 907c7a4214aSPierre Jolivet ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 908c7a4214aSPierre Jolivet A->ops->getdiagonal = MatGetDiagonal_Htool; 909c7a4214aSPierre Jolivet A->ops->getdiagonalblock = MatGetDiagonalBlock_Htool; 910c7a4214aSPierre Jolivet A->ops->mult = MatMult_Htool; 911c7a4214aSPierre Jolivet A->ops->multadd = MatMultAdd_Htool; 912c7a4214aSPierre Jolivet A->ops->increaseoverlap = MatIncreaseOverlap_Htool; 913c7a4214aSPierre Jolivet A->ops->createsubmatrices = MatCreateSubMatrices_Htool; 914c7a4214aSPierre Jolivet A->ops->transpose = MatTranspose_Htool; 915c7a4214aSPierre Jolivet A->ops->destroy = MatDestroy_Htool; 916c7a4214aSPierre Jolivet A->ops->view = MatView_Htool; 917c7a4214aSPierre Jolivet A->ops->setfromoptions = MatSetFromOptions_Htool; 918c7a4214aSPierre Jolivet A->ops->scale = MatScale_Htool; 919c7a4214aSPierre Jolivet A->ops->getrow = MatGetRow_Htool; 920c7a4214aSPierre Jolivet A->ops->restorerow = MatRestoreRow_Htool; 921c7a4214aSPierre Jolivet A->ops->assemblyend = MatAssemblyEnd_Htool; 922c7a4214aSPierre Jolivet a->dim = 0; 923c7a4214aSPierre Jolivet a->gcoords_target = NULL; 924c7a4214aSPierre Jolivet a->gcoords_source = NULL; 925c7a4214aSPierre Jolivet a->s = 1.0; 926c7a4214aSPierre Jolivet a->bs[0] = 10; 927c7a4214aSPierre Jolivet a->bs[1] = 1000000; 928c7a4214aSPierre Jolivet a->epsilon = PetscSqrtReal(PETSC_SMALL); 929c7a4214aSPierre Jolivet a->eta = 10.0; 930c7a4214aSPierre Jolivet a->depth[0] = 0; 931c7a4214aSPierre Jolivet a->depth[1] = 0; 932c7a4214aSPierre Jolivet a->compressor = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA; 933c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",MatProductSetFromOptions_Htool);CHKERRQ(ierr); 934c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",MatProductSetFromOptions_Htool);CHKERRQ(ierr); 935c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",MatConvert_Htool_Dense);CHKERRQ(ierr); 936c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",MatConvert_Htool_Dense);CHKERRQ(ierr); 937c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",MatHtoolGetHierarchicalMat_Htool);CHKERRQ(ierr); 938c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",MatHtoolSetKernel_Htool);CHKERRQ(ierr); 939*98e73e17SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",MatHtoolGetPermutationSource_Htool);CHKERRQ(ierr); 940*98e73e17SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",MatHtoolGetPermutationTarget_Htool);CHKERRQ(ierr); 941*98e73e17SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",MatHtoolUsePermutation_Htool);CHKERRQ(ierr); 942c7a4214aSPierre Jolivet PetscFunctionReturn(0); 943c7a4214aSPierre Jolivet } 944