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 */ 200*b94d7dedSBarry 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)); 205*b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 2069566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT)); 207c7a4214aSPierre Jolivet } else { 2089566063dSJacob Faibussowitsch PetscCall(MatTranspose(B,MAT_REUSE_MATRIX,&BT)); 209c7a4214aSPierre Jolivet } 2109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 212c7a4214aSPierre Jolivet } else { 213c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { /* block C from above */ 214c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(m,1,idxr,idxc+m+k,ptr+(m+k)*nrow); 215c7a4214aSPierre Jolivet } 216c7a4214aSPierre Jolivet } 217c7a4214aSPierre Jolivet } 218c7a4214aSPierre Jolivet if (m+A->rmap->n != nrow) { 219c7a4214aSPierre 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 */ 220c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 221*b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 2229566063dSJacob 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)); 2239566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B,nrow)); 2249566063dSJacob 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)); 2259566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT,nrow)); 226*b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 2279566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B,MAT_REUSE_MATRIX,&BT)); 228c7a4214aSPierre Jolivet } else { 2299566063dSJacob Faibussowitsch PetscCall(MatTranspose(B,MAT_REUSE_MATRIX,&BT)); 230c7a4214aSPierre Jolivet } 2319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 233c7a4214aSPierre Jolivet } else { 234c7a4214aSPierre Jolivet for (PetscInt k=0; k<A->rmap->n; ++k) { /* block F from above */ 235c7a4214aSPierre 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); 236c7a4214aSPierre Jolivet } 237c7a4214aSPierre Jolivet } 238c7a4214aSPierre Jolivet } 239c7a4214aSPierre Jolivet } /* complete local diagonal block not in IS */ 240c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */ 241c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* rmap->rstart not in IS */ 242c7a4214aSPierre Jolivet } /* unsorted IS */ 243c7a4214aSPierre Jolivet } 244c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* different row and column IS */ 245c7a4214aSPierre Jolivet if (!flg) a->wrapper->copy_submatrix(nrow,m,idxr,idxc,ptr); /* reassemble everything */ 2469566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irow[i],&idxr)); 2479566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icol[i],&idxc)); 2489566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite((*submat)[i],&ptr)); 2499566063dSJacob Faibussowitsch PetscCall(MatScale((*submat)[i],a->s)); 250c7a4214aSPierre Jolivet } 251c7a4214aSPierre Jolivet PetscFunctionReturn(0); 252c7a4214aSPierre Jolivet } 253c7a4214aSPierre Jolivet 254c7a4214aSPierre Jolivet static PetscErrorCode MatDestroy_Htool(Mat A) 255c7a4214aSPierre Jolivet { 256c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 257c7a4214aSPierre Jolivet PetscContainer container; 258c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 259c7a4214aSPierre Jolivet 260c7a4214aSPierre Jolivet PetscFunctionBegin; 2619566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A,NULL)); 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",NULL)); 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",NULL)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",NULL)); 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",NULL)); 2669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",NULL)); 2679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",NULL)); 2689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",NULL)); 2699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",NULL)); 2709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",NULL)); 2719566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A,"KernelTranspose",(PetscObject*)&container)); 272c7a4214aSPierre Jolivet if (container) { /* created in MatTranspose_Htool() */ 2739566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container,(void**)&kernelt)); 2749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 2759566063dSJacob Faibussowitsch PetscCall(PetscFree(kernelt)); 2769566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&container)); 2779566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A,"KernelTranspose",NULL)); 278c7a4214aSPierre Jolivet } 279c7a4214aSPierre Jolivet if (a->gcoords_source != a->gcoords_target) { 2809566063dSJacob Faibussowitsch PetscCall(PetscFree(a->gcoords_source)); 281c7a4214aSPierre Jolivet } 2829566063dSJacob Faibussowitsch PetscCall(PetscFree(a->gcoords_target)); 2839566063dSJacob Faibussowitsch PetscCall(PetscFree2(a->work_source,a->work_target)); 284c7a4214aSPierre Jolivet delete a->wrapper; 285c7a4214aSPierre Jolivet delete a->hmatrix; 2869566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 287c7a4214aSPierre Jolivet PetscFunctionReturn(0); 288c7a4214aSPierre Jolivet } 289c7a4214aSPierre Jolivet 290c7a4214aSPierre Jolivet static PetscErrorCode MatView_Htool(Mat A,PetscViewer pv) 291c7a4214aSPierre Jolivet { 292c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 293c7a4214aSPierre Jolivet PetscBool flg; 294c7a4214aSPierre Jolivet 295c7a4214aSPierre Jolivet PetscFunctionBegin; 2969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pv,PETSCVIEWERASCII,&flg)); 297c7a4214aSPierre Jolivet if (flg) { 2989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"symmetry: %c\n",a->hmatrix->get_symmetry_type())); 299c7a4214aSPierre Jolivet if (PetscAbsScalar(a->s-1.0) > PETSC_MACHINE_EPSILON) { 300c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX) 3019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"scaling: %g+%gi\n",(double)PetscRealPart(a->s),(double)PetscImaginaryPart(a->s))); 302c7a4214aSPierre Jolivet #else 3039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"scaling: %g\n",(double)a->s)); 304c7a4214aSPierre Jolivet #endif 305c7a4214aSPierre Jolivet } 3069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"minimum cluster size: %" PetscInt_FMT "\n",a->bs[0])); 3079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"maximum block size: %" PetscInt_FMT "\n",a->bs[1])); 3089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"epsilon: %g\n",(double)a->epsilon)); 3099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"eta: %g\n",(double)a->eta)); 3109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"minimum target depth: %" PetscInt_FMT "\n",a->depth[0])); 3119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"minimum source depth: %" PetscInt_FMT "\n",a->depth[1])); 3129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"compressor: %s\n",MatHtoolCompressorTypes[a->compressor])); 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"clustering: %s\n",MatHtoolClusteringTypes[a->clustering])); 3149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"compression ratio: %s\n",a->hmatrix->get_infos("Compression_ratio").c_str())); 3159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"space saving: %s\n",a->hmatrix->get_infos("Space_saving").c_str())); 3169566063dSJacob 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())); 3179566063dSJacob 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())); 3189566063dSJacob 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())); 3199566063dSJacob 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())); 320c7a4214aSPierre Jolivet } 321c7a4214aSPierre Jolivet PetscFunctionReturn(0); 322c7a4214aSPierre Jolivet } 323c7a4214aSPierre Jolivet 324c7a4214aSPierre Jolivet static PetscErrorCode MatScale_Htool(Mat A,PetscScalar s) 325c7a4214aSPierre Jolivet { 326c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 327c7a4214aSPierre Jolivet 328c7a4214aSPierre Jolivet PetscFunctionBegin; 329c7a4214aSPierre Jolivet a->s *= s; 330c7a4214aSPierre Jolivet PetscFunctionReturn(0); 331c7a4214aSPierre Jolivet } 332c7a4214aSPierre Jolivet 333c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */ 334c7a4214aSPierre Jolivet static PetscErrorCode MatGetRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 335c7a4214aSPierre Jolivet { 336c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 337c7a4214aSPierre Jolivet PetscInt *idxc; 338c7a4214aSPierre Jolivet PetscBLASInt one = 1,bn; 339c7a4214aSPierre Jolivet 340c7a4214aSPierre Jolivet PetscFunctionBegin; 341c7a4214aSPierre Jolivet if (nz) *nz = A->cmap->N; 342c7a4214aSPierre Jolivet if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */ 3439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N,&idxc)); 344c7a4214aSPierre Jolivet for (PetscInt i=0; i<A->cmap->N; ++i) idxc[i] = i; 345c7a4214aSPierre Jolivet } 346c7a4214aSPierre Jolivet if (idx) *idx = idxc; 347c7a4214aSPierre Jolivet if (v) { 3489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N,v)); 349c7a4214aSPierre Jolivet if (a->wrapper) a->wrapper->copy_submatrix(1,A->cmap->N,&row,idxc,*v); 350cab00e0dSPierre Jolivet else reinterpret_cast<htool::VirtualGenerator<PetscScalar>*>(a->kernelctx)->copy_submatrix(1,A->cmap->N,&row,idxc,*v); 3519566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(A->cmap->N,&bn)); 352c7a4214aSPierre Jolivet PetscStackCallBLAS("BLASscal",BLASscal_(&bn,&a->s,*v,&one)); 353c7a4214aSPierre Jolivet } 354c7a4214aSPierre Jolivet if (!idx) { 3559566063dSJacob Faibussowitsch PetscCall(PetscFree(idxc)); 356c7a4214aSPierre Jolivet } 357c7a4214aSPierre Jolivet PetscFunctionReturn(0); 358c7a4214aSPierre Jolivet } 359c7a4214aSPierre Jolivet 360c7a4214aSPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 361c7a4214aSPierre Jolivet { 362c7a4214aSPierre Jolivet PetscFunctionBegin; 363c7a4214aSPierre Jolivet if (nz) *nz = 0; 364c7a4214aSPierre Jolivet if (idx) { 3659566063dSJacob Faibussowitsch PetscCall(PetscFree(*idx)); 366c7a4214aSPierre Jolivet } 367c7a4214aSPierre Jolivet if (v) { 3689566063dSJacob Faibussowitsch PetscCall(PetscFree(*v)); 369c7a4214aSPierre Jolivet } 370c7a4214aSPierre Jolivet PetscFunctionReturn(0); 371c7a4214aSPierre Jolivet } 372c7a4214aSPierre Jolivet 373c7a4214aSPierre Jolivet static PetscErrorCode MatSetFromOptions_Htool(PetscOptionItems *PetscOptionsObject,Mat A) 374c7a4214aSPierre Jolivet { 375c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 376c7a4214aSPierre Jolivet PetscInt n; 377c7a4214aSPierre Jolivet PetscBool flg; 378c7a4214aSPierre Jolivet 379c7a4214aSPierre Jolivet PetscFunctionBegin; 380d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Htool options"); 3819566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_min_cluster_size","Minimal leaf size in cluster tree",NULL,a->bs[0],a->bs,NULL)); 3829566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_max_block_size","Maximal number of coefficients in a dense block",NULL,a->bs[1],a->bs + 1,NULL)); 3839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_htool_epsilon","Relative error in Frobenius norm when approximating a block",NULL,a->epsilon,&a->epsilon,NULL)); 3849566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_htool_eta","Admissibility condition tolerance",NULL,a->eta,&a->eta,NULL)); 3859566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_min_target_depth","Minimal cluster tree depth associated with the rows",NULL,a->depth[0],a->depth,NULL)); 3869566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_min_source_depth","Minimal cluster tree depth associated with the columns",NULL,a->depth[1],a->depth + 1,NULL)); 387c7a4214aSPierre Jolivet n = 0; 388dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_compressor","Type of compression","MatHtoolCompressorType",MatHtoolCompressorTypes,PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes),MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA],&n,&flg)); 389c7a4214aSPierre Jolivet if (flg) a->compressor = MatHtoolCompressorType(n); 39098e73e17SPierre Jolivet n = 0; 391dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_clustering","Type of clustering","MatHtoolClusteringType",MatHtoolClusteringTypes,PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes),MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR],&n,&flg)); 39298e73e17SPierre Jolivet if (flg) a->clustering = MatHtoolClusteringType(n); 393d0609cedSBarry Smith PetscOptionsHeadEnd(); 394c7a4214aSPierre Jolivet PetscFunctionReturn(0); 395c7a4214aSPierre Jolivet } 396c7a4214aSPierre Jolivet 397c7a4214aSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A,MatAssemblyType type) 398c7a4214aSPierre Jolivet { 399c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 400c7a4214aSPierre Jolivet const PetscInt *ranges; 401c7a4214aSPierre Jolivet PetscInt *offset; 402c7a4214aSPierre Jolivet PetscMPIInt size; 403*b94d7dedSBarry Smith char S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'),uplo = S == 'N' ? 'N' : 'U'; 404cab00e0dSPierre Jolivet htool::VirtualGenerator<PetscScalar> *generator = nullptr; 40598e73e17SPierre Jolivet std::shared_ptr<htool::VirtualCluster> t,s = nullptr; 4063b9338faSPierre Jolivet std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor = nullptr; 407c7a4214aSPierre Jolivet 408c7a4214aSPierre Jolivet PetscFunctionBegin; 4099566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(HtoolCitation,&HtoolCite)); 410c7a4214aSPierre Jolivet delete a->wrapper; 411c7a4214aSPierre Jolivet delete a->hmatrix; 4129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 4139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2*size,&offset)); 4149566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A,&ranges)); 415c7a4214aSPierre Jolivet for (PetscInt i=0; i<size; ++i) { 416c7a4214aSPierre Jolivet offset[2*i] = ranges[i]; 417c7a4214aSPierre Jolivet offset[2*i+1] = ranges[i+1] - ranges[i]; 418c7a4214aSPierre Jolivet } 41998e73e17SPierre Jolivet switch (a->clustering) { 42098e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 42198e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 42298e73e17SPierre Jolivet break; 42398e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 42498e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 42598e73e17SPierre Jolivet break; 42698e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 42798e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 42898e73e17SPierre Jolivet break; 42998e73e17SPierre Jolivet default: 43098e73e17SPierre Jolivet t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 43198e73e17SPierre Jolivet } 432c7a4214aSPierre Jolivet t->set_minclustersize(a->bs[0]); 43398e73e17SPierre Jolivet t->build(A->rmap->N,a->gcoords_target,offset); 434c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx); 435c7a4214aSPierre Jolivet else { 436c7a4214aSPierre Jolivet a->wrapper = NULL; 437cab00e0dSPierre Jolivet generator = reinterpret_cast<htool::VirtualGenerator<PetscScalar>*>(a->kernelctx); 438c7a4214aSPierre Jolivet } 439c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 4409566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(A,&ranges)); 441c7a4214aSPierre Jolivet for (PetscInt i=0; i<size; ++i) { 442c7a4214aSPierre Jolivet offset[2*i] = ranges[i]; 443c7a4214aSPierre Jolivet offset[2*i+1] = ranges[i+1] - ranges[i]; 444c7a4214aSPierre Jolivet } 44598e73e17SPierre Jolivet switch (a->clustering) { 44698e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 44798e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 44898e73e17SPierre Jolivet break; 44998e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 45098e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 45198e73e17SPierre Jolivet break; 45298e73e17SPierre Jolivet case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 45398e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 45498e73e17SPierre Jolivet break; 45598e73e17SPierre Jolivet default: 45698e73e17SPierre Jolivet s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 45798e73e17SPierre Jolivet } 458c7a4214aSPierre Jolivet s->set_minclustersize(a->bs[0]); 45998e73e17SPierre Jolivet s->build(A->cmap->N,a->gcoords_source,offset); 460c7a4214aSPierre Jolivet S = uplo = 'N'; 461c7a4214aSPierre Jolivet } 4629566063dSJacob Faibussowitsch PetscCall(PetscFree(offset)); 463c7a4214aSPierre Jolivet switch (a->compressor) { 464c7a4214aSPierre Jolivet case MAT_HTOOL_COMPRESSOR_FULL_ACA: 4653b9338faSPierre Jolivet compressor = std::make_shared<htool::fullACA<PetscScalar>>(); 466c7a4214aSPierre Jolivet break; 467c7a4214aSPierre Jolivet case MAT_HTOOL_COMPRESSOR_SVD: 4683b9338faSPierre Jolivet compressor = std::make_shared<htool::SVD<PetscScalar>>(); 469c7a4214aSPierre Jolivet break; 470c7a4214aSPierre Jolivet default: 4713b9338faSPierre Jolivet compressor = std::make_shared<htool::sympartialACA<PetscScalar>>(); 472c7a4214aSPierre Jolivet } 4733b9338faSPierre Jolivet a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar>*>(new htool::HMatrix<PetscScalar>(t,s ? s : t,a->epsilon,a->eta,S,uplo)); 4743b9338faSPierre Jolivet a->hmatrix->set_compression(compressor); 475c7a4214aSPierre Jolivet a->hmatrix->set_maxblocksize(a->bs[1]); 476c7a4214aSPierre Jolivet a->hmatrix->set_mintargetdepth(a->depth[0]); 477c7a4214aSPierre Jolivet a->hmatrix->set_minsourcedepth(a->depth[1]); 4783b9338faSPierre Jolivet if (s) a->hmatrix->build(a->wrapper ? *a->wrapper : *generator,a->gcoords_target,a->gcoords_source); 4793b9338faSPierre Jolivet else a->hmatrix->build(a->wrapper ? *a->wrapper : *generator,a->gcoords_target); 480c7a4214aSPierre Jolivet PetscFunctionReturn(0); 481c7a4214aSPierre Jolivet } 482c7a4214aSPierre Jolivet 483c7a4214aSPierre Jolivet static PetscErrorCode MatProductNumeric_Htool(Mat C) 484c7a4214aSPierre Jolivet { 485c7a4214aSPierre Jolivet Mat_Product *product = C->product; 486c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)product->A->data; 487c7a4214aSPierre Jolivet const PetscScalar *in; 488c7a4214aSPierre Jolivet PetscScalar *out; 4898b8ddd11SPierre Jolivet PetscInt N,lda; 490c7a4214aSPierre Jolivet 491c7a4214aSPierre Jolivet PetscFunctionBegin; 492c7a4214aSPierre Jolivet MatCheckProduct(C,1); 4939566063dSJacob Faibussowitsch PetscCall(MatGetSize(C,NULL,&N)); 4949566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C,&lda)); 4955f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")",lda,C->rmap->n); 4969566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(product->B,&in)); 4979566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C,&out)); 4988b8ddd11SPierre Jolivet switch (product->type) { 4998b8ddd11SPierre Jolivet case MATPRODUCT_AB: 500c7a4214aSPierre Jolivet a->hmatrix->mvprod_local_to_local(in,out,N); 5018b8ddd11SPierre Jolivet break; 5028b8ddd11SPierre Jolivet case MATPRODUCT_AtB: 5038b8ddd11SPierre Jolivet a->hmatrix->mvprod_transp_local_to_local(in,out,N); 504c7a4214aSPierre Jolivet break; 505c7a4214aSPierre Jolivet default: 50698921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType %s is not supported",MatProductTypes[product->type]); 507c7a4214aSPierre Jolivet } 5089566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C,&out)); 5099566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(product->B,&in)); 5109566063dSJacob Faibussowitsch PetscCall(MatScale(C,a->s)); 511c7a4214aSPierre Jolivet PetscFunctionReturn(0); 512c7a4214aSPierre Jolivet } 513c7a4214aSPierre Jolivet 514c7a4214aSPierre Jolivet static PetscErrorCode MatProductSymbolic_Htool(Mat C) 515c7a4214aSPierre Jolivet { 516c7a4214aSPierre Jolivet Mat_Product *product = C->product; 517c7a4214aSPierre Jolivet Mat A,B; 518c7a4214aSPierre Jolivet PetscBool flg; 519c7a4214aSPierre Jolivet 520c7a4214aSPierre Jolivet PetscFunctionBegin; 521c7a4214aSPierre Jolivet MatCheckProduct(C,1); 522c7a4214aSPierre Jolivet A = product->A; 523c7a4214aSPierre Jolivet B = product->B; 5249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,"")); 5255f80ce2aSJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MatProduct_AB not supported for %s",((PetscObject)product->B)->type_name); 526c7a4214aSPierre Jolivet switch (product->type) { 527c7a4214aSPierre Jolivet case MATPRODUCT_AB: 528c7a4214aSPierre Jolivet if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) { 5299566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,A->rmap->n,B->cmap->n,A->rmap->N,B->cmap->N)); 530c7a4214aSPierre Jolivet } 5318b8ddd11SPierre Jolivet break; 5328b8ddd11SPierre Jolivet case MATPRODUCT_AtB: 5338b8ddd11SPierre Jolivet if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) { 5349566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,A->cmap->n,B->cmap->n,A->cmap->N,B->cmap->N)); 5358b8ddd11SPierre Jolivet } 5368b8ddd11SPierre Jolivet break; 5378b8ddd11SPierre Jolivet default: 53898921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"ProductType %s is not supported",MatProductTypes[product->type]); 5398b8ddd11SPierre Jolivet } 5409566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATDENSE)); 5419566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 5429566063dSJacob Faibussowitsch PetscCall(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 5439566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 5449566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 545c7a4214aSPierre Jolivet C->ops->productsymbolic = NULL; 546c7a4214aSPierre Jolivet C->ops->productnumeric = MatProductNumeric_Htool; 547c7a4214aSPierre Jolivet PetscFunctionReturn(0); 548c7a4214aSPierre Jolivet } 549c7a4214aSPierre Jolivet 550c7a4214aSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Htool(Mat C) 551c7a4214aSPierre Jolivet { 552c7a4214aSPierre Jolivet PetscFunctionBegin; 553c7a4214aSPierre Jolivet MatCheckProduct(C,1); 5548b8ddd11SPierre Jolivet if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool; 555c7a4214aSPierre Jolivet PetscFunctionReturn(0); 556c7a4214aSPierre Jolivet } 557c7a4214aSPierre Jolivet 55898e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix) 559c7a4214aSPierre Jolivet { 560c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 561c7a4214aSPierre Jolivet 562c7a4214aSPierre Jolivet PetscFunctionBegin; 563c7a4214aSPierre Jolivet *hmatrix = a->hmatrix; 564c7a4214aSPierre Jolivet PetscFunctionReturn(0); 565c7a4214aSPierre Jolivet } 566c7a4214aSPierre Jolivet 567c7a4214aSPierre Jolivet /*@C 568c7a4214aSPierre Jolivet MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a MATHTOOL. 569c7a4214aSPierre Jolivet 570c7a4214aSPierre Jolivet Input Parameter: 571c7a4214aSPierre Jolivet . A - hierarchical matrix 572c7a4214aSPierre Jolivet 573c7a4214aSPierre Jolivet Output Parameter: 574c7a4214aSPierre Jolivet . hmatrix - opaque pointer to a Htool virtual matrix 575c7a4214aSPierre Jolivet 576c7a4214aSPierre Jolivet Level: advanced 577c7a4214aSPierre Jolivet 578db781477SPatrick Sanan .seealso: `MATHTOOL` 579c7a4214aSPierre Jolivet @*/ 58098e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A,const htool::VirtualHMatrix<PetscScalar> **hmatrix) 581c7a4214aSPierre Jolivet { 582c7a4214aSPierre Jolivet PetscFunctionBegin; 583c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 584c7a4214aSPierre Jolivet PetscValidPointer(hmatrix,2); 585cac4c232SBarry Smith PetscTryMethod(A,"MatHtoolGetHierarchicalMat_C",(Mat,const htool::VirtualHMatrix<PetscScalar>**),(A,hmatrix)); 586c7a4214aSPierre Jolivet PetscFunctionReturn(0); 587c7a4214aSPierre Jolivet } 588c7a4214aSPierre Jolivet 589c7a4214aSPierre Jolivet static PetscErrorCode MatHtoolSetKernel_Htool(Mat A,MatHtoolKernel kernel,void *kernelctx) 590c7a4214aSPierre Jolivet { 591c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 592c7a4214aSPierre Jolivet 593c7a4214aSPierre Jolivet PetscFunctionBegin; 594c7a4214aSPierre Jolivet a->kernel = kernel; 595c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 596c7a4214aSPierre Jolivet delete a->wrapper; 597c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N,A->cmap->N,a->dim,a->kernel,a->kernelctx); 598c7a4214aSPierre Jolivet PetscFunctionReturn(0); 599c7a4214aSPierre Jolivet } 600c7a4214aSPierre Jolivet 601c7a4214aSPierre Jolivet /*@C 60298e73e17SPierre Jolivet MatHtoolSetKernel - Sets the kernel and context used for the assembly of a MATHTOOL. 603c7a4214aSPierre Jolivet 604c7a4214aSPierre Jolivet Input Parameters: 605c7a4214aSPierre Jolivet + A - hierarchical matrix 606c7a4214aSPierre Jolivet . kernel - computational kernel (or NULL) 607cab00e0dSPierre Jolivet - kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 608c7a4214aSPierre Jolivet 609c7a4214aSPierre Jolivet Level: advanced 610c7a4214aSPierre Jolivet 611db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatCreateHtoolFromKernel()` 612c7a4214aSPierre Jolivet @*/ 613c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat A,MatHtoolKernel kernel,void *kernelctx) 614c7a4214aSPierre Jolivet { 615c7a4214aSPierre Jolivet PetscFunctionBegin; 616c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 617c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel,2); 618c7a4214aSPierre Jolivet if (!kernel) PetscValidPointer(kernelctx,3); 619cac4c232SBarry Smith PetscTryMethod(A,"MatHtoolSetKernel_C",(Mat,MatHtoolKernel,void*),(A,kernel,kernelctx)); 620c7a4214aSPierre Jolivet PetscFunctionReturn(0); 621c7a4214aSPierre Jolivet } 622c7a4214aSPierre Jolivet 62398e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A,IS* is) 62498e73e17SPierre Jolivet { 62598e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 62698e73e17SPierre Jolivet std::vector<PetscInt> source; 62798e73e17SPierre Jolivet 62898e73e17SPierre Jolivet PetscFunctionBegin; 629a9918087SPierre Jolivet source = a->hmatrix->get_source_cluster()->get_local_perm(); 6309566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A),source.size(),source.data(),PETSC_COPY_VALUES,is)); 6319566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 63298e73e17SPierre Jolivet PetscFunctionReturn(0); 63398e73e17SPierre Jolivet } 63498e73e17SPierre Jolivet 63598e73e17SPierre Jolivet /*@C 63698e73e17SPierre Jolivet MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster. 63798e73e17SPierre Jolivet 63898e73e17SPierre Jolivet Input Parameter: 63998e73e17SPierre Jolivet . A - hierarchical matrix 64098e73e17SPierre Jolivet 64198e73e17SPierre Jolivet Output Parameter: 64298e73e17SPierre Jolivet . is - permutation 64398e73e17SPierre Jolivet 64498e73e17SPierre Jolivet Level: advanced 64598e73e17SPierre Jolivet 646db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()` 64798e73e17SPierre Jolivet @*/ 64898e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat A,IS* is) 64998e73e17SPierre Jolivet { 65098e73e17SPierre Jolivet PetscFunctionBegin; 65198e73e17SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 65298e73e17SPierre Jolivet if (!is) PetscValidPointer(is,2); 653cac4c232SBarry Smith PetscTryMethod(A,"MatHtoolGetPermutationSource_C",(Mat,IS*),(A,is)); 65498e73e17SPierre Jolivet PetscFunctionReturn(0); 65598e73e17SPierre Jolivet } 65698e73e17SPierre Jolivet 65798e73e17SPierre Jolivet static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A,IS* is) 65898e73e17SPierre Jolivet { 65998e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 66098e73e17SPierre Jolivet std::vector<PetscInt> target; 66198e73e17SPierre Jolivet 66298e73e17SPierre Jolivet PetscFunctionBegin; 663a9918087SPierre Jolivet target = a->hmatrix->get_target_cluster()->get_local_perm(); 6649566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A),target.size(),target.data(),PETSC_COPY_VALUES,is)); 6659566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 66698e73e17SPierre Jolivet PetscFunctionReturn(0); 66798e73e17SPierre Jolivet } 66898e73e17SPierre Jolivet 66998e73e17SPierre Jolivet /*@C 67098e73e17SPierre Jolivet MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster. 67198e73e17SPierre Jolivet 67298e73e17SPierre Jolivet Input Parameter: 67398e73e17SPierre Jolivet . A - hierarchical matrix 67498e73e17SPierre Jolivet 67598e73e17SPierre Jolivet Output Parameter: 67698e73e17SPierre Jolivet . is - permutation 67798e73e17SPierre Jolivet 67898e73e17SPierre Jolivet Level: advanced 67998e73e17SPierre Jolivet 680db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()` 68198e73e17SPierre Jolivet @*/ 68298e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat A,IS* is) 68398e73e17SPierre Jolivet { 68498e73e17SPierre Jolivet PetscFunctionBegin; 68598e73e17SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 68698e73e17SPierre Jolivet if (!is) PetscValidPointer(is,2); 687cac4c232SBarry Smith PetscTryMethod(A,"MatHtoolGetPermutationTarget_C",(Mat,IS*),(A,is)); 68898e73e17SPierre Jolivet PetscFunctionReturn(0); 68998e73e17SPierre Jolivet } 69098e73e17SPierre Jolivet 69198e73e17SPierre Jolivet static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A,PetscBool use) 69298e73e17SPierre Jolivet { 69398e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 69498e73e17SPierre Jolivet 69598e73e17SPierre Jolivet PetscFunctionBegin; 69698e73e17SPierre Jolivet a->hmatrix->set_use_permutation(use); 69798e73e17SPierre Jolivet PetscFunctionReturn(0); 69898e73e17SPierre Jolivet } 69998e73e17SPierre Jolivet 70098e73e17SPierre Jolivet /*@C 70198e73e17SPierre Jolivet MatHtoolUsePermutation - Sets whether MATHTOOL should permute input (resp. output) vectors following its internal source (resp. target) permutation. 70298e73e17SPierre Jolivet 70398e73e17SPierre Jolivet Input Parameters: 70498e73e17SPierre Jolivet + A - hierarchical matrix 70598e73e17SPierre Jolivet - use - Boolean value 70698e73e17SPierre Jolivet 70798e73e17SPierre Jolivet Level: advanced 70898e73e17SPierre Jolivet 709db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()` 71098e73e17SPierre Jolivet @*/ 71198e73e17SPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat A,PetscBool use) 71298e73e17SPierre Jolivet { 71398e73e17SPierre Jolivet PetscFunctionBegin; 71498e73e17SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 71598e73e17SPierre Jolivet PetscValidLogicalCollectiveBool(A,use,2); 716cac4c232SBarry Smith PetscTryMethod(A,"MatHtoolUsePermutation_C",(Mat,PetscBool),(A,use)); 71798e73e17SPierre Jolivet PetscFunctionReturn(0); 71898e73e17SPierre Jolivet } 71998e73e17SPierre Jolivet 720c7a4214aSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B) 721c7a4214aSPierre Jolivet { 722c7a4214aSPierre Jolivet Mat C; 723c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data; 724c7a4214aSPierre Jolivet PetscInt lda; 725c7a4214aSPierre Jolivet PetscScalar *array; 726c7a4214aSPierre Jolivet 727c7a4214aSPierre Jolivet PetscFunctionBegin; 728c7a4214aSPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 729c7a4214aSPierre Jolivet C = *B; 7305f80ce2aSJacob Faibussowitsch PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible dimensions"); 7319566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C,&lda)); 7325f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")",lda,C->rmap->n); 733c7a4214aSPierre Jolivet } else { 7349566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C)); 7359566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N)); 7369566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATDENSE)); 7379566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7389566063dSJacob Faibussowitsch PetscCall(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 739c7a4214aSPierre Jolivet } 7409566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 7419566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 7429566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C,&array)); 743c7a4214aSPierre Jolivet a->hmatrix->copy_local_dense_perm(array); 7449566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C,&array)); 7459566063dSJacob Faibussowitsch PetscCall(MatScale(C,a->s)); 746c7a4214aSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 7479566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A,&C)); 748c7a4214aSPierre Jolivet } else *B = C; 749c7a4214aSPierre Jolivet PetscFunctionReturn(0); 750c7a4214aSPierre Jolivet } 751c7a4214aSPierre Jolivet 75298e73e17SPierre Jolivet static PetscErrorCode GenEntriesTranspose(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *rows,const PetscInt *cols,PetscScalar *ptr,void *ctx) 753c7a4214aSPierre Jolivet { 754c7a4214aSPierre Jolivet MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose*)ctx; 75598e73e17SPierre Jolivet PetscScalar *tmp; 75698e73e17SPierre Jolivet 75798e73e17SPierre Jolivet PetscFunctionBegin; 75898e73e17SPierre Jolivet generator->kernel(sdim,N,M,cols,rows,ptr,generator->kernelctx); 7599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M*N,&tmp)); 7609566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp,ptr,M*N)); 76198e73e17SPierre Jolivet for (PetscInt i=0; i<M; ++i) { 76298e73e17SPierre Jolivet for (PetscInt j=0; j<N; ++j) ptr[i+j*M] = tmp[j+i*N]; 76398e73e17SPierre Jolivet } 7649566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 76598e73e17SPierre Jolivet PetscFunctionReturn(0); 766c7a4214aSPierre Jolivet } 767c7a4214aSPierre Jolivet 768c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */ 769c7a4214aSPierre Jolivet static PetscErrorCode MatTranspose_Htool(Mat A,MatReuse reuse,Mat *B) 770c7a4214aSPierre Jolivet { 771c7a4214aSPierre Jolivet Mat C; 772c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool*)A->data,*c; 773c7a4214aSPierre Jolivet PetscInt M = A->rmap->N,N = A->cmap->N,m = A->rmap->n,n = A->cmap->n; 774c7a4214aSPierre Jolivet PetscContainer container; 775c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 776c7a4214aSPierre Jolivet 777c7a4214aSPierre Jolivet PetscFunctionBegin; 7785f80ce2aSJacob Faibussowitsch PetscCheck(reuse != MAT_INPLACE_MATRIX,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatTranspose() with MAT_INPLACE_MATRIX not supported"); 779c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 7809566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C)); 7819566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,n,m,N,M)); 7829566063dSJacob Faibussowitsch PetscCall(MatSetType(C,((PetscObject)A)->type_name)); 7839566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7849566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)C),&container)); 7859566063dSJacob Faibussowitsch PetscCall(PetscNew(&kernelt)); 7869566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(container,kernelt)); 7879566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)C,"KernelTranspose",(PetscObject)container)); 788c7a4214aSPierre Jolivet } else { 789c7a4214aSPierre Jolivet C = *B; 7909566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)C,"KernelTranspose",(PetscObject*)&container)); 7915f80ce2aSJacob Faibussowitsch PetscCheck(container,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatTranspose() with MAT_INITIAL_MATRIX first"); 7929566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container,(void**)&kernelt)); 793c7a4214aSPierre Jolivet } 794c7a4214aSPierre Jolivet c = (Mat_Htool*)C->data; 795c7a4214aSPierre Jolivet c->dim = a->dim; 796c7a4214aSPierre Jolivet c->s = a->s; 79798e73e17SPierre Jolivet c->kernel = GenEntriesTranspose; 798c7a4214aSPierre Jolivet if (kernelt->A != A) { 7999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 800c7a4214aSPierre Jolivet kernelt->A = A; 8019566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 802c7a4214aSPierre Jolivet } 803c7a4214aSPierre Jolivet kernelt->kernel = a->kernel; 804c7a4214aSPierre Jolivet kernelt->kernelctx = a->kernelctx; 805c7a4214aSPierre Jolivet c->kernelctx = kernelt; 806c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 8079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N*c->dim,&c->gcoords_target)); 8089566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_target,a->gcoords_source,N*c->dim)); 809c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 8109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M*c->dim,&c->gcoords_source)); 8119566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_source,a->gcoords_target,M*c->dim)); 812c7a4214aSPierre Jolivet } else c->gcoords_source = c->gcoords_target; 8139566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(M,&c->work_source,N,&c->work_target)); 814c7a4214aSPierre Jolivet } 8159566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 8169566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 817c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) *B = C; 818c7a4214aSPierre Jolivet PetscFunctionReturn(0); 819c7a4214aSPierre Jolivet } 820c7a4214aSPierre Jolivet 821c7a4214aSPierre Jolivet /*@C 822c7a4214aSPierre Jolivet MatCreateHtoolFromKernel - Creates a MATHTOOL from a user-supplied kernel. 823c7a4214aSPierre Jolivet 824c7a4214aSPierre Jolivet Input Parameters: 825c7a4214aSPierre Jolivet + comm - MPI communicator 826c7a4214aSPierre Jolivet . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 827c7a4214aSPierre Jolivet . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 828c7a4214aSPierre Jolivet . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 829c7a4214aSPierre Jolivet . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 830c7a4214aSPierre Jolivet . spacedim - dimension of the space coordinates 831c7a4214aSPierre Jolivet . coords_target - coordinates of the target 832c7a4214aSPierre Jolivet . coords_source - coordinates of the source 833c7a4214aSPierre Jolivet . kernel - computational kernel (or NULL) 834cab00e0dSPierre Jolivet - kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 835c7a4214aSPierre Jolivet 836c7a4214aSPierre Jolivet Output Parameter: 837c7a4214aSPierre Jolivet . B - matrix 838c7a4214aSPierre Jolivet 839c7a4214aSPierre Jolivet Options Database Keys: 840c7a4214aSPierre Jolivet + -mat_htool_min_cluster_size <PetscInt> - minimal leaf size in cluster tree 841c7a4214aSPierre Jolivet . -mat_htool_max_block_size <PetscInt> - maximal number of coefficients in a dense block 842c7a4214aSPierre Jolivet . -mat_htool_epsilon <PetscReal> - relative error in Frobenius norm when approximating a block 843c7a4214aSPierre Jolivet . -mat_htool_eta <PetscReal> - admissibility condition tolerance 844c7a4214aSPierre Jolivet . -mat_htool_min_target_depth <PetscInt> - minimal cluster tree depth associated with the rows 845c7a4214aSPierre Jolivet . -mat_htool_min_source_depth <PetscInt> - minimal cluster tree depth associated with the columns 84698e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression 84798e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering 848c7a4214aSPierre Jolivet 849c7a4214aSPierre Jolivet Level: intermediate 850c7a4214aSPierre Jolivet 851db781477SPatrick Sanan .seealso: `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()` 852c7a4214aSPierre Jolivet @*/ 853c7a4214aSPierre 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) 854c7a4214aSPierre Jolivet { 855c7a4214aSPierre Jolivet Mat A; 856c7a4214aSPierre Jolivet Mat_Htool *a; 857c7a4214aSPierre Jolivet 858c7a4214aSPierre Jolivet PetscFunctionBegin; 8599566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&A)); 860c7a4214aSPierre Jolivet PetscValidLogicalCollectiveInt(A,spacedim,6); 861c7a4214aSPierre Jolivet PetscValidRealPointer(coords_target,7); 862c7a4214aSPierre Jolivet PetscValidRealPointer(coords_source,8); 863c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel,9); 864c7a4214aSPierre Jolivet if (!kernel) PetscValidPointer(kernelctx,10); 8659566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,m,n,M,N)); 8669566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATHTOOL)); 8679566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 868c7a4214aSPierre Jolivet a = (Mat_Htool*)A->data; 869c7a4214aSPierre Jolivet a->dim = spacedim; 870c7a4214aSPierre Jolivet a->s = 1.0; 871c7a4214aSPierre Jolivet a->kernel = kernel; 872c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 8739566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->rmap->N*spacedim,&a->gcoords_target)); 8749566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_target+A->rmap->rstart*spacedim,coords_target,A->rmap->n*spacedim)); 8751c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE,a->gcoords_target,A->rmap->N*spacedim,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)A))); /* global target coordinates */ 876c7a4214aSPierre Jolivet if (coords_target != coords_source) { 8779566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->cmap->N*spacedim,&a->gcoords_source)); 8789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_source+A->cmap->rstart*spacedim,coords_source,A->cmap->n*spacedim)); 8791c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE,a->gcoords_source,A->cmap->N*spacedim,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)A))); /* global source coordinates */ 880c7a4214aSPierre Jolivet } else a->gcoords_source = a->gcoords_target; 8819566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(A->cmap->N,&a->work_source,A->rmap->N,&a->work_target)); 882c7a4214aSPierre Jolivet *B = A; 883c7a4214aSPierre Jolivet PetscFunctionReturn(0); 884c7a4214aSPierre Jolivet } 885c7a4214aSPierre Jolivet 886c7a4214aSPierre Jolivet /*MC 887c7a4214aSPierre Jolivet MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package. 888c7a4214aSPierre Jolivet 889c7a4214aSPierre Jolivet Use ./configure --download-htool to install PETSc to use Htool. 890c7a4214aSPierre Jolivet 891c7a4214aSPierre Jolivet Options Database Keys: 892e85a1b63SPatrick Sanan . -mat_type htool - matrix type to "htool" during a call to MatSetFromOptions() 893c7a4214aSPierre Jolivet 894c7a4214aSPierre Jolivet Level: beginner 895c7a4214aSPierre Jolivet 896db781477SPatrick Sanan .seealso: `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()` 897c7a4214aSPierre Jolivet M*/ 898c7a4214aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A) 899c7a4214aSPierre Jolivet { 900c7a4214aSPierre Jolivet Mat_Htool *a; 901c7a4214aSPierre Jolivet 902c7a4214aSPierre Jolivet PetscFunctionBegin; 9039566063dSJacob Faibussowitsch PetscCall(PetscNewLog(A,&a)); 904c7a4214aSPierre Jolivet A->data = (void*)a; 9059566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATHTOOL)); 9069566063dSJacob Faibussowitsch PetscCall(PetscMemzero(A->ops,sizeof(struct _MatOps))); 907c7a4214aSPierre Jolivet A->ops->getdiagonal = MatGetDiagonal_Htool; 908c7a4214aSPierre Jolivet A->ops->getdiagonalblock = MatGetDiagonalBlock_Htool; 909c7a4214aSPierre Jolivet A->ops->mult = MatMult_Htool; 910c7a4214aSPierre Jolivet A->ops->multadd = MatMultAdd_Htool; 9118b8ddd11SPierre Jolivet A->ops->multtranspose = MatMultTranspose_Htool; 9128b8ddd11SPierre Jolivet if (!PetscDefined(USE_COMPLEX)) A->ops->multhermitiantranspose = MatMultTranspose_Htool; 913c7a4214aSPierre Jolivet A->ops->increaseoverlap = MatIncreaseOverlap_Htool; 914c7a4214aSPierre Jolivet A->ops->createsubmatrices = MatCreateSubMatrices_Htool; 915c7a4214aSPierre Jolivet A->ops->transpose = MatTranspose_Htool; 916c7a4214aSPierre Jolivet A->ops->destroy = MatDestroy_Htool; 917c7a4214aSPierre Jolivet A->ops->view = MatView_Htool; 918c7a4214aSPierre Jolivet A->ops->setfromoptions = MatSetFromOptions_Htool; 919c7a4214aSPierre Jolivet A->ops->scale = MatScale_Htool; 920c7a4214aSPierre Jolivet A->ops->getrow = MatGetRow_Htool; 921c7a4214aSPierre Jolivet A->ops->restorerow = MatRestoreRow_Htool; 922c7a4214aSPierre Jolivet A->ops->assemblyend = MatAssemblyEnd_Htool; 923c7a4214aSPierre Jolivet a->dim = 0; 924c7a4214aSPierre Jolivet a->gcoords_target = NULL; 925c7a4214aSPierre Jolivet a->gcoords_source = NULL; 926c7a4214aSPierre Jolivet a->s = 1.0; 927c7a4214aSPierre Jolivet a->bs[0] = 10; 928c7a4214aSPierre Jolivet a->bs[1] = 1000000; 929c7a4214aSPierre Jolivet a->epsilon = PetscSqrtReal(PETSC_SMALL); 930c7a4214aSPierre Jolivet a->eta = 10.0; 931c7a4214aSPierre Jolivet a->depth[0] = 0; 932c7a4214aSPierre Jolivet a->depth[1] = 0; 933c7a4214aSPierre Jolivet a->compressor = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA; 9349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_seqdense_C",MatProductSetFromOptions_Htool)); 9359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_htool_mpidense_C",MatProductSetFromOptions_Htool)); 9369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_seqdense_C",MatConvert_Htool_Dense)); 9379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_htool_mpidense_C",MatConvert_Htool_Dense)); 9389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetHierarchicalMat_C",MatHtoolGetHierarchicalMat_Htool)); 9399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolSetKernel_C",MatHtoolSetKernel_Htool)); 9409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationSource_C",MatHtoolGetPermutationSource_Htool)); 9419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolGetPermutationTarget_C",MatHtoolGetPermutationTarget_Htool)); 9429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatHtoolUsePermutation_C",MatHtoolUsePermutation_Htool)); 943c7a4214aSPierre Jolivet PetscFunctionReturn(0); 944c7a4214aSPierre Jolivet } 945