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