1c7a4214aSPierre Jolivet #include <../src/mat/impls/htool/htool.hpp> /*I "petscmat.h" I*/ 2c7a4214aSPierre Jolivet #include <petscblaslapack.h> 3c7a4214aSPierre Jolivet #include <set> 4c7a4214aSPierre Jolivet 5c7a4214aSPierre Jolivet const char *const MatHtoolCompressorTypes[] = { "sympartialACA", "fullACA", "SVD" }; 6c7a4214aSPierre Jolivet const char HtoolCitation[] = "@article{marchand2020two,\n" 7c7a4214aSPierre Jolivet " Author = {Marchand, Pierre and Claeys, Xavier and Jolivet, Pierre and Nataf, Fr\\'ed\\'eric and Tournier, Pierre-Henri},\n" 8c7a4214aSPierre Jolivet " Title = {Two-level preconditioning for $h$-version boundary element approximation of hypersingular operator with {GenEO}},\n" 9c7a4214aSPierre Jolivet " Year = {2020},\n" 10c7a4214aSPierre Jolivet " Publisher = {Elsevier},\n" 11c7a4214aSPierre Jolivet " Journal = {Numerische Mathematik},\n" 12c7a4214aSPierre Jolivet " Volume = {146},\n" 13c7a4214aSPierre Jolivet " Pages = {597--628},\n" 14c7a4214aSPierre Jolivet " Url = {https://github.com/htool-ddm/htool}\n" 15c7a4214aSPierre Jolivet "}\n"; 16c7a4214aSPierre Jolivet static PetscBool HtoolCite = PETSC_FALSE; 17c7a4214aSPierre Jolivet 18c7a4214aSPierre Jolivet static PetscErrorCode MatGetDiagonal_Htool(Mat A,Vec v) 19c7a4214aSPierre Jolivet { 20c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 21c7a4214aSPierre Jolivet PetscScalar *x; 22c7a4214aSPierre Jolivet PetscBool flg; 23c7a4214aSPierre Jolivet PetscErrorCode ierr; 24c7a4214aSPierre Jolivet 25c7a4214aSPierre Jolivet PetscFunctionBegin; 26c7a4214aSPierre Jolivet ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr); 27c7a4214aSPierre Jolivet if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only congruent layouts supported"); 28c7a4214aSPierre Jolivet ierr = VecGetArrayWrite(v,&x);CHKERRQ(ierr); 29c7a4214aSPierre Jolivet a->hmatrix->copy_local_diagonal(x); 30c7a4214aSPierre Jolivet ierr = VecRestoreArrayWrite(v,&x);CHKERRQ(ierr); 31c7a4214aSPierre Jolivet ierr = VecScale(v,a->s);CHKERRQ(ierr); 32c7a4214aSPierre Jolivet PetscFunctionReturn(0); 33c7a4214aSPierre Jolivet } 34c7a4214aSPierre Jolivet 35c7a4214aSPierre Jolivet static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A,Mat *b) 36c7a4214aSPierre Jolivet { 37c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 38c7a4214aSPierre Jolivet Mat B; 39c7a4214aSPierre Jolivet PetscScalar *ptr; 40c7a4214aSPierre Jolivet PetscBool flg; 41c7a4214aSPierre Jolivet PetscErrorCode ierr; 42c7a4214aSPierre Jolivet 43c7a4214aSPierre Jolivet PetscFunctionBegin; 44c7a4214aSPierre Jolivet ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr); 45c7a4214aSPierre Jolivet if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only congruent layouts supported"); 46c7a4214aSPierre Jolivet ierr = PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B);CHKERRQ(ierr); /* same logic as in MatGetDiagonalBlock_MPIDense() */ 47c7a4214aSPierre Jolivet if (!B) { 48c7a4214aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_SELF,A->rmap->n,A->rmap->n,A->rmap->n,A->rmap->n,NULL,&B);CHKERRQ(ierr); 49c7a4214aSPierre Jolivet ierr = MatDenseGetArrayWrite(B,&ptr);CHKERRQ(ierr); 50c7a4214aSPierre Jolivet a->hmatrix->copy_local_diagonal_block(ptr); 51c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayWrite(B,&ptr);CHKERRQ(ierr); 52c7a4214aSPierre Jolivet ierr = MatPropagateSymmetryOptions(A,B);CHKERRQ(ierr); 53c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 54c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 55c7a4214aSPierre Jolivet ierr = MatScale(B,a->s);CHKERRQ(ierr); 56c7a4214aSPierre Jolivet ierr = PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B);CHKERRQ(ierr); 57c7a4214aSPierre Jolivet *b = B; 58c7a4214aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 59c7a4214aSPierre Jolivet } else *b = B; 60c7a4214aSPierre Jolivet PetscFunctionReturn(0); 61c7a4214aSPierre Jolivet } 62c7a4214aSPierre Jolivet 63c7a4214aSPierre Jolivet static PetscErrorCode MatMult_Htool(Mat A,Vec x,Vec y) 64c7a4214aSPierre Jolivet { 65c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 66c7a4214aSPierre Jolivet const PetscScalar *in; 67c7a4214aSPierre Jolivet PetscScalar *out; 68c7a4214aSPierre Jolivet PetscErrorCode ierr; 69c7a4214aSPierre Jolivet 70c7a4214aSPierre Jolivet PetscFunctionBegin; 71c7a4214aSPierre Jolivet ierr = VecGetArrayRead(x,&in);CHKERRQ(ierr); 72c7a4214aSPierre Jolivet ierr = VecGetArrayWrite(y,&out);CHKERRQ(ierr); 73c7a4214aSPierre Jolivet a->hmatrix->mvprod_local_to_local(in,out); 74c7a4214aSPierre Jolivet ierr = VecRestoreArrayRead(x,&in);CHKERRQ(ierr); 75c7a4214aSPierre Jolivet ierr = VecRestoreArrayWrite(y,&out);CHKERRQ(ierr); 76c7a4214aSPierre Jolivet ierr = VecScale(y,a->s);CHKERRQ(ierr); 77c7a4214aSPierre Jolivet PetscFunctionReturn(0); 78c7a4214aSPierre Jolivet } 79c7a4214aSPierre Jolivet 80c7a4214aSPierre Jolivet /* naive implementation of MatMultAdd() needed for FEM-BEM coupling via MATNEST */ 81c7a4214aSPierre Jolivet static PetscErrorCode MatMultAdd_Htool(Mat A,Vec v1,Vec v2,Vec v3) 82c7a4214aSPierre Jolivet { 83c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 84c7a4214aSPierre Jolivet Vec tmp; 85c7a4214aSPierre Jolivet const PetscScalar scale = a->s; 86c7a4214aSPierre Jolivet PetscErrorCode ierr; 87c7a4214aSPierre Jolivet 88c7a4214aSPierre Jolivet PetscFunctionBegin; 89c7a4214aSPierre Jolivet ierr = VecDuplicate(v2,&tmp);CHKERRQ(ierr); 90c7a4214aSPierre 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 */ 91c7a4214aSPierre Jolivet a->s = 1.0; /* set s to 1.0 since VecAXPY() may be used to scale the MatMult() output Vec */ 92c7a4214aSPierre Jolivet ierr = MatMult_Htool(A,v1,tmp);CHKERRQ(ierr); 93c7a4214aSPierre Jolivet ierr = VecAXPY(v3,scale,tmp);CHKERRQ(ierr); 94c7a4214aSPierre Jolivet ierr = VecDestroy(&tmp);CHKERRQ(ierr); 95c7a4214aSPierre Jolivet a->s = scale; /* set s back to its original value */ 96c7a4214aSPierre Jolivet PetscFunctionReturn(0); 97c7a4214aSPierre Jolivet } 98c7a4214aSPierre Jolivet 99c7a4214aSPierre Jolivet static PetscErrorCode MatIncreaseOverlap_Htool(Mat A,PetscInt is_max,IS is[],PetscInt ov) 100c7a4214aSPierre Jolivet { 101c7a4214aSPierre Jolivet std::set<PetscInt> set; 102c7a4214aSPierre Jolivet const PetscInt *idx; 103c7a4214aSPierre Jolivet PetscInt *oidx,size; 104c7a4214aSPierre Jolivet PetscMPIInt csize; 105c7a4214aSPierre Jolivet PetscErrorCode ierr; 106c7a4214aSPierre Jolivet 107c7a4214aSPierre Jolivet PetscFunctionBegin; 108c7a4214aSPierre Jolivet for (PetscInt i=0; i<is_max; ++i) { 109c7a4214aSPierre Jolivet /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */ 110c7a4214aSPierre Jolivet /* needed to avoid subdomain matrices to replicate A since it is dense */ 111c7a4214aSPierre Jolivet ierr = MPI_Comm_size(PetscObjectComm((PetscObject)is[i]),&csize);CHKERRMPI(ierr); 112c7a4214aSPierre Jolivet if (csize != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported parallel IS"); 113c7a4214aSPierre Jolivet ierr = ISGetSize(is[i],&size);CHKERRQ(ierr); 114c7a4214aSPierre Jolivet ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr); 115c7a4214aSPierre Jolivet for (PetscInt j=0; j<size; ++j) { 116c7a4214aSPierre Jolivet set.insert(idx[j]); 117c7a4214aSPierre Jolivet for (PetscInt k=1; k<=ov; ++k) { /* for each layer of overlap */ 118c7a4214aSPierre Jolivet if (idx[j] - k >= 0) set.insert(idx[j] - k); /* do not insert negative indices */ 119c7a4214aSPierre 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 */ 120c7a4214aSPierre Jolivet } 121c7a4214aSPierre Jolivet } 122c7a4214aSPierre Jolivet ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr); 123c7a4214aSPierre Jolivet ierr = ISDestroy(is+i);CHKERRQ(ierr); 124c7a4214aSPierre Jolivet size = set.size(); /* size with overlap */ 125c7a4214aSPierre Jolivet ierr = PetscMalloc1(size,&oidx);CHKERRQ(ierr); 126c7a4214aSPierre Jolivet for (const PetscInt j : set) *oidx++ = j; 127c7a4214aSPierre Jolivet oidx -= size; 128c7a4214aSPierre Jolivet ierr = ISCreateGeneral(PETSC_COMM_SELF,size,oidx,PETSC_OWN_POINTER,is+i);CHKERRQ(ierr); 129c7a4214aSPierre Jolivet } 130c7a4214aSPierre Jolivet PetscFunctionReturn(0); 131c7a4214aSPierre Jolivet } 132c7a4214aSPierre Jolivet 133c7a4214aSPierre Jolivet static PetscErrorCode MatCreateSubMatrices_Htool(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 134c7a4214aSPierre Jolivet { 135c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 136c7a4214aSPierre Jolivet Mat D,B,BT; 137c7a4214aSPierre Jolivet const PetscScalar *copy; 138c7a4214aSPierre Jolivet PetscScalar *ptr; 139c7a4214aSPierre Jolivet const PetscInt *idxr,*idxc,*it; 140c7a4214aSPierre Jolivet PetscInt nrow,m,i; 141c7a4214aSPierre Jolivet PetscBool flg; 142c7a4214aSPierre Jolivet PetscErrorCode ierr; 143c7a4214aSPierre Jolivet 144c7a4214aSPierre Jolivet PetscFunctionBegin; 145c7a4214aSPierre Jolivet if (scall != MAT_REUSE_MATRIX) { 146c7a4214aSPierre Jolivet ierr = PetscCalloc1(n,submat);CHKERRQ(ierr); 147c7a4214aSPierre Jolivet } 148c7a4214aSPierre Jolivet for (i=0; i<n; ++i) { 149c7a4214aSPierre Jolivet ierr = ISGetLocalSize(irow[i],&nrow);CHKERRQ(ierr); 150c7a4214aSPierre Jolivet ierr = ISGetLocalSize(icol[i],&m);CHKERRQ(ierr); 151c7a4214aSPierre Jolivet ierr = ISGetIndices(irow[i],&idxr);CHKERRQ(ierr); 152c7a4214aSPierre Jolivet ierr = ISGetIndices(icol[i],&idxc);CHKERRQ(ierr); 153c7a4214aSPierre Jolivet if (scall != MAT_REUSE_MATRIX) { 154c7a4214aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_SELF,nrow,m,nrow,m,NULL,(*submat)+i);CHKERRQ(ierr); 155c7a4214aSPierre Jolivet } 156c7a4214aSPierre Jolivet ierr = MatDenseGetArrayWrite((*submat)[i],&ptr);CHKERRQ(ierr); 157c7a4214aSPierre Jolivet if (irow[i] == icol[i]) { /* same row and column IS? */ 158c7a4214aSPierre Jolivet ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr); 159c7a4214aSPierre Jolivet if (flg) { 160c7a4214aSPierre Jolivet ierr = ISSorted(irow[i],&flg);CHKERRQ(ierr); 161c7a4214aSPierre Jolivet if (flg) { /* sorted IS? */ 162c7a4214aSPierre Jolivet it = std::lower_bound(idxr,idxr+nrow,A->rmap->rstart); 163c7a4214aSPierre Jolivet if (it != idxr+nrow && *it == A->rmap->rstart) { /* rmap->rstart in IS? */ 164c7a4214aSPierre Jolivet if (std::distance(idxr,it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */ 165c7a4214aSPierre Jolivet for (PetscInt j=0; j<A->rmap->n && flg; ++j) if (PetscUnlikely(it[j] != A->rmap->rstart+j)) flg = PETSC_FALSE; 166c7a4214aSPierre Jolivet if (flg) { /* complete local diagonal block in IS? */ 167c7a4214aSPierre Jolivet /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM 168c7a4214aSPierre Jolivet * [ B C E ] 169c7a4214aSPierre Jolivet * A = [ B D E ] 170c7a4214aSPierre Jolivet * [ B F E ] 171c7a4214aSPierre Jolivet */ 172c7a4214aSPierre Jolivet m = std::distance(idxr,it); /* shift of the coefficient (0,0) of block D from above */ 173c7a4214aSPierre Jolivet ierr = MatGetDiagonalBlock_Htool(A,&D);CHKERRQ(ierr); 174c7a4214aSPierre Jolivet ierr = MatDenseGetArrayRead(D,©);CHKERRQ(ierr); 175c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { 176c7a4214aSPierre Jolivet ierr = PetscArraycpy(ptr+(m+k)*nrow+m,copy+k*A->rmap->n,A->rmap->n);CHKERRQ(ierr); /* block D from above */ 177c7a4214aSPierre Jolivet } 178c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayRead(D,©);CHKERRQ(ierr); 179c7a4214aSPierre Jolivet if (m) { 180c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* vertical block B from above */ 181c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 182c7a4214aSPierre Jolivet if (A->symmetric || A->hermitian) { 183c7a4214aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_SELF,A->rmap->n,m,A->rmap->n,m,ptr+m,&B);CHKERRQ(ierr); 184c7a4214aSPierre Jolivet ierr = MatDenseSetLDA(B,nrow);CHKERRQ(ierr); 185c7a4214aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_SELF,m,A->rmap->n,m,A->rmap->n,ptr+m*nrow,&BT);CHKERRQ(ierr); 186c7a4214aSPierre Jolivet ierr = MatDenseSetLDA(BT,nrow);CHKERRQ(ierr); 187c7a4214aSPierre Jolivet if (A->hermitian && PetscDefined(USE_COMPLEX)) { 188c7a4214aSPierre Jolivet ierr = MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr); 189c7a4214aSPierre Jolivet } else { 190c7a4214aSPierre Jolivet ierr = MatTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr); 191c7a4214aSPierre Jolivet } 192c7a4214aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 193c7a4214aSPierre Jolivet ierr = MatDestroy(&BT);CHKERRQ(ierr); 194c7a4214aSPierre Jolivet } else { 195c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { /* block C from above */ 196c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(m,1,idxr,idxc+m+k,ptr+(m+k)*nrow); 197c7a4214aSPierre Jolivet } 198c7a4214aSPierre Jolivet } 199c7a4214aSPierre Jolivet } 200c7a4214aSPierre Jolivet if (m+A->rmap->n != nrow) { 201c7a4214aSPierre 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 */ 202c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 203c7a4214aSPierre Jolivet if (A->symmetric || A->hermitian) { 204c7a4214aSPierre 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); 205c7a4214aSPierre Jolivet ierr = MatDenseSetLDA(B,nrow);CHKERRQ(ierr); 206c7a4214aSPierre 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); 207c7a4214aSPierre Jolivet ierr = MatDenseSetLDA(BT,nrow);CHKERRQ(ierr); 208c7a4214aSPierre Jolivet if (A->hermitian && PetscDefined(USE_COMPLEX)) { 209c7a4214aSPierre Jolivet ierr = MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr); 210c7a4214aSPierre Jolivet } else { 211c7a4214aSPierre Jolivet ierr = MatTranspose(B,MAT_REUSE_MATRIX,&BT);CHKERRQ(ierr); 212c7a4214aSPierre Jolivet } 213c7a4214aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 214c7a4214aSPierre Jolivet ierr = MatDestroy(&BT);CHKERRQ(ierr); 215c7a4214aSPierre Jolivet } else { 216c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { /* block F from above */ 217c7a4214aSPierre 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); 218c7a4214aSPierre Jolivet } 219c7a4214aSPierre Jolivet } 220c7a4214aSPierre Jolivet } 221c7a4214aSPierre Jolivet } /* complete local diagonal block not in IS */ 222c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */ 223c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* rmap->rstart not in IS */ 224c7a4214aSPierre Jolivet } /* unsorted IS */ 225c7a4214aSPierre Jolivet } 226c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* different row and column IS */ 227c7a4214aSPierre Jolivet if (!flg) a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* reassemble everything */ 228c7a4214aSPierre Jolivet ierr = ISRestoreIndices(irow[i],&idxr);CHKERRQ(ierr); 229c7a4214aSPierre Jolivet ierr = ISRestoreIndices(icol[i],&idxc);CHKERRQ(ierr); 230c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayWrite((*submat)[i],&ptr);CHKERRQ(ierr); 231c7a4214aSPierre Jolivet ierr = MatAssemblyBegin((*submat)[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 232c7a4214aSPierre Jolivet ierr = MatAssemblyEnd((*submat)[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 233c7a4214aSPierre Jolivet ierr = MatScale((*submat)[i],a->s);CHKERRQ(ierr); 234c7a4214aSPierre Jolivet } 235c7a4214aSPierre Jolivet PetscFunctionReturn(0); 236c7a4214aSPierre Jolivet } 237c7a4214aSPierre Jolivet 238c7a4214aSPierre Jolivet static PetscErrorCode MatDestroy_Htool(Mat A) 239c7a4214aSPierre Jolivet { 240c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 241c7a4214aSPierre Jolivet PetscContainer container; 242c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 243c7a4214aSPierre Jolivet PetscErrorCode ierr; 244c7a4214aSPierre Jolivet 245c7a4214aSPierre Jolivet PetscFunctionBegin; 246c7a4214aSPierre Jolivet ierr = PetscObjectChangeTypeName((PetscObject)A,NULL);CHKERRQ(ierr); 247c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",NULL);CHKERRQ(ierr); 248c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",NULL);CHKERRQ(ierr); 249c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",NULL);CHKERRQ(ierr); 250c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",NULL);CHKERRQ(ierr); 251c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",NULL);CHKERRQ(ierr); 252c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",NULL);CHKERRQ(ierr); 253c7a4214aSPierre Jolivet ierr = PetscObjectQuery((PetscObject)A,"KernelTranspose",(PetscObject*)&container);CHKERRQ(ierr); 254c7a4214aSPierre Jolivet if (container) { /* created in MatTranspose_Htool() */ 255c7a4214aSPierre Jolivet ierr = PetscContainerGetPointer(container,(void**)&kernelt);CHKERRQ(ierr); 256c7a4214aSPierre Jolivet ierr = MatDestroy(&kernelt->A);CHKERRQ(ierr); 257c7a4214aSPierre Jolivet ierr = PetscFree(kernelt);CHKERRQ(ierr); 258c7a4214aSPierre Jolivet ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 259c7a4214aSPierre Jolivet ierr = PetscObjectCompose((PetscObject)A,"KernelTranspose",NULL);CHKERRQ(ierr); 260c7a4214aSPierre Jolivet } 261c7a4214aSPierre Jolivet if (a->gcoords_source != a->gcoords_target) { 262c7a4214aSPierre Jolivet ierr = PetscFree(a->gcoords_source);CHKERRQ(ierr); 263c7a4214aSPierre Jolivet } 264c7a4214aSPierre Jolivet ierr = PetscFree(a->gcoords_target);CHKERRQ(ierr); 265c7a4214aSPierre Jolivet ierr = PetscFree2(a->work_source,a->work_target);CHKERRQ(ierr); 266c7a4214aSPierre Jolivet delete a->wrapper; 267c7a4214aSPierre Jolivet delete a->hmatrix; 268c7a4214aSPierre Jolivet ierr = PetscFree(A->data);CHKERRQ(ierr); 269c7a4214aSPierre Jolivet PetscFunctionReturn(0); 270c7a4214aSPierre Jolivet } 271c7a4214aSPierre Jolivet 272c7a4214aSPierre Jolivet static PetscErrorCode MatView_Htool(Mat A,PetscViewer pv) 273c7a4214aSPierre Jolivet { 274c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 275c7a4214aSPierre Jolivet PetscBool flg; 276c7a4214aSPierre Jolivet PetscErrorCode ierr; 277c7a4214aSPierre Jolivet 278c7a4214aSPierre Jolivet PetscFunctionBegin; 279c7a4214aSPierre Jolivet ierr = PetscObjectTypeCompare((PetscObject)pv,PETSCVIEWERASCII,&flg);CHKERRQ(ierr); 280c7a4214aSPierre Jolivet if (flg) { 281c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"symmetry: %c\n",a->hmatrix->get_symmetry_type());CHKERRQ(ierr); 282c7a4214aSPierre Jolivet if (PetscAbsScalar(a->s-1.0) > PETSC_MACHINE_EPSILON) { 283c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX) 284c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"scaling: %g+%gi\n",(double)PetscRealPart(a->s),(double)PetscImaginaryPart(a->s));CHKERRQ(ierr); 285c7a4214aSPierre Jolivet #else 286c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"scaling: %g\n",(double)a->s);CHKERRQ(ierr); 287c7a4214aSPierre Jolivet #endif 288c7a4214aSPierre Jolivet } 289c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"minimum cluster size: %D\n",a->bs[0]);CHKERRQ(ierr); 290c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"maximum block size: %D\n",a->bs[1]);CHKERRQ(ierr); 291c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"epsilon: %g\n",(double)a->epsilon);CHKERRQ(ierr); 292c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"eta: %g\n",(double)a->eta);CHKERRQ(ierr); 293c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"minimum target depth: %D\n",a->depth[0]);CHKERRQ(ierr); 294c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"minimum source depth: %D\n",a->depth[1]);CHKERRQ(ierr); 295c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"compressor: %s\n",MatHtoolCompressorTypes[a->compressor]);CHKERRQ(ierr); 296c7a4214aSPierre Jolivet ierr = PetscViewerASCIIPrintf(pv,"compression: %s\n",a->hmatrix->get_infos("Compression").c_str());CHKERRQ(ierr); 297c7a4214aSPierre 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); 298c7a4214aSPierre 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); 299c7a4214aSPierre 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); 300c7a4214aSPierre 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); 301c7a4214aSPierre Jolivet } 302c7a4214aSPierre Jolivet PetscFunctionReturn(0); 303c7a4214aSPierre Jolivet } 304c7a4214aSPierre Jolivet 305c7a4214aSPierre Jolivet static PetscErrorCode MatScale_Htool(Mat A,PetscScalar s) 306c7a4214aSPierre Jolivet { 307c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 308c7a4214aSPierre Jolivet 309c7a4214aSPierre Jolivet PetscFunctionBegin; 310c7a4214aSPierre Jolivet a->s *= s; 311c7a4214aSPierre Jolivet PetscFunctionReturn(0); 312c7a4214aSPierre Jolivet } 313c7a4214aSPierre Jolivet 314c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */ 315c7a4214aSPierre Jolivet static PetscErrorCode MatGetRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 316c7a4214aSPierre Jolivet { 317c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 318c7a4214aSPierre Jolivet PetscInt *idxc; 319c7a4214aSPierre Jolivet PetscBLASInt one = 1,bn; 320c7a4214aSPierre Jolivet PetscErrorCode ierr; 321c7a4214aSPierre Jolivet 322c7a4214aSPierre Jolivet PetscFunctionBegin; 323c7a4214aSPierre Jolivet if (nz) *nz = A->cmap->N; 324c7a4214aSPierre Jolivet if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */ 325c7a4214aSPierre Jolivet ierr = PetscMalloc1(A->cmap->N,&idxc);CHKERRQ(ierr); 326c7a4214aSPierre Jolivet for (PetscInt i=0; i<A->cmap->N; ++i) idxc[i] = i; 327c7a4214aSPierre Jolivet } 328c7a4214aSPierre Jolivet if (idx) *idx = idxc; 329c7a4214aSPierre Jolivet if (v) { 330c7a4214aSPierre Jolivet ierr = PetscMalloc1(A->cmap->N,v);CHKERRQ(ierr); 331c7a4214aSPierre Jolivet if (a->wrapper) a->wrapper->copy_submatrix(1,A->cmap->N,&row,idxc,*v); 332c7a4214aSPierre Jolivet else reinterpret_cast<htool::IMatrix<PetscScalar>*>(a->kernelctx)->copy_submatrix(1,A->cmap->N,&row,idxc,*v); 333c7a4214aSPierre Jolivet ierr = PetscBLASIntCast(A->cmap->N,&bn);CHKERRQ(ierr); 334c7a4214aSPierre Jolivet PetscStackCallBLAS("BLASscal",BLASscal_(&bn,&a->s,*v,&one)); 335c7a4214aSPierre Jolivet } 336c7a4214aSPierre Jolivet if (!idx) { 337c7a4214aSPierre Jolivet ierr = PetscFree(idxc);CHKERRQ(ierr); 338c7a4214aSPierre Jolivet } 339c7a4214aSPierre Jolivet PetscFunctionReturn(0); 340c7a4214aSPierre Jolivet } 341c7a4214aSPierre Jolivet 342c7a4214aSPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 343c7a4214aSPierre Jolivet { 344c7a4214aSPierre Jolivet PetscErrorCode ierr; 345c7a4214aSPierre Jolivet 346c7a4214aSPierre Jolivet PetscFunctionBegin; 347c7a4214aSPierre Jolivet if (nz) *nz = 0; 348c7a4214aSPierre Jolivet if (idx) { 349c7a4214aSPierre Jolivet ierr = PetscFree(*idx);CHKERRQ(ierr); 350c7a4214aSPierre Jolivet } 351c7a4214aSPierre Jolivet if (v) { 352c7a4214aSPierre Jolivet ierr = PetscFree(*v);CHKERRQ(ierr); 353c7a4214aSPierre Jolivet } 354c7a4214aSPierre Jolivet PetscFunctionReturn(0); 355c7a4214aSPierre Jolivet } 356c7a4214aSPierre Jolivet 357c7a4214aSPierre Jolivet static PetscErrorCode MatSetFromOptions_Htool(PetscOptionItems *PetscOptionsObject,Mat A) 358c7a4214aSPierre Jolivet { 359c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 360c7a4214aSPierre Jolivet PetscInt n; 361c7a4214aSPierre Jolivet PetscBool flg; 362c7a4214aSPierre Jolivet PetscErrorCode ierr; 363c7a4214aSPierre Jolivet 364c7a4214aSPierre Jolivet PetscFunctionBegin; 365c7a4214aSPierre Jolivet ierr = PetscOptionsHead(PetscOptionsObject,"Htool options");CHKERRQ(ierr); 366c7a4214aSPierre Jolivet ierr = PetscOptionsInt("-mat_htool_min_cluster_size","Minimal leaf size in cluster tree",NULL,a->bs[0],a->bs,NULL);CHKERRQ(ierr); 367c7a4214aSPierre 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); 368c7a4214aSPierre Jolivet ierr = PetscOptionsReal("-mat_htool_epsilon","Relative error in Frobenius norm when approximating a block",NULL,a->epsilon,&a->epsilon,NULL);CHKERRQ(ierr); 369c7a4214aSPierre Jolivet ierr = PetscOptionsReal("-mat_htool_eta","Admissibility condition tolerance",NULL,a->eta,&a->eta,NULL);CHKERRQ(ierr); 370c7a4214aSPierre 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); 371c7a4214aSPierre 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); 372c7a4214aSPierre Jolivet n = 0; 373c7a4214aSPierre Jolivet ierr = PetscOptionsEList("-mat_htool_compressor","Type of compression","MatHtoolCompressorType",MatHtoolCompressorTypes,3,MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA],&n,&flg);CHKERRQ(ierr); 374c7a4214aSPierre Jolivet if (flg) a->compressor = MatHtoolCompressorType(n); 375c7a4214aSPierre Jolivet ierr = PetscOptionsTail();CHKERRQ(ierr); 376c7a4214aSPierre Jolivet PetscFunctionReturn(0); 377c7a4214aSPierre Jolivet } 378c7a4214aSPierre Jolivet 379c7a4214aSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A,MatAssemblyType type) 380c7a4214aSPierre Jolivet { 381c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 382c7a4214aSPierre Jolivet const PetscInt *ranges; 383c7a4214aSPierre Jolivet PetscInt *offset; 384c7a4214aSPierre Jolivet PetscMPIInt size; 385c7a4214aSPierre Jolivet char S = PetscDefined(USE_COMPLEX) && A->hermitian ? 'H' : (A->symmetric ? 'S' : 'N'),uplo = S == 'N' ? 'N' : 'U'; 386c7a4214aSPierre Jolivet htool::IMatrix<PetscScalar> *generator = nullptr; 387c7a4214aSPierre Jolivet std::shared_ptr<htool::RegularClustering> t,s = nullptr; 388c7a4214aSPierre Jolivet PetscErrorCode ierr; 389c7a4214aSPierre Jolivet 390c7a4214aSPierre Jolivet PetscFunctionBegin; 391c7a4214aSPierre Jolivet ierr = PetscCitationsRegister(HtoolCitation,&HtoolCite);CHKERRQ(ierr); 392c7a4214aSPierre Jolivet delete a->wrapper; 393c7a4214aSPierre Jolivet delete a->hmatrix; 394c7a4214aSPierre Jolivet ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr); 395c7a4214aSPierre Jolivet ierr = PetscMalloc1(2*size,&offset);CHKERRQ(ierr); 396c7a4214aSPierre Jolivet ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 397c7a4214aSPierre Jolivet for (PetscInt i=0; i<size; ++i) { 398c7a4214aSPierre Jolivet offset[2*i] = ranges[i]; 399c7a4214aSPierre Jolivet offset[2*i+1] = ranges[i+1] - ranges[i]; 400c7a4214aSPierre Jolivet } 401c7a4214aSPierre Jolivet t = std::make_shared<htool::RegularClustering>(a->dim); 402c7a4214aSPierre Jolivet t->set_minclustersize(a->bs[0]); 403c7a4214aSPierre Jolivet t->build_local_auto(A->rmap->N,a->gcoords_target,offset); 404c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx); 405c7a4214aSPierre Jolivet else { 406c7a4214aSPierre Jolivet a->wrapper = NULL; 407c7a4214aSPierre Jolivet generator = reinterpret_cast<htool::IMatrix<PetscScalar>*>(a->kernelctx); 408c7a4214aSPierre Jolivet } 409c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 410c7a4214aSPierre Jolivet ierr = MatGetOwnershipRangesColumn(A,&ranges);CHKERRQ(ierr); 411c7a4214aSPierre Jolivet for (PetscInt i=0; i<size; ++i) { 412c7a4214aSPierre Jolivet offset[2*i] = ranges[i]; 413c7a4214aSPierre Jolivet offset[2*i+1] = ranges[i+1] - ranges[i]; 414c7a4214aSPierre Jolivet } 415c7a4214aSPierre Jolivet s = std::make_shared<htool::RegularClustering>(a->dim); 416c7a4214aSPierre Jolivet s->set_minclustersize(a->bs[0]); 417c7a4214aSPierre Jolivet s->build_local_auto(A->cmap->N,a->gcoords_source,offset); 418c7a4214aSPierre Jolivet S = uplo = 'N'; 419c7a4214aSPierre Jolivet } 420c7a4214aSPierre Jolivet ierr = PetscFree(offset);CHKERRQ(ierr); 421c7a4214aSPierre Jolivet switch (a->compressor) { 422c7a4214aSPierre Jolivet case MAT_HTOOL_COMPRESSOR_FULL_ACA: 423c7a4214aSPierre Jolivet a->hmatrix = dynamic_cast<htool::HMatrixVirtual<PetscScalar>*>(new htool::HMatrix<PetscScalar,htool::fullACA,htool::RegularClustering,htool::RjasanowSteinbach>(t,s?s:t,a->epsilon,a->eta,S,uplo)); 424c7a4214aSPierre Jolivet break; 425c7a4214aSPierre Jolivet case MAT_HTOOL_COMPRESSOR_SVD: 426c7a4214aSPierre Jolivet a->hmatrix = dynamic_cast<htool::HMatrixVirtual<PetscScalar>*>(new htool::HMatrix<PetscScalar,htool::SVD,htool::RegularClustering,htool::RjasanowSteinbach>(t,s?s:t,a->epsilon,a->eta,S,uplo)); 427c7a4214aSPierre Jolivet break; 428c7a4214aSPierre Jolivet default: 429c7a4214aSPierre Jolivet a->hmatrix = dynamic_cast<htool::HMatrixVirtual<PetscScalar>*>(new htool::HMatrix<PetscScalar,htool::sympartialACA,htool::RegularClustering,htool::RjasanowSteinbach>(t,s?s:t,a->epsilon,a->eta,S,uplo)); 430c7a4214aSPierre Jolivet } 431c7a4214aSPierre Jolivet a->hmatrix->set_maxblocksize(a->bs[1]); 432c7a4214aSPierre Jolivet a->hmatrix->set_mintargetdepth(a->depth[0]); 433c7a4214aSPierre Jolivet a->hmatrix->set_minsourcedepth(a->depth[1]); 434c7a4214aSPierre Jolivet if (s) a->hmatrix->build_auto(a->wrapper ? *a->wrapper : *generator,a->gcoords_target,a->gcoords_source); 435c7a4214aSPierre Jolivet else a->hmatrix->build_auto_sym(a->wrapper ? *a->wrapper : *generator,a->gcoords_target); 436c7a4214aSPierre Jolivet PetscFunctionReturn(0); 437c7a4214aSPierre Jolivet } 438c7a4214aSPierre Jolivet 439c7a4214aSPierre Jolivet static PetscErrorCode MatProductNumeric_Htool(Mat C) 440c7a4214aSPierre Jolivet { 441c7a4214aSPierre Jolivet Mat_Product *product = C->product; 442c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)product->A->data; 443c7a4214aSPierre Jolivet const PetscScalar *in; 444c7a4214aSPierre Jolivet PetscScalar *out; 445c7a4214aSPierre Jolivet PetscInt lda; 446c7a4214aSPierre Jolivet PetscErrorCode ierr; 447c7a4214aSPierre Jolivet 448c7a4214aSPierre Jolivet PetscFunctionBegin; 449c7a4214aSPierre Jolivet MatCheckProduct(C,1); 450c7a4214aSPierre Jolivet switch (product->type) { 451c7a4214aSPierre Jolivet case MATPRODUCT_AB: 452c7a4214aSPierre Jolivet PetscInt N; 453c7a4214aSPierre Jolivet ierr = MatGetSize(C,NULL,&N);CHKERRQ(ierr); 454c7a4214aSPierre Jolivet ierr = MatDenseGetLDA(C,&lda);CHKERRQ(ierr); 455c7a4214aSPierre Jolivet if (lda != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%D != %D)",lda,C->rmap->n); 456c7a4214aSPierre Jolivet ierr = MatDenseGetArrayRead(product->B,&in);CHKERRQ(ierr); 457c7a4214aSPierre Jolivet ierr = MatDenseGetArrayWrite(C,&out);CHKERRQ(ierr); 458c7a4214aSPierre Jolivet a->hmatrix->mvprod_local_to_local(in,out,N); 459c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayWrite(C,&out);CHKERRQ(ierr); 460c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayRead(product->B,&in);CHKERRQ(ierr); 461c7a4214aSPierre Jolivet ierr = MatScale(C,a->s);CHKERRQ(ierr); 462c7a4214aSPierre Jolivet break; 463c7a4214aSPierre Jolivet default: 464c7a4214aSPierre Jolivet SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType %s is not supported",MatProductTypes[product->type]); 465c7a4214aSPierre Jolivet } 466c7a4214aSPierre Jolivet PetscFunctionReturn(0); 467c7a4214aSPierre Jolivet } 468c7a4214aSPierre Jolivet 469c7a4214aSPierre Jolivet static PetscErrorCode MatProductSymbolic_Htool(Mat C) 470c7a4214aSPierre Jolivet { 471c7a4214aSPierre Jolivet Mat_Product *product = C->product; 472c7a4214aSPierre Jolivet Mat A,B; 473c7a4214aSPierre Jolivet PetscBool flg; 474c7a4214aSPierre Jolivet PetscErrorCode ierr; 475c7a4214aSPierre Jolivet 476c7a4214aSPierre Jolivet PetscFunctionBegin; 477c7a4214aSPierre Jolivet MatCheckProduct(C,1); 478c7a4214aSPierre Jolivet A = product->A; 479c7a4214aSPierre Jolivet B = product->B; 480c7a4214aSPierre Jolivet ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 481c7a4214aSPierre Jolivet if (!flg) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MatProduct_AB not supported for %s",((PetscObject)product->B)->type_name); 482c7a4214aSPierre Jolivet switch (product->type) { 483c7a4214aSPierre Jolivet case MATPRODUCT_AB: 484c7a4214aSPierre Jolivet if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) { 485c7a4214aSPierre Jolivet ierr = MatSetSizes(C,A->rmap->n,B->cmap->n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 486c7a4214aSPierre Jolivet } 487c7a4214aSPierre Jolivet ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr); 488c7a4214aSPierre Jolivet ierr = MatSetUp(C);CHKERRQ(ierr); 489c7a4214aSPierre Jolivet ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 490c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 491c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 492c7a4214aSPierre Jolivet break; 493c7a4214aSPierre Jolivet default: 494c7a4214aSPierre Jolivet SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"ProductType %s is not supported",MatProductTypes[product->type]); 495c7a4214aSPierre Jolivet } 496c7a4214aSPierre Jolivet C->ops->productsymbolic = NULL; 497c7a4214aSPierre Jolivet C->ops->productnumeric = MatProductNumeric_Htool; 498c7a4214aSPierre Jolivet PetscFunctionReturn(0); 499c7a4214aSPierre Jolivet } 500c7a4214aSPierre Jolivet 501c7a4214aSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Htool(Mat C) 502c7a4214aSPierre Jolivet { 503c7a4214aSPierre Jolivet PetscFunctionBegin; 504c7a4214aSPierre Jolivet MatCheckProduct(C,1); 505c7a4214aSPierre Jolivet if (C->product->type == MATPRODUCT_AB) C->ops->productsymbolic = MatProductSymbolic_Htool; 506c7a4214aSPierre Jolivet PetscFunctionReturn(0); 507c7a4214aSPierre Jolivet } 508c7a4214aSPierre Jolivet 509c7a4214aSPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A,const htool::HMatrixVirtual<PetscScalar> **hmatrix) 510c7a4214aSPierre Jolivet { 511c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 512c7a4214aSPierre Jolivet 513c7a4214aSPierre Jolivet PetscFunctionBegin; 514c7a4214aSPierre Jolivet *hmatrix = a->hmatrix; 515c7a4214aSPierre Jolivet PetscFunctionReturn(0); 516c7a4214aSPierre Jolivet } 517c7a4214aSPierre Jolivet 518c7a4214aSPierre Jolivet /*@C 519c7a4214aSPierre Jolivet MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a MATHTOOL. 520c7a4214aSPierre Jolivet 521c7a4214aSPierre Jolivet Input Parameter: 522c7a4214aSPierre Jolivet . A - hierarchical matrix 523c7a4214aSPierre Jolivet 524c7a4214aSPierre Jolivet Output Parameter: 525c7a4214aSPierre Jolivet . hmatrix - opaque pointer to a Htool virtual matrix 526c7a4214aSPierre Jolivet 527c7a4214aSPierre Jolivet Level: advanced 528c7a4214aSPierre Jolivet 529c7a4214aSPierre Jolivet .seealso: MATHTOOL 530c7a4214aSPierre Jolivet @*/ 531c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A,const htool::HMatrixVirtual<PetscScalar> **hmatrix) 532c7a4214aSPierre Jolivet { 533c7a4214aSPierre Jolivet PetscErrorCode ierr; 534c7a4214aSPierre Jolivet 535c7a4214aSPierre Jolivet PetscFunctionBegin; 536c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 537c7a4214aSPierre Jolivet PetscValidPointer(hmatrix,2); 538c7a4214aSPierre Jolivet ierr = PetscTryMethod(A,"MatHtoolGetHierarchicalMat_C",(Mat,const htool::HMatrixVirtual<PetscScalar>**),(A,hmatrix));CHKERRQ(ierr); 539c7a4214aSPierre Jolivet PetscFunctionReturn(0); 540c7a4214aSPierre Jolivet } 541c7a4214aSPierre Jolivet 542c7a4214aSPierre Jolivet static PetscErrorCode MatHtoolSetKernel_Htool(Mat A,MatHtoolKernel kernel,void *kernelctx) 543c7a4214aSPierre Jolivet { 544c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 545c7a4214aSPierre Jolivet 546c7a4214aSPierre Jolivet PetscFunctionBegin; 547c7a4214aSPierre Jolivet a->kernel = kernel; 548c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 549c7a4214aSPierre Jolivet delete a->wrapper; 550c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx); 551c7a4214aSPierre Jolivet PetscFunctionReturn(0); 552c7a4214aSPierre Jolivet } 553c7a4214aSPierre Jolivet 554c7a4214aSPierre Jolivet /*@C 555c7a4214aSPierre Jolivet MatHtoolSetKernel - Sets the kernel and context used for the assembly of a MATHTOOL 556c7a4214aSPierre Jolivet 557c7a4214aSPierre Jolivet Input Parameters: 558c7a4214aSPierre Jolivet + A - hierarchical matrix 559c7a4214aSPierre Jolivet . kernel - computational kernel (or NULL) 560c7a4214aSPierre Jolivet - kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::IMatrix<PetscScalar>*) 561c7a4214aSPierre Jolivet 562c7a4214aSPierre Jolivet Level: advanced 563c7a4214aSPierre Jolivet 564c7a4214aSPierre Jolivet .seealso: MATHTOOL, MatCreateHtoolFromKernel() 565c7a4214aSPierre Jolivet @*/ 566c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat A,MatHtoolKernel kernel,void *kernelctx) 567c7a4214aSPierre Jolivet { 568c7a4214aSPierre Jolivet PetscErrorCode ierr; 569c7a4214aSPierre Jolivet 570c7a4214aSPierre Jolivet PetscFunctionBegin; 571c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 572c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel,2); 573c7a4214aSPierre Jolivet if (!kernel) PetscValidPointer(kernelctx,3); 574c7a4214aSPierre Jolivet ierr = PetscTryMethod(A,"MatHtoolSetKernel_C",(Mat,MatHtoolKernel,void*),(A,kernel,kernelctx));CHKERRQ(ierr); 575c7a4214aSPierre Jolivet PetscFunctionReturn(0); 576c7a4214aSPierre Jolivet } 577c7a4214aSPierre Jolivet 578c7a4214aSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B) 579c7a4214aSPierre Jolivet { 580c7a4214aSPierre Jolivet Mat C; 581c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 582c7a4214aSPierre Jolivet PetscInt lda; 583c7a4214aSPierre Jolivet PetscScalar *array; 584c7a4214aSPierre Jolivet PetscErrorCode ierr; 585c7a4214aSPierre Jolivet 586c7a4214aSPierre Jolivet PetscFunctionBegin; 587c7a4214aSPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 588c7a4214aSPierre Jolivet C = *B; 589c7a4214aSPierre Jolivet if (C->rmap->n != A->rmap->n || C->cmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible dimensions"); 590c7a4214aSPierre Jolivet ierr = MatDenseGetLDA(C,&lda);CHKERRQ(ierr); 591c7a4214aSPierre Jolivet if (lda != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%D != %D)",lda,C->rmap->n); 592c7a4214aSPierre Jolivet } else { 593c7a4214aSPierre Jolivet ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 594c7a4214aSPierre Jolivet ierr = MatSetSizes(C,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 595c7a4214aSPierre Jolivet ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr); 596c7a4214aSPierre Jolivet ierr = MatSetUp(C);CHKERRQ(ierr); 597c7a4214aSPierre Jolivet ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 598c7a4214aSPierre Jolivet } 599c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 600c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 601c7a4214aSPierre Jolivet ierr = MatDenseGetArrayWrite(C,&array);CHKERRQ(ierr); 602c7a4214aSPierre Jolivet a->hmatrix->copy_local_dense_perm(array); 603c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayWrite(C,&array);CHKERRQ(ierr); 604c7a4214aSPierre Jolivet ierr = MatScale(C,a->s);CHKERRQ(ierr); 605c7a4214aSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 606c7a4214aSPierre Jolivet ierr = MatHeaderReplace(A,&C);CHKERRQ(ierr); 607c7a4214aSPierre Jolivet } else *B = C; 608c7a4214aSPierre Jolivet PetscFunctionReturn(0); 609c7a4214aSPierre Jolivet } 610c7a4214aSPierre Jolivet 611c7a4214aSPierre Jolivet static PetscScalar GenEntryTranspose(PetscInt sdim,PetscInt i,PetscInt j,void *ctx) 612c7a4214aSPierre Jolivet { 613c7a4214aSPierre Jolivet MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose*)ctx; 614c7a4214aSPierre Jolivet return generator->kernel(sdim,j,i,generator->kernelctx); 615c7a4214aSPierre Jolivet } 616c7a4214aSPierre Jolivet 617c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */ 618c7a4214aSPierre Jolivet static PetscErrorCode MatTranspose_Htool(Mat A,MatReuse reuse,Mat *B) 619c7a4214aSPierre Jolivet { 620c7a4214aSPierre Jolivet Mat C; 621c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data,*c; 622c7a4214aSPierre Jolivet PetscInt M = A->rmap->N,N = A->cmap->N,m = A->rmap->n,n = A->cmap->n; 623c7a4214aSPierre Jolivet PetscContainer container; 624c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 625c7a4214aSPierre Jolivet PetscErrorCode ierr; 626c7a4214aSPierre Jolivet 627c7a4214aSPierre Jolivet PetscFunctionBegin; 628c7a4214aSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatTranspose() with MAT_INPLACE_MATRIX not supported"); 629c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 630c7a4214aSPierre Jolivet ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 631c7a4214aSPierre Jolivet ierr = MatSetSizes(C,n,m,N,M);CHKERRQ(ierr); 632c7a4214aSPierre Jolivet ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr); 633c7a4214aSPierre Jolivet ierr = MatSetUp(C);CHKERRQ(ierr); 634c7a4214aSPierre Jolivet ierr = PetscContainerCreate(PetscObjectComm((PetscObject)C),&container);CHKERRQ(ierr); 635c7a4214aSPierre Jolivet ierr = PetscNew(&kernelt);CHKERRQ(ierr); 636c7a4214aSPierre Jolivet ierr = PetscContainerSetPointer(container,kernelt);CHKERRQ(ierr); 637c7a4214aSPierre Jolivet ierr = PetscObjectCompose((PetscObject)C,"KernelTranspose",(PetscObject)container);CHKERRQ(ierr); 638c7a4214aSPierre Jolivet } else { 639c7a4214aSPierre Jolivet C = *B; 640c7a4214aSPierre Jolivet ierr = PetscObjectQuery((PetscObject)C,"KernelTranspose",(PetscObject*)&container);CHKERRQ(ierr); 641c7a4214aSPierre Jolivet if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatTranspose() with MAT_INITIAL_MATRIX first"); 642c7a4214aSPierre Jolivet ierr = PetscContainerGetPointer(container,(void**)&kernelt);CHKERRQ(ierr); 643c7a4214aSPierre Jolivet } 644c7a4214aSPierre Jolivet c = (Mat_Htool*)C->data; 645c7a4214aSPierre Jolivet c->dim = a->dim; 646c7a4214aSPierre Jolivet c->s = a->s; 647c7a4214aSPierre Jolivet c->kernel = GenEntryTranspose; 648c7a4214aSPierre Jolivet if (kernelt->A != A) { 649c7a4214aSPierre Jolivet ierr = MatDestroy(&kernelt->A);CHKERRQ(ierr); 650c7a4214aSPierre Jolivet kernelt->A = A; 651c7a4214aSPierre Jolivet ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 652c7a4214aSPierre Jolivet } 653c7a4214aSPierre Jolivet kernelt->kernel = a->kernel; 654c7a4214aSPierre Jolivet kernelt->kernelctx = a->kernelctx; 655c7a4214aSPierre Jolivet c->kernelctx = kernelt; 656c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 657c7a4214aSPierre Jolivet ierr = PetscMalloc1(N*c->dim,&c->gcoords_target);CHKERRQ(ierr); 658c7a4214aSPierre Jolivet ierr = PetscArraycpy(c->gcoords_target,a->gcoords_source,N*c->dim);CHKERRQ(ierr); 659c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 660c7a4214aSPierre Jolivet ierr = PetscMalloc1(M*c->dim,&c->gcoords_source);CHKERRQ(ierr); 661c7a4214aSPierre Jolivet ierr = PetscArraycpy(c->gcoords_source,a->gcoords_target,M*c->dim);CHKERRQ(ierr); 662c7a4214aSPierre Jolivet } else c->gcoords_source = c->gcoords_target; 663c7a4214aSPierre Jolivet ierr = PetscCalloc2(M,&c->work_source,N,&c->work_target);CHKERRQ(ierr); 664c7a4214aSPierre Jolivet } 665c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 666c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 667c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) *B = C; 668c7a4214aSPierre Jolivet PetscFunctionReturn(0); 669c7a4214aSPierre Jolivet } 670c7a4214aSPierre Jolivet 671c7a4214aSPierre Jolivet /*@C 672c7a4214aSPierre Jolivet MatCreateHtoolFromKernel - Creates a MATHTOOL from a user-supplied kernel. 673c7a4214aSPierre Jolivet 674c7a4214aSPierre Jolivet Input Parameters: 675c7a4214aSPierre Jolivet + comm - MPI communicator 676c7a4214aSPierre Jolivet . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 677c7a4214aSPierre Jolivet . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 678c7a4214aSPierre Jolivet . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 679c7a4214aSPierre Jolivet . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 680c7a4214aSPierre Jolivet . spacedim - dimension of the space coordinates 681c7a4214aSPierre Jolivet . coords_target - coordinates of the target 682c7a4214aSPierre Jolivet . coords_source - coordinates of the source 683c7a4214aSPierre Jolivet . kernel - computational kernel (or NULL) 684c7a4214aSPierre Jolivet - kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::IMatrix<PetscScalar>*) 685c7a4214aSPierre Jolivet 686c7a4214aSPierre Jolivet Output Parameter: 687c7a4214aSPierre Jolivet . B - matrix 688c7a4214aSPierre Jolivet 689c7a4214aSPierre Jolivet Options Database Keys: 690c7a4214aSPierre Jolivet + -mat_htool_min_cluster_size <PetscInt> - minimal leaf size in cluster tree 691c7a4214aSPierre Jolivet . -mat_htool_max_block_size <PetscInt> - maximal number of coefficients in a dense block 692c7a4214aSPierre Jolivet . -mat_htool_epsilon <PetscReal> - relative error in Frobenius norm when approximating a block 693c7a4214aSPierre Jolivet . -mat_htool_eta <PetscReal> - admissibility condition tolerance 694c7a4214aSPierre Jolivet . -mat_htool_min_target_depth <PetscInt> - minimal cluster tree depth associated with the rows 695c7a4214aSPierre Jolivet . -mat_htool_min_source_depth <PetscInt> - minimal cluster tree depth associated with the columns 696c7a4214aSPierre Jolivet - -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression 697c7a4214aSPierre Jolivet 698c7a4214aSPierre Jolivet Level: intermediate 699c7a4214aSPierre Jolivet 700c7a4214aSPierre Jolivet .seealso: MatCreate(), MATHTOOL, PCSetCoordinates(), MatHtoolSetKernel(), MatHtoolCompressorType, MATHARA, MatCreateHaraFromKernel() 701c7a4214aSPierre Jolivet @*/ 702c7a4214aSPierre 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) 703c7a4214aSPierre Jolivet { 704c7a4214aSPierre Jolivet Mat A; 705c7a4214aSPierre Jolivet Mat_Htool *a; 706c7a4214aSPierre Jolivet PetscErrorCode ierr; 707c7a4214aSPierre Jolivet 708c7a4214aSPierre Jolivet PetscFunctionBegin; 709c7a4214aSPierre Jolivet ierr = MatCreate(comm,&A);CHKERRQ(ierr); 710c7a4214aSPierre Jolivet PetscValidLogicalCollectiveInt(A,spacedim,6); 711c7a4214aSPierre Jolivet PetscValidRealPointer(coords_target,7); 712c7a4214aSPierre Jolivet PetscValidRealPointer(coords_source,8); 713c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel,9); 714c7a4214aSPierre Jolivet if (!kernel) PetscValidPointer(kernelctx,10); 715c7a4214aSPierre Jolivet ierr = MatSetSizes(A,m,n,M,N);CHKERRQ(ierr); 716c7a4214aSPierre Jolivet ierr = MatSetType(A,MATHTOOL);CHKERRQ(ierr); 717c7a4214aSPierre Jolivet ierr = MatSetUp(A);CHKERRQ(ierr); 718c7a4214aSPierre Jolivet a = (Mat_Htool*)A->data; 719c7a4214aSPierre Jolivet a->dim = spacedim; 720c7a4214aSPierre Jolivet a->s = 1.0; 721c7a4214aSPierre Jolivet a->kernel = kernel; 722c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 723c7a4214aSPierre Jolivet ierr = PetscCalloc1(A->rmap->N*spacedim,&a->gcoords_target);CHKERRQ(ierr); 724c7a4214aSPierre Jolivet ierr = PetscArraycpy(a->gcoords_target+A->rmap->rstart*spacedim,coords_target,A->rmap->n*spacedim);CHKERRQ(ierr); 725c7a4214aSPierre 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 */ 726c7a4214aSPierre Jolivet if (coords_target != coords_source) { 727c7a4214aSPierre Jolivet ierr = PetscCalloc1(A->cmap->N*spacedim,&a->gcoords_source);CHKERRQ(ierr); 728c7a4214aSPierre Jolivet ierr = PetscArraycpy(a->gcoords_source+A->cmap->rstart*spacedim,coords_source,A->cmap->n*spacedim);CHKERRQ(ierr); 729c7a4214aSPierre 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 */ 730c7a4214aSPierre Jolivet } else a->gcoords_source = a->gcoords_target; 731c7a4214aSPierre Jolivet ierr = PetscCalloc2(A->cmap->N,&a->work_source,A->rmap->N,&a->work_target);CHKERRQ(ierr); 732c7a4214aSPierre Jolivet *B = A; 733c7a4214aSPierre Jolivet PetscFunctionReturn(0); 734c7a4214aSPierre Jolivet } 735c7a4214aSPierre Jolivet 736c7a4214aSPierre Jolivet /*MC 737c7a4214aSPierre Jolivet MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package. 738c7a4214aSPierre Jolivet 739c7a4214aSPierre Jolivet Use ./configure --download-htool to install PETSc to use Htool. 740c7a4214aSPierre Jolivet 741c7a4214aSPierre Jolivet Options Database Keys: 742*e85a1b63SPatrick Sanan . -mat_type htool - matrix type to "htool" during a call to MatSetFromOptions() 743c7a4214aSPierre Jolivet 744c7a4214aSPierre Jolivet Level: beginner 745c7a4214aSPierre Jolivet 746c7a4214aSPierre Jolivet .seealso: MATHARA, MATDENSE, MatCreateHtoolFromKernel(), MatHtoolSetKernel() 747c7a4214aSPierre Jolivet M*/ 748c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A) 749c7a4214aSPierre Jolivet { 750c7a4214aSPierre Jolivet Mat_Htool *a; 751c7a4214aSPierre Jolivet PetscErrorCode ierr; 752c7a4214aSPierre Jolivet 753c7a4214aSPierre Jolivet PetscFunctionBegin; 754c7a4214aSPierre Jolivet ierr = PetscNewLog(A,&a);CHKERRQ(ierr); 755c7a4214aSPierre Jolivet A->data = (void*)a; 756c7a4214aSPierre Jolivet ierr = PetscObjectChangeTypeName((PetscObject)A,MATHTOOL);CHKERRQ(ierr); 757c7a4214aSPierre Jolivet ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 758c7a4214aSPierre Jolivet A->ops->getdiagonal = MatGetDiagonal_Htool; 759c7a4214aSPierre Jolivet A->ops->getdiagonalblock = MatGetDiagonalBlock_Htool; 760c7a4214aSPierre Jolivet A->ops->mult = MatMult_Htool; 761c7a4214aSPierre Jolivet A->ops->multadd = MatMultAdd_Htool; 762c7a4214aSPierre Jolivet A->ops->increaseoverlap = MatIncreaseOverlap_Htool; 763c7a4214aSPierre Jolivet A->ops->createsubmatrices = MatCreateSubMatrices_Htool; 764c7a4214aSPierre Jolivet A->ops->transpose = MatTranspose_Htool; 765c7a4214aSPierre Jolivet A->ops->destroy = MatDestroy_Htool; 766c7a4214aSPierre Jolivet A->ops->view = MatView_Htool; 767c7a4214aSPierre Jolivet A->ops->setfromoptions = MatSetFromOptions_Htool; 768c7a4214aSPierre Jolivet A->ops->scale = MatScale_Htool; 769c7a4214aSPierre Jolivet A->ops->getrow = MatGetRow_Htool; 770c7a4214aSPierre Jolivet A->ops->restorerow = MatRestoreRow_Htool; 771c7a4214aSPierre Jolivet A->ops->assemblyend = MatAssemblyEnd_Htool; 772c7a4214aSPierre Jolivet a->dim = 0; 773c7a4214aSPierre Jolivet a->gcoords_target = NULL; 774c7a4214aSPierre Jolivet a->gcoords_source = NULL; 775c7a4214aSPierre Jolivet a->s = 1.0; 776c7a4214aSPierre Jolivet a->bs[0] = 10; 777c7a4214aSPierre Jolivet a->bs[1] = 1000000; 778c7a4214aSPierre Jolivet a->epsilon = PetscSqrtReal(PETSC_SMALL); 779c7a4214aSPierre Jolivet a->eta = 10.0; 780c7a4214aSPierre Jolivet a->depth[0] = 0; 781c7a4214aSPierre Jolivet a->depth[1] = 0; 782c7a4214aSPierre Jolivet a->compressor = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA; 783c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",MatProductSetFromOptions_Htool);CHKERRQ(ierr); 784c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",MatProductSetFromOptions_Htool);CHKERRQ(ierr); 785c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",MatConvert_Htool_Dense);CHKERRQ(ierr); 786c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",MatConvert_Htool_Dense);CHKERRQ(ierr); 787c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",MatHtoolGetHierarchicalMat_Htool);CHKERRQ(ierr); 788c7a4214aSPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",MatHtoolSetKernel_Htool);CHKERRQ(ierr); 789c7a4214aSPierre Jolivet PetscFunctionReturn(0); 790c7a4214aSPierre Jolivet } 791