1c7a4214aSPierre Jolivet #include <../src/mat/impls/htool/htool.hpp> /*I "petscmat.h" I*/ 2c7a4214aSPierre Jolivet #include <petscblaslapack.h> 3c7a4214aSPierre Jolivet #include <set> 4c7a4214aSPierre Jolivet 5c7a4214aSPierre Jolivet const char *const MatHtoolCompressorTypes[] = { "sympartialACA", "fullACA", "SVD" }; 698e73e17SPierre Jolivet const char *const MatHtoolClusteringTypes[] = { "PCARegular", "PCAGeometric", "BoundingBox1Regular", "BoundingBox1Geometric" }; 7c7a4214aSPierre Jolivet const char HtoolCitation[] = "@article{marchand2020two,\n" 8c7a4214aSPierre Jolivet " Author = {Marchand, Pierre and Claeys, Xavier and Jolivet, Pierre and Nataf, Fr\\'ed\\'eric and Tournier, Pierre-Henri},\n" 9c7a4214aSPierre Jolivet " Title = {Two-level preconditioning for $h$-version boundary element approximation of hypersingular operator with {GenEO}},\n" 10c7a4214aSPierre Jolivet " Year = {2020},\n" 11c7a4214aSPierre Jolivet " Publisher = {Elsevier},\n" 12c7a4214aSPierre Jolivet " Journal = {Numerische Mathematik},\n" 13c7a4214aSPierre Jolivet " Volume = {146},\n" 14c7a4214aSPierre Jolivet " Pages = {597--628},\n" 15c7a4214aSPierre Jolivet " Url = {https://github.com/htool-ddm/htool}\n" 16c7a4214aSPierre Jolivet "}\n"; 17c7a4214aSPierre Jolivet static PetscBool HtoolCite = PETSC_FALSE; 18c7a4214aSPierre Jolivet 19c7a4214aSPierre Jolivet static PetscErrorCode MatGetDiagonal_Htool(Mat A,Vec v) 20c7a4214aSPierre Jolivet { 21c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 22c7a4214aSPierre Jolivet PetscScalar *x; 23c7a4214aSPierre Jolivet PetscBool flg; 24c7a4214aSPierre Jolivet 25c7a4214aSPierre Jolivet PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A,&flg)); 275f80ce2aSJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only congruent layouts supported"); 289566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(v,&x)); 29c7a4214aSPierre Jolivet a->hmatrix->copy_local_diagonal(x); 309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(v,&x)); 319566063dSJacob Faibussowitsch PetscCall(VecScale(v,a->s)); 32c7a4214aSPierre Jolivet PetscFunctionReturn(0); 33c7a4214aSPierre Jolivet } 34c7a4214aSPierre Jolivet 35c7a4214aSPierre Jolivet static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A,Mat *b) 36c7a4214aSPierre Jolivet { 37c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 38c7a4214aSPierre Jolivet Mat B; 39c7a4214aSPierre Jolivet PetscScalar *ptr; 40c7a4214aSPierre Jolivet PetscBool flg; 41c7a4214aSPierre Jolivet 42c7a4214aSPierre Jolivet PetscFunctionBegin; 439566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A,&flg)); 445f80ce2aSJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only congruent layouts supported"); 459566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B)); /* same logic as in MatGetDiagonalBlock_MPIDense() */ 46c7a4214aSPierre Jolivet if (!B) { 479566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF,A->rmap->n,A->rmap->n,A->rmap->n,A->rmap->n,NULL,&B)); 489566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(B,&ptr)); 49c7a4214aSPierre Jolivet a->hmatrix->copy_local_diagonal_block(ptr); 509566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(B,&ptr)); 519566063dSJacob Faibussowitsch PetscCall(MatPropagateSymmetryOptions(A,B)); 529566063dSJacob Faibussowitsch PetscCall(MatScale(B,a->s)); 539566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B)); 54c7a4214aSPierre Jolivet *b = B; 559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 56c7a4214aSPierre Jolivet } else *b = B; 57c7a4214aSPierre Jolivet PetscFunctionReturn(0); 58c7a4214aSPierre Jolivet } 59c7a4214aSPierre Jolivet 60c7a4214aSPierre Jolivet static PetscErrorCode MatMult_Htool(Mat A,Vec x,Vec y) 61c7a4214aSPierre Jolivet { 62c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 63c7a4214aSPierre Jolivet const PetscScalar *in; 64c7a4214aSPierre Jolivet PetscScalar *out; 65c7a4214aSPierre Jolivet 66c7a4214aSPierre Jolivet PetscFunctionBegin; 679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&in)); 689566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y,&out)); 69c7a4214aSPierre Jolivet a->hmatrix->mvprod_local_to_local(in,out); 709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&in)); 719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y,&out)); 729566063dSJacob Faibussowitsch PetscCall(VecScale(y,a->s)); 73c7a4214aSPierre Jolivet PetscFunctionReturn(0); 74c7a4214aSPierre Jolivet } 75c7a4214aSPierre Jolivet 76c7a4214aSPierre Jolivet /* naive implementation of MatMultAdd() needed for FEM-BEM coupling via MATNEST */ 77c7a4214aSPierre Jolivet static PetscErrorCode MatMultAdd_Htool(Mat A,Vec v1,Vec v2,Vec v3) 78c7a4214aSPierre Jolivet { 79c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 80c7a4214aSPierre Jolivet Vec tmp; 81c7a4214aSPierre Jolivet const PetscScalar scale = a->s; 82c7a4214aSPierre Jolivet 83c7a4214aSPierre Jolivet PetscFunctionBegin; 849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v2,&tmp)); 859566063dSJacob Faibussowitsch PetscCall(VecCopy(v2,v3)); /* no-op in MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]) since VecCopy() checks for x == y */ 86c7a4214aSPierre Jolivet a->s = 1.0; /* set s to 1.0 since VecAXPY() may be used to scale the MatMult() output Vec */ 879566063dSJacob Faibussowitsch PetscCall(MatMult_Htool(A,v1,tmp)); 889566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3,scale,tmp)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp)); 90c7a4214aSPierre Jolivet a->s = scale; /* set s back to its original value */ 91c7a4214aSPierre Jolivet PetscFunctionReturn(0); 92c7a4214aSPierre Jolivet } 93c7a4214aSPierre Jolivet 948b8ddd11SPierre Jolivet static PetscErrorCode MatMultTranspose_Htool(Mat A,Vec x,Vec y) 958b8ddd11SPierre Jolivet { 968b8ddd11SPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 978b8ddd11SPierre Jolivet const PetscScalar *in; 988b8ddd11SPierre Jolivet PetscScalar *out; 998b8ddd11SPierre Jolivet 1008b8ddd11SPierre Jolivet PetscFunctionBegin; 1019566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&in)); 1029566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y,&out)); 1038b8ddd11SPierre Jolivet a->hmatrix->mvprod_transp_local_to_local(in,out); 1049566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&in)); 1059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y,&out)); 1069566063dSJacob Faibussowitsch PetscCall(VecScale(y,a->s)); 1078b8ddd11SPierre Jolivet PetscFunctionReturn(0); 1088b8ddd11SPierre Jolivet } 1098b8ddd11SPierre Jolivet 110c7a4214aSPierre Jolivet static PetscErrorCode MatIncreaseOverlap_Htool(Mat A,PetscInt is_max,IS is[],PetscInt ov) 111c7a4214aSPierre Jolivet { 112c7a4214aSPierre Jolivet std::set<PetscInt> set; 113c7a4214aSPierre Jolivet const PetscInt *idx; 1148f308287SPierre Jolivet PetscInt *oidx,size,bs[2]; 115c7a4214aSPierre Jolivet PetscMPIInt csize; 116c7a4214aSPierre Jolivet 117c7a4214aSPierre Jolivet PetscFunctionBegin; 1189566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A,bs,bs+1)); 1198f308287SPierre Jolivet if (bs[0] != bs[1]) bs[0] = 1; 120c7a4214aSPierre Jolivet for (PetscInt i=0; i<is_max; ++i) { 121c7a4214aSPierre Jolivet /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */ 122c7a4214aSPierre Jolivet /* needed to avoid subdomain matrices to replicate A since it is dense */ 1239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]),&csize)); 1245f80ce2aSJacob Faibussowitsch PetscCheck(csize == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported parallel IS"); 1259566063dSJacob Faibussowitsch PetscCall(ISGetSize(is[i],&size)); 1269566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is[i],&idx)); 127c7a4214aSPierre Jolivet for (PetscInt j=0; j<size; ++j) { 128c7a4214aSPierre Jolivet set.insert(idx[j]); 129c7a4214aSPierre Jolivet for (PetscInt k=1; k<=ov; ++k) { /* for each layer of overlap */ 130c7a4214aSPierre Jolivet if (idx[j] - k >= 0) set.insert(idx[j] - k); /* do not insert negative indices */ 131c7a4214aSPierre 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 */ 132c7a4214aSPierre Jolivet } 133c7a4214aSPierre Jolivet } 1349566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is[i],&idx)); 1359566063dSJacob Faibussowitsch PetscCall(ISDestroy(is+i)); 1368f308287SPierre Jolivet if (bs[0] > 1) { 1378f308287SPierre Jolivet for (std::set<PetscInt>::iterator it=set.cbegin(); it!=set.cend(); it++) { 1388f308287SPierre Jolivet std::vector<PetscInt> block(bs[0]); 1398f308287SPierre Jolivet std::iota(block.begin(),block.end(),(*it/bs[0])*bs[0]); 1408f308287SPierre Jolivet set.insert(block.cbegin(),block.cend()); 1418f308287SPierre Jolivet } 1428f308287SPierre Jolivet } 143c7a4214aSPierre Jolivet size = set.size(); /* size with overlap */ 1449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&oidx)); 145c7a4214aSPierre Jolivet for (const PetscInt j : set) *oidx++ = j; 146c7a4214aSPierre Jolivet oidx -= size; 1479566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,size,oidx,PETSC_OWN_POINTER,is+i)); 148c7a4214aSPierre Jolivet } 149c7a4214aSPierre Jolivet PetscFunctionReturn(0); 150c7a4214aSPierre Jolivet } 151c7a4214aSPierre Jolivet 152c7a4214aSPierre Jolivet static PetscErrorCode MatCreateSubMatrices_Htool(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 153c7a4214aSPierre Jolivet { 154c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 155c7a4214aSPierre Jolivet Mat D,B,BT; 156c7a4214aSPierre Jolivet const PetscScalar *copy; 157c7a4214aSPierre Jolivet PetscScalar *ptr; 158c7a4214aSPierre Jolivet const PetscInt *idxr,*idxc,*it; 159c7a4214aSPierre Jolivet PetscInt nrow,m,i; 160c7a4214aSPierre Jolivet PetscBool flg; 161c7a4214aSPierre Jolivet 162c7a4214aSPierre Jolivet PetscFunctionBegin; 163c7a4214aSPierre Jolivet if (scall != MAT_REUSE_MATRIX) { 1649566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(n,submat)); 165c7a4214aSPierre Jolivet } 166c7a4214aSPierre Jolivet for (i=0; i<n; ++i) { 1679566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irow[i],&nrow)); 1689566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icol[i],&m)); 1699566063dSJacob Faibussowitsch PetscCall(ISGetIndices(irow[i],&idxr)); 1709566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icol[i],&idxc)); 171c7a4214aSPierre Jolivet if (scall != MAT_REUSE_MATRIX) { 1729566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF,nrow,m,nrow,m,NULL,(*submat)+i)); 173c7a4214aSPierre Jolivet } 1749566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite((*submat)[i],&ptr)); 175c7a4214aSPierre Jolivet if (irow[i] == icol[i]) { /* same row and column IS? */ 1769566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A,&flg)); 177c7a4214aSPierre Jolivet if (flg) { 1789566063dSJacob Faibussowitsch PetscCall(ISSorted(irow[i],&flg)); 179c7a4214aSPierre Jolivet if (flg) { /* sorted IS? */ 180c7a4214aSPierre Jolivet it = std::lower_bound(idxr,idxr+nrow,A->rmap->rstart); 181c7a4214aSPierre Jolivet if (it != idxr+nrow && *it == A->rmap->rstart) { /* rmap->rstart in IS? */ 182c7a4214aSPierre Jolivet if (std::distance(idxr,it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */ 183c7a4214aSPierre Jolivet for (PetscInt j=0; j<A->rmap->n && flg; ++j) if (PetscUnlikely(it[j] != A->rmap->rstart+j)) flg = PETSC_FALSE; 184c7a4214aSPierre Jolivet if (flg) { /* complete local diagonal block in IS? */ 185c7a4214aSPierre Jolivet /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM 186c7a4214aSPierre Jolivet * [ B C E ] 187c7a4214aSPierre Jolivet * A = [ B D E ] 188c7a4214aSPierre Jolivet * [ B F E ] 189c7a4214aSPierre Jolivet */ 190c7a4214aSPierre Jolivet m = std::distance(idxr,it); /* shift of the coefficient (0,0) of block D from above */ 1919566063dSJacob Faibussowitsch PetscCall(MatGetDiagonalBlock_Htool(A,&D)); 1929566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(D,©)); 193c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { 1949566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ptr+(m+k)*nrow+m,copy+k*A->rmap->n,A->rmap->n)); /* block D from above */ 195c7a4214aSPierre Jolivet } 1969566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(D,©)); 197c7a4214aSPierre Jolivet if (m) { 198c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* vertical block B from above */ 199c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 200b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 2019566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF,A->rmap->n,m,A->rmap->n,m,ptr+m,&B)); 2029566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B,nrow)); 2039566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF,m,A->rmap->n,m,A->rmap->n,ptr+m*nrow,&BT)); 2049566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT,nrow)); 205b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 2069566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT)); 207c7a4214aSPierre Jolivet } else { 208*7fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(B,BT)); 2099566063dSJacob Faibussowitsch PetscCall(MatTranspose(B,MAT_REUSE_MATRIX,&BT)); 210c7a4214aSPierre Jolivet } 2119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 213c7a4214aSPierre Jolivet } else { 214c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { /* block C from above */ 215c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(m,1,idxr,idxc+m+k,ptr+(m+k)*nrow); 216c7a4214aSPierre Jolivet } 217c7a4214aSPierre Jolivet } 218c7a4214aSPierre Jolivet } 219c7a4214aSPierre Jolivet if (m+A->rmap->n != nrow) { 220c7a4214aSPierre 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 */ 221c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 222b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 2239566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF,A->rmap->n,nrow-(m+A->rmap->n),A->rmap->n,nrow-(m+A->rmap->n),ptr+(m+A->rmap->n)*nrow+m,&B)); 2249566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B,nrow)); 2259566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF,nrow-(m+A->rmap->n),A->rmap->n,nrow-(m+A->rmap->n),A->rmap->n,ptr+m*nrow+m+A->rmap->n,&BT)); 2269566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT,nrow)); 227b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 2289566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT)); 229c7a4214aSPierre Jolivet } else { 230*7fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(B,BT)); 2319566063dSJacob Faibussowitsch PetscCall(MatTranspose(B,MAT_REUSE_MATRIX,&BT)); 232c7a4214aSPierre Jolivet } 2339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 235c7a4214aSPierre Jolivet } else { 236c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { /* block F from above */ 237c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(std::distance(it+A->rmap->n,idxr+nrow),1,it+A->rmap->n,idxc+m+k,ptr+(m+k)*nrow+m+A->rmap->n); 238c7a4214aSPierre Jolivet } 239c7a4214aSPierre Jolivet } 240c7a4214aSPierre Jolivet } 241c7a4214aSPierre Jolivet } /* complete local diagonal block not in IS */ 242c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */ 243c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* rmap->rstart not in IS */ 244c7a4214aSPierre Jolivet } /* unsorted IS */ 245c7a4214aSPierre Jolivet } 246c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* different row and column IS */ 247c7a4214aSPierre Jolivet if (!flg) a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* reassemble everything */ 2489566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irow[i],&idxr)); 2499566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icol[i],&idxc)); 2509566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite((*submat)[i],&ptr)); 2519566063dSJacob Faibussowitsch PetscCall(MatScale((*submat)[i],a->s)); 252c7a4214aSPierre Jolivet } 253c7a4214aSPierre Jolivet PetscFunctionReturn(0); 254c7a4214aSPierre Jolivet } 255c7a4214aSPierre Jolivet 256c7a4214aSPierre Jolivet static PetscErrorCode MatDestroy_Htool(Mat A) 257c7a4214aSPierre Jolivet { 258c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 259c7a4214aSPierre Jolivet PetscContainer container; 260c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 261c7a4214aSPierre Jolivet 262c7a4214aSPierre Jolivet PetscFunctionBegin; 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A,NULL)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",NULL)); 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",NULL)); 2669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",NULL)); 2679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",NULL)); 2689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",NULL)); 2699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",NULL)); 2709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",NULL)); 2719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",NULL)); 2729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",NULL)); 2739566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A,"KernelTranspose",(PetscObject*)&container)); 274c7a4214aSPierre Jolivet if (container) { /* created in MatTranspose_Htool() */ 2759566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container,(void**)&kernelt)); 2769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 2779566063dSJacob Faibussowitsch PetscCall(PetscFree(kernelt)); 2789566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&container)); 2799566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A,"KernelTranspose",NULL)); 280c7a4214aSPierre Jolivet } 281c7a4214aSPierre Jolivet if (a->gcoords_source != a->gcoords_target) { 2829566063dSJacob Faibussowitsch PetscCall(PetscFree(a->gcoords_source)); 283c7a4214aSPierre Jolivet } 2849566063dSJacob Faibussowitsch PetscCall(PetscFree(a->gcoords_target)); 2859566063dSJacob Faibussowitsch PetscCall(PetscFree2(a->work_source,a->work_target)); 286c7a4214aSPierre Jolivet delete a->wrapper; 287c7a4214aSPierre Jolivet delete a->hmatrix; 2889566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 289c7a4214aSPierre Jolivet PetscFunctionReturn(0); 290c7a4214aSPierre Jolivet } 291c7a4214aSPierre Jolivet 292c7a4214aSPierre Jolivet static PetscErrorCode MatView_Htool(Mat A,PetscViewer pv) 293c7a4214aSPierre Jolivet { 294c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 295c7a4214aSPierre Jolivet PetscBool flg; 296c7a4214aSPierre Jolivet 297c7a4214aSPierre Jolivet PetscFunctionBegin; 2989566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pv,PETSCVIEWERASCII,&flg)); 299c7a4214aSPierre Jolivet if (flg) { 3009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"symmetry: %c\n",a->hmatrix->get_symmetry_type())); 301c7a4214aSPierre Jolivet if (PetscAbsScalar(a->s-1.0) > PETSC_MACHINE_EPSILON) { 302c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX) 3039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"scaling: %g+%gi\n",(double)PetscRealPart(a->s),(double)PetscImaginaryPart(a->s))); 304c7a4214aSPierre Jolivet #else 3059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"scaling: %g\n",(double)a->s)); 306c7a4214aSPierre Jolivet #endif 307c7a4214aSPierre Jolivet } 3089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"minimum cluster size: %" PetscInt_FMT "\n",a->bs[0])); 3099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"maximum block size: %" PetscInt_FMT "\n",a->bs[1])); 3109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"epsilon: %g\n",(double)a->epsilon)); 3119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"eta: %g\n",(double)a->eta)); 3129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"minimum target depth: %" PetscInt_FMT "\n",a->depth[0])); 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"minimum source depth: %" PetscInt_FMT "\n",a->depth[1])); 3149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"compressor: %s\n",MatHtoolCompressorTypes[a->compressor])); 3159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"clustering: %s\n",MatHtoolClusteringTypes[a->clustering])); 3169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"compression ratio: %s\n",a->hmatrix->get_infos("Compression_ratio").c_str())); 3179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"space saving: %s\n",a->hmatrix->get_infos("Space_saving").c_str())); 3189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"number of dense (resp. low rank) matrices: %s (resp. %s)\n",a->hmatrix->get_infos("Number_of_dmat").c_str(),a->hmatrix->get_infos("Number_of_lrmat").c_str())); 3199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"(minimum, mean, maximum) dense block sizes: (%s, %s, %s)\n",a->hmatrix->get_infos("Dense_block_size_min").c_str(),a->hmatrix->get_infos("Dense_block_size_mean").c_str(),a->hmatrix->get_infos("Dense_block_size_max").c_str())); 3209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"(minimum, mean, maximum) low rank block sizes: (%s, %s, %s)\n",a->hmatrix->get_infos("Low_rank_block_size_min").c_str(),a->hmatrix->get_infos("Low_rank_block_size_mean").c_str(),a->hmatrix->get_infos("Low_rank_block_size_max").c_str())); 3219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"(minimum, mean, maximum) ranks: (%s, %s, %s)\n",a->hmatrix->get_infos("Rank_min").c_str(),a->hmatrix->get_infos("Rank_mean").c_str(),a->hmatrix->get_infos("Rank_max").c_str())); 322c7a4214aSPierre Jolivet } 323c7a4214aSPierre Jolivet PetscFunctionReturn(0); 324c7a4214aSPierre Jolivet } 325c7a4214aSPierre Jolivet 326c7a4214aSPierre Jolivet static PetscErrorCode MatScale_Htool(Mat A,PetscScalar s) 327c7a4214aSPierre Jolivet { 328c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 329c7a4214aSPierre Jolivet 330c7a4214aSPierre Jolivet PetscFunctionBegin; 331c7a4214aSPierre Jolivet a->s *= s; 332c7a4214aSPierre Jolivet PetscFunctionReturn(0); 333c7a4214aSPierre Jolivet } 334c7a4214aSPierre Jolivet 335c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */ 336c7a4214aSPierre Jolivet static PetscErrorCode MatGetRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 337c7a4214aSPierre Jolivet { 338c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 339c7a4214aSPierre Jolivet PetscInt *idxc; 340c7a4214aSPierre Jolivet PetscBLASInt one = 1,bn; 341c7a4214aSPierre Jolivet 342c7a4214aSPierre Jolivet PetscFunctionBegin; 343c7a4214aSPierre Jolivet if (nz) *nz = A->cmap->N; 344c7a4214aSPierre Jolivet if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */ 3459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N,&idxc)); 346c7a4214aSPierre Jolivet for (PetscInt i=0; i<A->cmap->N; ++i) idxc[i] = i; 347c7a4214aSPierre Jolivet } 348c7a4214aSPierre Jolivet if (idx) *idx = idxc; 349c7a4214aSPierre Jolivet if (v) { 3509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N,v)); 351c7a4214aSPierre Jolivet if (a->wrapper) a->wrapper->copy_submatrix(1,A->cmap->N,&row,idxc,*v); 352cab00e0dSPierre Jolivet else reinterpret_cast<htool::VirtualGenerator<PetscScalar>*>(a->kernelctx)->copy_submatrix(1,A->cmap->N,&row,idxc,*v); 3539566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(A->cmap->N,&bn)); 354792fecdfSBarry Smith PetscCallBLAS("BLASscal",BLASscal_(&bn,&a->s,*v,&one)); 355c7a4214aSPierre Jolivet } 356c7a4214aSPierre Jolivet if (!idx) { 3579566063dSJacob Faibussowitsch PetscCall(PetscFree(idxc)); 358c7a4214aSPierre Jolivet } 359c7a4214aSPierre Jolivet PetscFunctionReturn(0); 360c7a4214aSPierre Jolivet } 361c7a4214aSPierre Jolivet 362c7a4214aSPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 363c7a4214aSPierre Jolivet { 364c7a4214aSPierre Jolivet PetscFunctionBegin; 365c7a4214aSPierre Jolivet if (nz) *nz = 0; 366c7a4214aSPierre Jolivet if (idx) { 3679566063dSJacob Faibussowitsch PetscCall(PetscFree(*idx)); 368c7a4214aSPierre Jolivet } 369c7a4214aSPierre Jolivet if (v) { 3709566063dSJacob Faibussowitsch PetscCall(PetscFree(*v)); 371c7a4214aSPierre Jolivet } 372c7a4214aSPierre Jolivet PetscFunctionReturn(0); 373c7a4214aSPierre Jolivet } 374c7a4214aSPierre Jolivet 375c7a4214aSPierre Jolivet static PetscErrorCode MatSetFromOptions_Htool(PetscOptionItems *PetscOptionsObject,Mat A) 376c7a4214aSPierre Jolivet { 377c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 378c7a4214aSPierre Jolivet PetscInt n; 379c7a4214aSPierre Jolivet PetscBool flg; 380c7a4214aSPierre Jolivet 381c7a4214aSPierre Jolivet PetscFunctionBegin; 382d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Htool options"); 3839566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_min_cluster_size","Minimal leaf size in cluster tree",NULL,a->bs[0],a->bs,NULL)); 3849566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_max_block_size","Maximal number of coefficients in a dense block",NULL,a->bs[1],a->bs + 1,NULL)); 3859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_htool_epsilon","Relative error in Frobenius norm when approximating a block",NULL,a->epsilon,&a->epsilon,NULL)); 3869566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_htool_eta","Admissibility condition tolerance",NULL,a->eta,&a->eta,NULL)); 3879566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_min_target_depth","Minimal cluster tree depth associated with the rows",NULL,a->depth[0],a->depth,NULL)); 3889566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_min_source_depth","Minimal cluster tree depth associated with the columns",NULL,a->depth[1],a->depth + 1,NULL)); 389c7a4214aSPierre Jolivet n = 0; 390dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_compressor","Type of compression","MatHtoolCompressorType",MatHtoolCompressorTypes,PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes),MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA],&n,&flg)); 391c7a4214aSPierre Jolivet if (flg) a->compressor = MatHtoolCompressorType(n); 39298e73e17SPierre Jolivet n = 0; 393dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_clustering","Type of clustering","MatHtoolClusteringType",MatHtoolClusteringTypes,PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes),MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR],&n,&flg)); 39498e73e17SPierre Jolivet if (flg) a->clustering = MatHtoolClusteringType(n); 395d0609cedSBarry Smith PetscOptionsHeadEnd(); 396c7a4214aSPierre Jolivet PetscFunctionReturn(0); 397c7a4214aSPierre Jolivet } 398c7a4214aSPierre Jolivet 399c7a4214aSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A,MatAssemblyType type) 400c7a4214aSPierre Jolivet { 401c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 402c7a4214aSPierre Jolivet const PetscInt *ranges; 403c7a4214aSPierre Jolivet PetscInt *offset; 404c7a4214aSPierre Jolivet PetscMPIInt size; 405b94d7dedSBarry Smith char S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'),uplo = S == 'N' ? 'N' : 'U'; 406cab00e0dSPierre Jolivet htool::VirtualGenerator<PetscScalar> *generator = nullptr; 40798e73e17SPierre Jolivet std::shared_ptr<htool::VirtualCluster> t,s = nullptr; 4083b9338faSPierre Jolivet std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor = nullptr; 409c7a4214aSPierre Jolivet 410c7a4214aSPierre Jolivet PetscFunctionBegin; 4119566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(HtoolCitation,&HtoolCite)); 412c7a4214aSPierre Jolivet delete a->wrapper; 413c7a4214aSPierre Jolivet delete a->hmatrix; 4149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 4159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2*size,&offset)); 4169566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A,&ranges)); 417c7a4214aSPierre Jolivet for (PetscInt i=0; i<size; ++i) { 418c7a4214aSPierre Jolivet offset[2*i] = ranges[i]; 419c7a4214aSPierre Jolivet offset[2*i+1] = ranges[i+1] - ranges[i]; 420c7a4214aSPierre Jolivet } 42198e73e17SPierre Jolivet switch (a->clustering) { 42298e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 42398e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 42498e73e17SPierre Jolivet break; 42598e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 42698e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 42798e73e17SPierre Jolivet break; 42898e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 42998e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 43098e73e17SPierre Jolivet break; 43198e73e17SPierre Jolivet default: 43298e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 43398e73e17SPierre Jolivet } 434c7a4214aSPierre Jolivet t->set_minclustersize(a->bs[0]); 43598e73e17SPierre Jolivet t->build(A->rmap->N,a->gcoords_target,offset); 436c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx); 437c7a4214aSPierre Jolivet else { 438c7a4214aSPierre Jolivet a->wrapper = NULL; 439cab00e0dSPierre Jolivet generator = reinterpret_cast<htool::VirtualGenerator<PetscScalar>*>(a->kernelctx); 440c7a4214aSPierre Jolivet } 441c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 4429566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(A,&ranges)); 443c7a4214aSPierre Jolivet for (PetscInt i=0; i<size; ++i) { 444c7a4214aSPierre Jolivet offset[2*i] = ranges[i]; 445c7a4214aSPierre Jolivet offset[2*i+1] = ranges[i+1] - ranges[i]; 446c7a4214aSPierre Jolivet } 44798e73e17SPierre Jolivet switch (a->clustering) { 44898e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 44998e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 45098e73e17SPierre Jolivet break; 45198e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 45298e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 45398e73e17SPierre Jolivet break; 45498e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 45598e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 45698e73e17SPierre Jolivet break; 45798e73e17SPierre Jolivet default: 45898e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 45998e73e17SPierre Jolivet } 460c7a4214aSPierre Jolivet s->set_minclustersize(a->bs[0]); 46198e73e17SPierre Jolivet s->build(A->cmap->N,a->gcoords_source,offset); 462c7a4214aSPierre Jolivet S = uplo = 'N'; 463c7a4214aSPierre Jolivet } 4649566063dSJacob Faibussowitsch PetscCall(PetscFree(offset)); 465c7a4214aSPierre Jolivet switch (a->compressor) { 466c7a4214aSPierre Jolivet case MAT_HTOOL_COMPRESSOR_FULL_ACA: 4673b9338faSPierre Jolivet compressor = std::make_shared<htool::fullACA<PetscScalar>>(); 468c7a4214aSPierre Jolivet break; 469c7a4214aSPierre Jolivet case MAT_HTOOL_COMPRESSOR_SVD: 4703b9338faSPierre Jolivet compressor = std::make_shared<htool::SVD<PetscScalar>>(); 471c7a4214aSPierre Jolivet break; 472c7a4214aSPierre Jolivet default: 4733b9338faSPierre Jolivet compressor = std::make_shared<htool::sympartialACA<PetscScalar>>(); 474c7a4214aSPierre Jolivet } 4753b9338faSPierre Jolivet a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar>*>(new htool::HMatrix<PetscScalar>(t,s ? s : t,a->epsilon,a->eta,S,uplo)); 4763b9338faSPierre Jolivet a->hmatrix->set_compression(compressor); 477c7a4214aSPierre Jolivet a->hmatrix->set_maxblocksize(a->bs[1]); 478c7a4214aSPierre Jolivet a->hmatrix->set_mintargetdepth(a->depth[0]); 479c7a4214aSPierre Jolivet a->hmatrix->set_minsourcedepth(a->depth[1]); 4803b9338faSPierre Jolivet if (s) a->hmatrix->build(a->wrapper ? *a->wrapper : *generator,a->gcoords_target,a->gcoords_source); 4813b9338faSPierre Jolivet else a->hmatrix->build(a->wrapper ? *a->wrapper : *generator,a->gcoords_target); 482c7a4214aSPierre Jolivet PetscFunctionReturn(0); 483c7a4214aSPierre Jolivet } 484c7a4214aSPierre Jolivet 485c7a4214aSPierre Jolivet static PetscErrorCode MatProductNumeric_Htool(Mat C) 486c7a4214aSPierre Jolivet { 487c7a4214aSPierre Jolivet Mat_Product *product = C->product; 488c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)product->A->data; 489c7a4214aSPierre Jolivet const PetscScalar *in; 490c7a4214aSPierre Jolivet PetscScalar *out; 4918b8ddd11SPierre Jolivet PetscInt N,lda; 492c7a4214aSPierre Jolivet 493c7a4214aSPierre Jolivet PetscFunctionBegin; 494c7a4214aSPierre Jolivet MatCheckProduct(C,1); 4959566063dSJacob Faibussowitsch PetscCall(MatGetSize(C,NULL,&N)); 4969566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C,&lda)); 4975f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")",lda,C->rmap->n); 4989566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(product->B,&in)); 4999566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C,&out)); 5008b8ddd11SPierre Jolivet switch (product->type) { 5018b8ddd11SPierre Jolivet case MATPRODUCT_AB: 502c7a4214aSPierre Jolivet a->hmatrix->mvprod_local_to_local(in,out,N); 5038b8ddd11SPierre Jolivet break; 5048b8ddd11SPierre Jolivet case MATPRODUCT_AtB: 5058b8ddd11SPierre Jolivet a->hmatrix->mvprod_transp_local_to_local(in,out,N); 506c7a4214aSPierre Jolivet break; 507c7a4214aSPierre Jolivet default: 50898921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType %s is not supported",MatProductTypes[product->type]); 509c7a4214aSPierre Jolivet } 5109566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C,&out)); 5119566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(product->B,&in)); 5129566063dSJacob Faibussowitsch PetscCall(MatScale(C,a->s)); 513c7a4214aSPierre Jolivet PetscFunctionReturn(0); 514c7a4214aSPierre Jolivet } 515c7a4214aSPierre Jolivet 516c7a4214aSPierre Jolivet static PetscErrorCode MatProductSymbolic_Htool(Mat C) 517c7a4214aSPierre Jolivet { 518c7a4214aSPierre Jolivet Mat_Product *product = C->product; 519c7a4214aSPierre Jolivet Mat A,B; 520c7a4214aSPierre Jolivet PetscBool flg; 521c7a4214aSPierre Jolivet 522c7a4214aSPierre Jolivet PetscFunctionBegin; 523c7a4214aSPierre Jolivet MatCheckProduct(C,1); 524c7a4214aSPierre Jolivet A = product->A; 525c7a4214aSPierre Jolivet B = product->B; 5269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,"")); 5275f80ce2aSJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MatProduct_AB not supported for %s",((PetscObject)product->B)->type_name); 528c7a4214aSPierre Jolivet switch (product->type) { 529c7a4214aSPierre Jolivet case MATPRODUCT_AB: 530c7a4214aSPierre Jolivet if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) { 5319566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,A->rmap->n,B->cmap->n,A->rmap->N,B->cmap->N)); 532c7a4214aSPierre Jolivet } 5338b8ddd11SPierre Jolivet break; 5348b8ddd11SPierre Jolivet case MATPRODUCT_AtB: 5358b8ddd11SPierre Jolivet if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) { 5369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,A->cmap->n,B->cmap->n,A->cmap->N,B->cmap->N)); 5378b8ddd11SPierre Jolivet } 5388b8ddd11SPierre Jolivet break; 5398b8ddd11SPierre Jolivet default: 54098921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"ProductType %s is not supported",MatProductTypes[product->type]); 5418b8ddd11SPierre Jolivet } 5429566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATDENSE)); 5439566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 5449566063dSJacob Faibussowitsch PetscCall(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 5459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 5469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 547c7a4214aSPierre Jolivet C->ops->productsymbolic = NULL; 548c7a4214aSPierre Jolivet C->ops->productnumeric = MatProductNumeric_Htool; 549c7a4214aSPierre Jolivet PetscFunctionReturn(0); 550c7a4214aSPierre Jolivet } 551c7a4214aSPierre Jolivet 552c7a4214aSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Htool(Mat C) 553c7a4214aSPierre Jolivet { 554c7a4214aSPierre Jolivet PetscFunctionBegin; 555c7a4214aSPierre Jolivet MatCheckProduct(C,1); 5568b8ddd11SPierre Jolivet if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool; 557c7a4214aSPierre Jolivet PetscFunctionReturn(0); 558c7a4214aSPierre Jolivet } 559c7a4214aSPierre Jolivet 56098e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix) 561c7a4214aSPierre Jolivet { 562c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 563c7a4214aSPierre Jolivet 564c7a4214aSPierre Jolivet PetscFunctionBegin; 565c7a4214aSPierre Jolivet *hmatrix = a->hmatrix; 566c7a4214aSPierre Jolivet PetscFunctionReturn(0); 567c7a4214aSPierre Jolivet } 568c7a4214aSPierre Jolivet 569c7a4214aSPierre Jolivet /*@C 570c7a4214aSPierre Jolivet MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a MATHTOOL. 571c7a4214aSPierre Jolivet 572c7a4214aSPierre Jolivet Input Parameter: 573c7a4214aSPierre Jolivet . A - hierarchical matrix 574c7a4214aSPierre Jolivet 575c7a4214aSPierre Jolivet Output Parameter: 576c7a4214aSPierre Jolivet . hmatrix - opaque pointer to a Htool virtual matrix 577c7a4214aSPierre Jolivet 578c7a4214aSPierre Jolivet Level: advanced 579c7a4214aSPierre Jolivet 580db781477SPatrick Sanan .seealso: `MATHTOOL` 581c7a4214aSPierre Jolivet @*/ 58298e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix) 583c7a4214aSPierre Jolivet { 584c7a4214aSPierre Jolivet PetscFunctionBegin; 585c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 586c7a4214aSPierre Jolivet PetscValidPointer(hmatrix,2); 587cac4c232SBarry Smith PetscTryMethod(A,"MatHtoolGetHierarchicalMat_C",(Mat,const htool::VirtualHMatrix<PetscScalar>**),(A,hmatrix)); 588c7a4214aSPierre Jolivet PetscFunctionReturn(0); 589c7a4214aSPierre Jolivet } 590c7a4214aSPierre Jolivet 591c7a4214aSPierre Jolivet static PetscErrorCode MatHtoolSetKernel_Htool(Mat A,MatHtoolKernel kernel,void *kernelctx) 592c7a4214aSPierre Jolivet { 593c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 594c7a4214aSPierre Jolivet 595c7a4214aSPierre Jolivet PetscFunctionBegin; 596c7a4214aSPierre Jolivet a->kernel = kernel; 597c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 598c7a4214aSPierre Jolivet delete a->wrapper; 599c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx); 600c7a4214aSPierre Jolivet PetscFunctionReturn(0); 601c7a4214aSPierre Jolivet } 602c7a4214aSPierre Jolivet 603c7a4214aSPierre Jolivet /*@C 60498e73e17SPierre Jolivet MatHtoolSetKernel - Sets the kernel and context used for the assembly of a MATHTOOL. 605c7a4214aSPierre Jolivet 606c7a4214aSPierre Jolivet Input Parameters: 607c7a4214aSPierre Jolivet + A - hierarchical matrix 608c7a4214aSPierre Jolivet . kernel - computational kernel (or NULL) 609cab00e0dSPierre Jolivet - kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 610c7a4214aSPierre Jolivet 611c7a4214aSPierre Jolivet Level: advanced 612c7a4214aSPierre Jolivet 613db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatCreateHtoolFromKernel()` 614c7a4214aSPierre Jolivet @*/ 615c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat A,MatHtoolKernel kernel,void *kernelctx) 616c7a4214aSPierre Jolivet { 617c7a4214aSPierre Jolivet PetscFunctionBegin; 618c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 619c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel,2); 620c7a4214aSPierre Jolivet if (!kernel) PetscValidPointer(kernelctx,3); 621cac4c232SBarry Smith PetscTryMethod(A,"MatHtoolSetKernel_C",(Mat,MatHtoolKernel,void*),(A,kernel,kernelctx)); 622c7a4214aSPierre Jolivet PetscFunctionReturn(0); 623c7a4214aSPierre Jolivet } 624c7a4214aSPierre Jolivet 62598e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A,IS* is) 62698e73e17SPierre Jolivet { 62798e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 62898e73e17SPierre Jolivet std::vector<PetscInt> source; 62998e73e17SPierre Jolivet 63098e73e17SPierre Jolivet PetscFunctionBegin; 631a9918087SPierre Jolivet source = a->hmatrix->get_source_cluster()->get_local_perm(); 6329566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A),source.size(),source.data(),PETSC_COPY_VALUES,is)); 6339566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 63498e73e17SPierre Jolivet PetscFunctionReturn(0); 63598e73e17SPierre Jolivet } 63698e73e17SPierre Jolivet 63798e73e17SPierre Jolivet /*@C 63898e73e17SPierre Jolivet MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster. 63998e73e17SPierre Jolivet 64098e73e17SPierre Jolivet Input Parameter: 64198e73e17SPierre Jolivet . A - hierarchical matrix 64298e73e17SPierre Jolivet 64398e73e17SPierre Jolivet Output Parameter: 64498e73e17SPierre Jolivet . is - permutation 64598e73e17SPierre Jolivet 64698e73e17SPierre Jolivet Level: advanced 64798e73e17SPierre Jolivet 648db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()` 64998e73e17SPierre Jolivet @*/ 65098e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat A,IS* is) 65198e73e17SPierre Jolivet { 65298e73e17SPierre Jolivet PetscFunctionBegin; 65398e73e17SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 65498e73e17SPierre Jolivet if (!is) PetscValidPointer(is,2); 655cac4c232SBarry Smith PetscTryMethod(A,"MatHtoolGetPermutationSource_C",(Mat,IS*),(A,is)); 65698e73e17SPierre Jolivet PetscFunctionReturn(0); 65798e73e17SPierre Jolivet } 65898e73e17SPierre Jolivet 65998e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A,IS* is) 66098e73e17SPierre Jolivet { 66198e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 66298e73e17SPierre Jolivet std::vector<PetscInt> target; 66398e73e17SPierre Jolivet 66498e73e17SPierre Jolivet PetscFunctionBegin; 665a9918087SPierre Jolivet target = a->hmatrix->get_target_cluster()->get_local_perm(); 6669566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A),target.size(),target.data(),PETSC_COPY_VALUES,is)); 6679566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 66898e73e17SPierre Jolivet PetscFunctionReturn(0); 66998e73e17SPierre Jolivet } 67098e73e17SPierre Jolivet 67198e73e17SPierre Jolivet /*@C 67298e73e17SPierre Jolivet MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster. 67398e73e17SPierre Jolivet 67498e73e17SPierre Jolivet Input Parameter: 67598e73e17SPierre Jolivet . A - hierarchical matrix 67698e73e17SPierre Jolivet 67798e73e17SPierre Jolivet Output Parameter: 67898e73e17SPierre Jolivet . is - permutation 67998e73e17SPierre Jolivet 68098e73e17SPierre Jolivet Level: advanced 68198e73e17SPierre Jolivet 682db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()` 68398e73e17SPierre Jolivet @*/ 68498e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat A,IS* is) 68598e73e17SPierre Jolivet { 68698e73e17SPierre Jolivet PetscFunctionBegin; 68798e73e17SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 68898e73e17SPierre Jolivet if (!is) PetscValidPointer(is,2); 689cac4c232SBarry Smith PetscTryMethod(A,"MatHtoolGetPermutationTarget_C",(Mat,IS*),(A,is)); 69098e73e17SPierre Jolivet PetscFunctionReturn(0); 69198e73e17SPierre Jolivet } 69298e73e17SPierre Jolivet 69398e73e17SPierre Jolivet static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A,PetscBool use) 69498e73e17SPierre Jolivet { 69598e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 69698e73e17SPierre Jolivet 69798e73e17SPierre Jolivet PetscFunctionBegin; 69898e73e17SPierre Jolivet a->hmatrix->set_use_permutation(use); 69998e73e17SPierre Jolivet PetscFunctionReturn(0); 70098e73e17SPierre Jolivet } 70198e73e17SPierre Jolivet 70298e73e17SPierre Jolivet /*@C 70398e73e17SPierre Jolivet MatHtoolUsePermutation - Sets whether MATHTOOL should permute input (resp. output) vectors following its internal source (resp. target) permutation. 70498e73e17SPierre Jolivet 70598e73e17SPierre Jolivet Input Parameters: 70698e73e17SPierre Jolivet + A - hierarchical matrix 70798e73e17SPierre Jolivet - use - Boolean value 70898e73e17SPierre Jolivet 70998e73e17SPierre Jolivet Level: advanced 71098e73e17SPierre Jolivet 711db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()` 71298e73e17SPierre Jolivet @*/ 71398e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat A,PetscBool use) 71498e73e17SPierre Jolivet { 71598e73e17SPierre Jolivet PetscFunctionBegin; 71698e73e17SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 71798e73e17SPierre Jolivet PetscValidLogicalCollectiveBool(A,use,2); 718cac4c232SBarry Smith PetscTryMethod(A,"MatHtoolUsePermutation_C",(Mat,PetscBool),(A,use)); 71998e73e17SPierre Jolivet PetscFunctionReturn(0); 72098e73e17SPierre Jolivet } 72198e73e17SPierre Jolivet 722c7a4214aSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B) 723c7a4214aSPierre Jolivet { 724c7a4214aSPierre Jolivet Mat C; 725c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 726c7a4214aSPierre Jolivet PetscInt lda; 727c7a4214aSPierre Jolivet PetscScalar *array; 728c7a4214aSPierre Jolivet 729c7a4214aSPierre Jolivet PetscFunctionBegin; 730c7a4214aSPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 731c7a4214aSPierre Jolivet C = *B; 7325f80ce2aSJacob Faibussowitsch PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible dimensions"); 7339566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C,&lda)); 7345f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")",lda,C->rmap->n); 735c7a4214aSPierre Jolivet } else { 7369566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C)); 7379566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N)); 7389566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATDENSE)); 7399566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7409566063dSJacob Faibussowitsch PetscCall(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 741c7a4214aSPierre Jolivet } 7429566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 7439566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 7449566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C,&array)); 745c7a4214aSPierre Jolivet a->hmatrix->copy_local_dense_perm(array); 7469566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C,&array)); 7479566063dSJacob Faibussowitsch PetscCall(MatScale(C,a->s)); 748c7a4214aSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 7499566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A,&C)); 750c7a4214aSPierre Jolivet } else *B = C; 751c7a4214aSPierre Jolivet PetscFunctionReturn(0); 752c7a4214aSPierre Jolivet } 753c7a4214aSPierre Jolivet 75498e73e17SPierre Jolivet static PetscErrorCode GenEntriesTranspose(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *rows,const PetscInt *cols,PetscScalar *ptr,void *ctx) 755c7a4214aSPierre Jolivet { 756c7a4214aSPierre Jolivet MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose*)ctx; 75798e73e17SPierre Jolivet PetscScalar *tmp; 75898e73e17SPierre Jolivet 75998e73e17SPierre Jolivet PetscFunctionBegin; 76098e73e17SPierre Jolivet generator->kernel(sdim,N,M,cols,rows,ptr,generator->kernelctx); 7619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M*N,&tmp)); 7629566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp,ptr,M*N)); 76398e73e17SPierre Jolivet for (PetscInt i=0; i<M; ++i) { 76498e73e17SPierre Jolivet for (PetscInt j=0; j<N; ++j) ptr[i+j*M] = tmp[j+i*N]; 76598e73e17SPierre Jolivet } 7669566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 76798e73e17SPierre Jolivet PetscFunctionReturn(0); 768c7a4214aSPierre Jolivet } 769c7a4214aSPierre Jolivet 770c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */ 771c7a4214aSPierre Jolivet static PetscErrorCode MatTranspose_Htool(Mat A,MatReuse reuse,Mat *B) 772c7a4214aSPierre Jolivet { 773c7a4214aSPierre Jolivet Mat C; 774c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data,*c; 775c7a4214aSPierre Jolivet PetscInt M = A->rmap->N,N = A->cmap->N,m = A->rmap->n,n = A->cmap->n; 776c7a4214aSPierre Jolivet PetscContainer container; 777c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 778c7a4214aSPierre Jolivet 779c7a4214aSPierre Jolivet PetscFunctionBegin; 780*7fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A,*B)); 7815f80ce2aSJacob Faibussowitsch PetscCheck(reuse != MAT_INPLACE_MATRIX,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatTranspose() with MAT_INPLACE_MATRIX not supported"); 782c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 7839566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C)); 7849566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,n,m,N,M)); 7859566063dSJacob Faibussowitsch PetscCall(MatSetType(C,((PetscObject)A)->type_name)); 7869566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7879566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)C),&container)); 7889566063dSJacob Faibussowitsch PetscCall(PetscNew(&kernelt)); 7899566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(container,kernelt)); 7909566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)C,"KernelTranspose",(PetscObject)container)); 791c7a4214aSPierre Jolivet } else { 792c7a4214aSPierre Jolivet C = *B; 7939566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)C,"KernelTranspose",(PetscObject*)&container)); 7945f80ce2aSJacob Faibussowitsch PetscCheck(container,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatTranspose() with MAT_INITIAL_MATRIX first"); 7959566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container,(void**)&kernelt)); 796c7a4214aSPierre Jolivet } 797c7a4214aSPierre Jolivet c = (Mat_Htool*)C->data; 798c7a4214aSPierre Jolivet c->dim = a->dim; 799c7a4214aSPierre Jolivet c->s = a->s; 80098e73e17SPierre Jolivet c->kernel = GenEntriesTranspose; 801c7a4214aSPierre Jolivet if (kernelt->A != A) { 8029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 803c7a4214aSPierre Jolivet kernelt->A = A; 8049566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 805c7a4214aSPierre Jolivet } 806c7a4214aSPierre Jolivet kernelt->kernel = a->kernel; 807c7a4214aSPierre Jolivet kernelt->kernelctx = a->kernelctx; 808c7a4214aSPierre Jolivet c->kernelctx = kernelt; 809c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 8109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N*c->dim,&c->gcoords_target)); 8119566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_target,a->gcoords_source,N*c->dim)); 812c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 8139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M*c->dim,&c->gcoords_source)); 8149566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_source,a->gcoords_target,M*c->dim)); 815c7a4214aSPierre Jolivet } else c->gcoords_source = c->gcoords_target; 8169566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(M,&c->work_source,N,&c->work_target)); 817c7a4214aSPierre Jolivet } 8189566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 8199566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 820c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) *B = C; 821c7a4214aSPierre Jolivet PetscFunctionReturn(0); 822c7a4214aSPierre Jolivet } 823c7a4214aSPierre Jolivet 824c7a4214aSPierre Jolivet /*@C 825c7a4214aSPierre Jolivet MatCreateHtoolFromKernel - Creates a MATHTOOL from a user-supplied kernel. 826c7a4214aSPierre Jolivet 827c7a4214aSPierre Jolivet Input Parameters: 828c7a4214aSPierre Jolivet + comm - MPI communicator 829c7a4214aSPierre Jolivet . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 830c7a4214aSPierre Jolivet . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 831c7a4214aSPierre Jolivet . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 832c7a4214aSPierre Jolivet . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 833c7a4214aSPierre Jolivet . spacedim - dimension of the space coordinates 834c7a4214aSPierre Jolivet . coords_target - coordinates of the target 835c7a4214aSPierre Jolivet . coords_source - coordinates of the source 836c7a4214aSPierre Jolivet . kernel - computational kernel (or NULL) 837cab00e0dSPierre Jolivet - kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 838c7a4214aSPierre Jolivet 839c7a4214aSPierre Jolivet Output Parameter: 840c7a4214aSPierre Jolivet . B - matrix 841c7a4214aSPierre Jolivet 842c7a4214aSPierre Jolivet Options Database Keys: 843c7a4214aSPierre Jolivet + -mat_htool_min_cluster_size <PetscInt> - minimal leaf size in cluster tree 844c7a4214aSPierre Jolivet . -mat_htool_max_block_size <PetscInt> - maximal number of coefficients in a dense block 845c7a4214aSPierre Jolivet . -mat_htool_epsilon <PetscReal> - relative error in Frobenius norm when approximating a block 846c7a4214aSPierre Jolivet . -mat_htool_eta <PetscReal> - admissibility condition tolerance 847c7a4214aSPierre Jolivet . -mat_htool_min_target_depth <PetscInt> - minimal cluster tree depth associated with the rows 848c7a4214aSPierre Jolivet . -mat_htool_min_source_depth <PetscInt> - minimal cluster tree depth associated with the columns 84998e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression 85098e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering 851c7a4214aSPierre Jolivet 852c7a4214aSPierre Jolivet Level: intermediate 853c7a4214aSPierre Jolivet 854db781477SPatrick Sanan .seealso: `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()` 855c7a4214aSPierre Jolivet @*/ 856c7a4214aSPierre 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) 857c7a4214aSPierre Jolivet { 858c7a4214aSPierre Jolivet Mat A; 859c7a4214aSPierre Jolivet Mat_Htool *a; 860c7a4214aSPierre Jolivet 861c7a4214aSPierre Jolivet PetscFunctionBegin; 8629566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&A)); 863c7a4214aSPierre Jolivet PetscValidLogicalCollectiveInt(A,spacedim,6); 864c7a4214aSPierre Jolivet PetscValidRealPointer(coords_target,7); 865c7a4214aSPierre Jolivet PetscValidRealPointer(coords_source,8); 866c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel,9); 867c7a4214aSPierre Jolivet if (!kernel) PetscValidPointer(kernelctx,10); 8689566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,m,n,M,N)); 8699566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATHTOOL)); 8709566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 871c7a4214aSPierre Jolivet a = (Mat_Htool*)A->data; 872c7a4214aSPierre Jolivet a->dim = spacedim; 873c7a4214aSPierre Jolivet a->s = 1.0; 874c7a4214aSPierre Jolivet a->kernel = kernel; 875c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 8769566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->rmap->N*spacedim,&a->gcoords_target)); 8779566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_target+A->rmap->rstart*spacedim,coords_target,A->rmap->n*spacedim)); 8781c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE,a->gcoords_target,A->rmap->N*spacedim,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)A))); /* global target coordinates */ 879c7a4214aSPierre Jolivet if (coords_target != coords_source) { 8809566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->cmap->N*spacedim,&a->gcoords_source)); 8819566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_source+A->cmap->rstart*spacedim,coords_source,A->cmap->n*spacedim)); 8821c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE,a->gcoords_source,A->cmap->N*spacedim,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)A))); /* global source coordinates */ 883c7a4214aSPierre Jolivet } else a->gcoords_source = a->gcoords_target; 8849566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(A->cmap->N,&a->work_source,A->rmap->N,&a->work_target)); 885c7a4214aSPierre Jolivet *B = A; 886c7a4214aSPierre Jolivet PetscFunctionReturn(0); 887c7a4214aSPierre Jolivet } 888c7a4214aSPierre Jolivet 889c7a4214aSPierre Jolivet /*MC 890c7a4214aSPierre Jolivet MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package. 891c7a4214aSPierre Jolivet 892c7a4214aSPierre Jolivet Use ./configure --download-htool to install PETSc to use Htool. 893c7a4214aSPierre Jolivet 894c7a4214aSPierre Jolivet Options Database Keys: 895e85a1b63SPatrick Sanan . -mat_type htool - matrix type to "htool" during a call to MatSetFromOptions() 896c7a4214aSPierre Jolivet 897c7a4214aSPierre Jolivet Level: beginner 898c7a4214aSPierre Jolivet 899db781477SPatrick Sanan .seealso: `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()` 900c7a4214aSPierre Jolivet M*/ 901c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A) 902c7a4214aSPierre Jolivet { 903c7a4214aSPierre Jolivet Mat_Htool *a; 904c7a4214aSPierre Jolivet 905c7a4214aSPierre Jolivet PetscFunctionBegin; 9069566063dSJacob Faibussowitsch PetscCall(PetscNewLog(A,&a)); 907c7a4214aSPierre Jolivet A->data = (void*)a; 9089566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATHTOOL)); 9099566063dSJacob Faibussowitsch PetscCall(PetscMemzero(A->ops,sizeof(struct _MatOps))); 910c7a4214aSPierre Jolivet A->ops->getdiagonal = MatGetDiagonal_Htool; 911c7a4214aSPierre Jolivet A->ops->getdiagonalblock = MatGetDiagonalBlock_Htool; 912c7a4214aSPierre Jolivet A->ops->mult = MatMult_Htool; 913c7a4214aSPierre Jolivet A->ops->multadd = MatMultAdd_Htool; 9148b8ddd11SPierre Jolivet A->ops->multtranspose = MatMultTranspose_Htool; 9158b8ddd11SPierre Jolivet if (!PetscDefined(USE_COMPLEX)) A->ops->multhermitiantranspose = MatMultTranspose_Htool; 916c7a4214aSPierre Jolivet A->ops->increaseoverlap = MatIncreaseOverlap_Htool; 917c7a4214aSPierre Jolivet A->ops->createsubmatrices = MatCreateSubMatrices_Htool; 918c7a4214aSPierre Jolivet A->ops->transpose = MatTranspose_Htool; 919c7a4214aSPierre Jolivet A->ops->destroy = MatDestroy_Htool; 920c7a4214aSPierre Jolivet A->ops->view = MatView_Htool; 921c7a4214aSPierre Jolivet A->ops->setfromoptions = MatSetFromOptions_Htool; 922c7a4214aSPierre Jolivet A->ops->scale = MatScale_Htool; 923c7a4214aSPierre Jolivet A->ops->getrow = MatGetRow_Htool; 924c7a4214aSPierre Jolivet A->ops->restorerow = MatRestoreRow_Htool; 925c7a4214aSPierre Jolivet A->ops->assemblyend = MatAssemblyEnd_Htool; 926c7a4214aSPierre Jolivet a->dim = 0; 927c7a4214aSPierre Jolivet a->gcoords_target = NULL; 928c7a4214aSPierre Jolivet a->gcoords_source = NULL; 929c7a4214aSPierre Jolivet a->s = 1.0; 930c7a4214aSPierre Jolivet a->bs[0] = 10; 931c7a4214aSPierre Jolivet a->bs[1] = 1000000; 932c7a4214aSPierre Jolivet a->epsilon = PetscSqrtReal(PETSC_SMALL); 933c7a4214aSPierre Jolivet a->eta = 10.0; 934c7a4214aSPierre Jolivet a->depth[0] = 0; 935c7a4214aSPierre Jolivet a->depth[1] = 0; 936c7a4214aSPierre Jolivet a->compressor = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA; 9379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",MatProductSetFromOptions_Htool)); 9389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",MatProductSetFromOptions_Htool)); 9399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",MatConvert_Htool_Dense)); 9409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",MatConvert_Htool_Dense)); 9419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",MatHtoolGetHierarchicalMat_Htool)); 9429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",MatHtoolSetKernel_Htool)); 9439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",MatHtoolGetPermutationSource_Htool)); 9449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",MatHtoolGetPermutationTarget_Htool)); 9459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",MatHtoolUsePermutation_Htool)); 946c7a4214aSPierre Jolivet PetscFunctionReturn(0); 947c7a4214aSPierre Jolivet } 948