xref: /petsc/src/ksp/pc/impls/h2opus/pch2opus.c (revision 8db5126349b744ef72bea061699bcc188cb355a4)
153022affSStefano Zampini #include <petsc/private/pcimpl.h>
253022affSStefano Zampini #include <petsc/private/matimpl.h>
353022affSStefano Zampini #include <h2opusconf.h>
453022affSStefano Zampini 
553022affSStefano Zampini /* Use GPU only if H2OPUS is configured for GPU */
653022affSStefano Zampini #if defined(PETSC_HAVE_CUDA) && defined(H2OPUS_USE_GPU)
753022affSStefano Zampini #define PETSC_H2OPUS_USE_GPU
853022affSStefano Zampini #endif
953022affSStefano Zampini 
1053022affSStefano Zampini typedef struct {
1153022affSStefano Zampini   Mat         A;
1253022affSStefano Zampini   Mat         M;
1353022affSStefano Zampini   PetscScalar s0;
1453022affSStefano Zampini 
1553022affSStefano Zampini   /* sampler for Newton-Schultz */
1653022affSStefano Zampini   Mat      S;
1753022affSStefano Zampini   PetscInt hyperorder;
1853022affSStefano Zampini   Vec      wns[4];
1953022affSStefano Zampini   Mat      wnsmat[4];
2053022affSStefano Zampini 
2153022affSStefano Zampini   /* convergence testing */
2253022affSStefano Zampini   Mat       T;
2353022affSStefano Zampini   Vec       w;
2453022affSStefano Zampini 
2553022affSStefano Zampini   /* Support for PCSetCoordinates */
2653022affSStefano Zampini   PetscInt  sdim;
2753022affSStefano Zampini   PetscInt  nlocc;
2853022affSStefano Zampini   PetscReal *coords;
2953022affSStefano Zampini 
3053022affSStefano Zampini   /* Newton-Schultz customization */
3153022affSStefano Zampini   PetscInt  maxits;
3253022affSStefano Zampini   PetscReal rtol,atol;
3353022affSStefano Zampini   PetscBool monitor;
3453022affSStefano Zampini   PetscBool useapproximatenorms;
3553022affSStefano Zampini   NormType  normtype;
3653022affSStefano Zampini 
3753022affSStefano Zampini   /* Used when pmat != MATH2OPUS */
3853022affSStefano Zampini   PetscReal eta;
3953022affSStefano Zampini   PetscInt  leafsize;
4053022affSStefano Zampini   PetscInt  max_rank;
4153022affSStefano Zampini   PetscInt  bs;
4253022affSStefano Zampini   PetscReal mrtol;
4353022affSStefano Zampini 
44*8db51263SStefano Zampini   /* CPU/GPU */
45*8db51263SStefano Zampini   PetscBool forcecpu;
4653022affSStefano Zampini   PetscBool boundtocpu;
4753022affSStefano Zampini } PC_H2OPUS;
4853022affSStefano Zampini 
4953022affSStefano Zampini PETSC_EXTERN PetscErrorCode MatNorm_H2OPUS(Mat,NormType,PetscReal*);
5053022affSStefano Zampini 
5153022affSStefano Zampini static PetscErrorCode PCReset_H2OPUS(PC pc)
5253022affSStefano Zampini {
5353022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
5453022affSStefano Zampini 
5553022affSStefano Zampini   PetscFunctionBegin;
5653022affSStefano Zampini   pch2opus->sdim  = 0;
5753022affSStefano Zampini   pch2opus->nlocc = 0;
589566063dSJacob Faibussowitsch   PetscCall(PetscFree(pch2opus->coords));
599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->A));
609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->M));
619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->T));
629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->w));
639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->S));
649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[0]));
659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[1]));
669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[2]));
679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[3]));
689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[0]));
699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[1]));
709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[2]));
719566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[3]));
7253022affSStefano Zampini   PetscFunctionReturn(0);
7353022affSStefano Zampini }
7453022affSStefano Zampini 
7553022affSStefano Zampini static PetscErrorCode PCSetCoordinates_H2OPUS(PC pc, PetscInt sdim, PetscInt nlocc, PetscReal *coords)
7653022affSStefano Zampini {
7753022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
7853022affSStefano Zampini   PetscBool  reset    = PETSC_TRUE;
7953022affSStefano Zampini 
8053022affSStefano Zampini   PetscFunctionBegin;
8153022affSStefano Zampini   if (pch2opus->sdim && sdim == pch2opus->sdim && nlocc == pch2opus->nlocc) {
829566063dSJacob Faibussowitsch     PetscCall(PetscArraycmp(pch2opus->coords,coords,sdim*nlocc,&reset));
8353022affSStefano Zampini     reset = (PetscBool)!reset;
8453022affSStefano Zampini   }
851c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(MPI_IN_PLACE,&reset,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc)));
8653022affSStefano Zampini   if (reset) {
879566063dSJacob Faibussowitsch     PetscCall(PCReset_H2OPUS(pc));
889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(sdim*nlocc,&pch2opus->coords));
899566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pch2opus->coords,coords,sdim*nlocc));
9053022affSStefano Zampini     pch2opus->sdim  = sdim;
9153022affSStefano Zampini     pch2opus->nlocc = nlocc;
9253022affSStefano Zampini   }
9353022affSStefano Zampini   PetscFunctionReturn(0);
9453022affSStefano Zampini }
9553022affSStefano Zampini 
9653022affSStefano Zampini static PetscErrorCode PCDestroy_H2OPUS(PC pc)
9753022affSStefano Zampini {
9853022affSStefano Zampini   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(PCReset_H2OPUS(pc));
1009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL));
1019566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
10253022affSStefano Zampini   PetscFunctionReturn(0);
10353022affSStefano Zampini }
10453022affSStefano Zampini 
10553022affSStefano Zampini static PetscErrorCode PCSetFromOptions_H2OPUS(PetscOptionItems *PetscOptionsObject,PC pc)
10653022affSStefano Zampini {
10753022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
10853022affSStefano Zampini 
10953022affSStefano Zampini   PetscFunctionBegin;
110d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"H2OPUS options");
1119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_h2opus_maxits","Maximum number of iterations for Newton-Schultz",NULL,pch2opus->maxits,&pch2opus->maxits,NULL));
1129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_h2opus_monitor","Monitor Newton-Schultz convergence",NULL,pch2opus->monitor,&pch2opus->monitor,NULL));
1139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_h2opus_atol","Absolute tolerance",NULL,pch2opus->atol,&pch2opus->atol,NULL));
1149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_h2opus_rtol","Relative tolerance",NULL,pch2opus->rtol,&pch2opus->rtol,NULL));
1159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_h2opus_norm_type","Norm type for convergence monitoring",NULL,NormTypes,(PetscEnum)pch2opus->normtype,(PetscEnum*)&pch2opus->normtype,NULL));
1169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_h2opus_hyperorder","Hyper power order of sampling",NULL,pch2opus->hyperorder,&pch2opus->hyperorder,NULL));
1179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_h2opus_leafsize","Leaf size when constructed from kernel",NULL,pch2opus->leafsize,&pch2opus->leafsize,NULL));
1189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_h2opus_eta","Admissibility condition tolerance",NULL,pch2opus->eta,&pch2opus->eta,NULL));
1199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_h2opus_maxrank","Maximum rank when constructed from matvecs",NULL,pch2opus->max_rank,&pch2opus->max_rank,NULL));
1209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_h2opus_samples","Number of samples to be taken concurrently when constructing from matvecs",NULL,pch2opus->bs,&pch2opus->bs,NULL));
1219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_h2opus_mrtol","Relative tolerance for construction from sampling",NULL,pch2opus->mrtol,&pch2opus->mrtol,NULL));
122*8db51263SStefano Zampini   PetscCall(PetscOptionsBool("-pc_h2opus_forcecpu","Force construction on CPU",NULL,pch2opus->forcecpu,&pch2opus->forcecpu,NULL));
123d0609cedSBarry Smith   PetscOptionsHeadEnd();
12453022affSStefano Zampini   PetscFunctionReturn(0);
12553022affSStefano Zampini }
12653022affSStefano Zampini 
12753022affSStefano Zampini typedef struct {
12853022affSStefano Zampini   Mat A;
12953022affSStefano Zampini   Mat M;
13053022affSStefano Zampini   Vec w;
13153022affSStefano Zampini } AAtCtx;
13253022affSStefano Zampini 
13353022affSStefano Zampini static PetscErrorCode MatMult_AAt(Mat A, Vec x, Vec y)
13453022affSStefano Zampini {
13553022affSStefano Zampini   AAtCtx *aat;
13653022affSStefano Zampini 
13753022affSStefano Zampini   PetscFunctionBegin;
1389566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A,(void*)&aat));
1399566063dSJacob Faibussowitsch   /* PetscCall(MatMultTranspose(aat->M,x,aat->w)); */
1409566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(aat->A,x,aat->w));
1419566063dSJacob Faibussowitsch   PetscCall(MatMult(aat->A,aat->w,y));
14253022affSStefano Zampini   PetscFunctionReturn(0);
14353022affSStefano Zampini }
14453022affSStefano Zampini 
14553022affSStefano Zampini static PetscErrorCode PCH2OpusSetUpInit(PC pc)
14653022affSStefano Zampini {
14753022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
14853022affSStefano Zampini   Mat        A        = pc->useAmat ? pc->mat : pc->pmat, AAt;
14953022affSStefano Zampini   PetscInt   M,m;
15053022affSStefano Zampini   VecType    vtype;
15153022affSStefano Zampini   PetscReal  n;
15253022affSStefano Zampini   AAtCtx     aat;
15353022affSStefano Zampini 
15453022affSStefano Zampini   PetscFunctionBegin;
15553022affSStefano Zampini   aat.A = A;
15653022affSStefano Zampini   aat.M = pch2opus->M; /* unused so far */
1579566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,NULL,&aat.w));
1589566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&M,NULL));
1599566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,NULL));
1609566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)A),m,m,M,M,&aat,&AAt));
1619566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(AAt,pch2opus->boundtocpu));
1629566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(AAt,MATOP_MULT,(void (*)(void))MatMult_AAt));
1639566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(AAt,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMult_AAt));
1649566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(AAt,MATOP_NORM,(void (*)(void))MatNorm_H2OPUS));
1659566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A,&vtype));
1669566063dSJacob Faibussowitsch   PetscCall(MatShellSetVecType(AAt,vtype));
1679566063dSJacob Faibussowitsch   PetscCall(MatNorm(AAt,NORM_1,&n));
16853022affSStefano Zampini   pch2opus->s0 = 1./n;
1699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aat.w));
1709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&AAt));
17153022affSStefano Zampini   PetscFunctionReturn(0);
17253022affSStefano Zampini }
17353022affSStefano Zampini 
17453022affSStefano Zampini static PetscErrorCode PCApplyKernel_H2OPUS(PC pc, Vec x, Vec y, PetscBool t)
17553022affSStefano Zampini {
17653022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
17753022affSStefano Zampini 
17853022affSStefano Zampini   PetscFunctionBegin;
1799566063dSJacob Faibussowitsch   if (t) PetscCall(MatMultTranspose(pch2opus->M,x,y));
1809566063dSJacob Faibussowitsch   else  PetscCall(MatMult(pch2opus->M,x,y));
18153022affSStefano Zampini   PetscFunctionReturn(0);
18253022affSStefano Zampini }
18353022affSStefano Zampini 
18453022affSStefano Zampini static PetscErrorCode PCApplyMatKernel_H2OPUS(PC pc, Mat X, Mat Y, PetscBool t)
18553022affSStefano Zampini {
18653022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
18753022affSStefano Zampini 
18853022affSStefano Zampini   PetscFunctionBegin;
1899566063dSJacob Faibussowitsch   if (t) PetscCall(MatTransposeMatMult(pch2opus->M,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y));
1909566063dSJacob Faibussowitsch   else   PetscCall(MatMatMult(pch2opus->M,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y));
19153022affSStefano Zampini   PetscFunctionReturn(0);
19253022affSStefano Zampini }
19353022affSStefano Zampini 
19453022affSStefano Zampini static PetscErrorCode PCApplyMat_H2OPUS(PC pc, Mat X, Mat Y)
19553022affSStefano Zampini {
19653022affSStefano Zampini   PetscFunctionBegin;
1979566063dSJacob Faibussowitsch   PetscCall(PCApplyMatKernel_H2OPUS(pc,X,Y,PETSC_FALSE));
19853022affSStefano Zampini   PetscFunctionReturn(0);
19953022affSStefano Zampini }
20053022affSStefano Zampini 
20153022affSStefano Zampini static PetscErrorCode PCApplyTransposeMat_H2OPUS(PC pc, Mat X, Mat Y)
20253022affSStefano Zampini {
20353022affSStefano Zampini   PetscFunctionBegin;
2049566063dSJacob Faibussowitsch   PetscCall(PCApplyMatKernel_H2OPUS(pc,X,Y,PETSC_TRUE));
20553022affSStefano Zampini   PetscFunctionReturn(0);
20653022affSStefano Zampini }
20753022affSStefano Zampini 
20853022affSStefano Zampini static PetscErrorCode PCApply_H2OPUS(PC pc, Vec x, Vec y)
20953022affSStefano Zampini {
21053022affSStefano Zampini   PetscFunctionBegin;
2119566063dSJacob Faibussowitsch   PetscCall(PCApplyKernel_H2OPUS(pc,x,y,PETSC_FALSE));
21253022affSStefano Zampini   PetscFunctionReturn(0);
21353022affSStefano Zampini }
21453022affSStefano Zampini 
21553022affSStefano Zampini static PetscErrorCode PCApplyTranspose_H2OPUS(PC pc, Vec x, Vec y)
21653022affSStefano Zampini {
21753022affSStefano Zampini   PetscFunctionBegin;
2189566063dSJacob Faibussowitsch   PetscCall(PCApplyKernel_H2OPUS(pc,x,y,PETSC_TRUE));
21953022affSStefano Zampini   PetscFunctionReturn(0);
22053022affSStefano Zampini }
22153022affSStefano Zampini 
22253022affSStefano Zampini /* used to test the norm of (M^-1 A - I) */
22353022affSStefano Zampini static PetscErrorCode MatMultKernel_MAmI(Mat M, Vec x, Vec y, PetscBool t)
22453022affSStefano Zampini {
22553022affSStefano Zampini   PC         pc;
22653022affSStefano Zampini   Mat        A;
22753022affSStefano Zampini   PC_H2OPUS *pch2opus;
22853022affSStefano Zampini   PetscBool  sideleft = PETSC_TRUE;
22953022affSStefano Zampini 
23053022affSStefano Zampini   PetscFunctionBegin;
2319566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M,(void*)&pc));
23253022affSStefano Zampini   pch2opus = (PC_H2OPUS*)pc->data;
2339566063dSJacob Faibussowitsch   if (!pch2opus->w) PetscCall(MatCreateVecs(pch2opus->M,&pch2opus->w,NULL));
23453022affSStefano Zampini   A = pch2opus->A;
2359566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->w,pch2opus->boundtocpu));
23653022affSStefano Zampini   if (t) {
23753022affSStefano Zampini     if (sideleft) {
2389566063dSJacob Faibussowitsch       PetscCall(PCApplyTranspose_H2OPUS(pc,x,pch2opus->w));
2399566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A,pch2opus->w,y));
24053022affSStefano Zampini     } else {
2419566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A,x,pch2opus->w));
2429566063dSJacob Faibussowitsch       PetscCall(PCApplyTranspose_H2OPUS(pc,pch2opus->w,y));
24353022affSStefano Zampini     }
24453022affSStefano Zampini   } else {
24553022affSStefano Zampini     if (sideleft) {
2469566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,pch2opus->w));
2479566063dSJacob Faibussowitsch       PetscCall(PCApply_H2OPUS(pc,pch2opus->w,y));
24853022affSStefano Zampini     } else {
2499566063dSJacob Faibussowitsch       PetscCall(PCApply_H2OPUS(pc,x,pch2opus->w));
2509566063dSJacob Faibussowitsch       PetscCall(MatMult(A,pch2opus->w,y));
25153022affSStefano Zampini     }
25253022affSStefano Zampini   }
253*8db51263SStefano Zampini   PetscCall(VecAXPY(y,-1.0,x));
25453022affSStefano Zampini   PetscFunctionReturn(0);
25553022affSStefano Zampini }
25653022affSStefano Zampini 
25753022affSStefano Zampini static PetscErrorCode MatMult_MAmI(Mat A, Vec x, Vec y)
25853022affSStefano Zampini {
25953022affSStefano Zampini   PetscFunctionBegin;
2609566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_MAmI(A,x,y,PETSC_FALSE));
26153022affSStefano Zampini   PetscFunctionReturn(0);
26253022affSStefano Zampini }
26353022affSStefano Zampini 
26453022affSStefano Zampini static PetscErrorCode MatMultTranspose_MAmI(Mat A, Vec x, Vec y)
26553022affSStefano Zampini {
26653022affSStefano Zampini   PetscFunctionBegin;
2679566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_MAmI(A,x,y,PETSC_TRUE));
26853022affSStefano Zampini   PetscFunctionReturn(0);
26953022affSStefano Zampini }
27053022affSStefano Zampini 
27153022affSStefano Zampini /* HyperPower kernel:
27253022affSStefano Zampini Y = R = x
27353022affSStefano Zampini for i = 1 . . . l - 1 do
27453022affSStefano Zampini   R = (I - A * Xk) * R
27553022affSStefano Zampini   Y = Y + R
27653022affSStefano Zampini Y = Xk * Y
27753022affSStefano Zampini */
27853022affSStefano Zampini static PetscErrorCode MatMultKernel_Hyper(Mat M, Vec x, Vec y, PetscBool t)
27953022affSStefano Zampini {
28053022affSStefano Zampini   PC         pc;
28153022affSStefano Zampini   Mat        A;
28253022affSStefano Zampini   PC_H2OPUS *pch2opus;
28353022affSStefano Zampini 
28453022affSStefano Zampini   PetscFunctionBegin;
2859566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M,(void*)&pc));
28653022affSStefano Zampini   pch2opus = (PC_H2OPUS*)pc->data;
28753022affSStefano Zampini   A = pch2opus->A;
2889566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pch2opus->M,pch2opus->wns[0] ? NULL : &pch2opus->wns[0],pch2opus->wns[1] ? NULL : &pch2opus->wns[1]));
2899566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pch2opus->M,pch2opus->wns[2] ? NULL : &pch2opus->wns[2],pch2opus->wns[3] ? NULL : &pch2opus->wns[3]));
2909566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[0],pch2opus->boundtocpu));
2919566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[1],pch2opus->boundtocpu));
2929566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[2],pch2opus->boundtocpu));
2939566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[3],pch2opus->boundtocpu));
2949566063dSJacob Faibussowitsch   PetscCall(VecCopy(x,pch2opus->wns[0]));
2959566063dSJacob Faibussowitsch   PetscCall(VecCopy(x,pch2opus->wns[3]));
29653022affSStefano Zampini   if (t) {
2975f80ce2aSJacob Faibussowitsch     for (PetscInt i=0;i<pch2opus->hyperorder-1;i++) {
2989566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A,pch2opus->wns[0],pch2opus->wns[1]));
2999566063dSJacob Faibussowitsch       PetscCall(PCApplyTranspose_H2OPUS(pc,pch2opus->wns[1],pch2opus->wns[2]));
3009566063dSJacob Faibussowitsch       PetscCall(VecAXPY(pch2opus->wns[0],-1.,pch2opus->wns[2]));
3019566063dSJacob Faibussowitsch       PetscCall(VecAXPY(pch2opus->wns[3], 1.,pch2opus->wns[0]));
30253022affSStefano Zampini     }
3039566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose_H2OPUS(pc,pch2opus->wns[3],y));
30453022affSStefano Zampini   } else {
3055f80ce2aSJacob Faibussowitsch     for (PetscInt i=0;i<pch2opus->hyperorder-1;i++) {
3069566063dSJacob Faibussowitsch       PetscCall(PCApply_H2OPUS(pc,pch2opus->wns[0],pch2opus->wns[1]));
3079566063dSJacob Faibussowitsch       PetscCall(MatMult(A,pch2opus->wns[1],pch2opus->wns[2]));
3089566063dSJacob Faibussowitsch       PetscCall(VecAXPY(pch2opus->wns[0],-1.,pch2opus->wns[2]));
3099566063dSJacob Faibussowitsch       PetscCall(VecAXPY(pch2opus->wns[3], 1.,pch2opus->wns[0]));
31053022affSStefano Zampini     }
3119566063dSJacob Faibussowitsch     PetscCall(PCApply_H2OPUS(pc,pch2opus->wns[3],y));
31253022affSStefano Zampini   }
31353022affSStefano Zampini   PetscFunctionReturn(0);
31453022affSStefano Zampini }
31553022affSStefano Zampini 
31653022affSStefano Zampini static PetscErrorCode MatMult_Hyper(Mat M, Vec x, Vec y)
31753022affSStefano Zampini {
31853022affSStefano Zampini   PetscFunctionBegin;
3199566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_Hyper(M,x,y,PETSC_FALSE));
32053022affSStefano Zampini   PetscFunctionReturn(0);
32153022affSStefano Zampini }
32253022affSStefano Zampini 
32353022affSStefano Zampini static PetscErrorCode MatMultTranspose_Hyper(Mat M, Vec x, Vec y)
32453022affSStefano Zampini {
32553022affSStefano Zampini   PetscFunctionBegin;
3269566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_Hyper(M,x,y,PETSC_TRUE));
32753022affSStefano Zampini   PetscFunctionReturn(0);
32853022affSStefano Zampini }
32953022affSStefano Zampini 
33053022affSStefano Zampini /* Hyper power kernel, MatMat version */
33153022affSStefano Zampini static PetscErrorCode MatMatMultKernel_Hyper(Mat M, Mat X, Mat Y, PetscBool t)
33253022affSStefano Zampini {
33353022affSStefano Zampini   PC         pc;
33453022affSStefano Zampini   Mat        A;
33553022affSStefano Zampini   PC_H2OPUS *pch2opus;
33653022affSStefano Zampini   PetscInt   i;
33753022affSStefano Zampini 
33853022affSStefano Zampini   PetscFunctionBegin;
3399566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M,(void*)&pc));
34053022affSStefano Zampini   pch2opus = (PC_H2OPUS*)pc->data;
34153022affSStefano Zampini   A = pch2opus->A;
34253022affSStefano Zampini   if (pch2opus->wnsmat[0] && pch2opus->wnsmat[0]->cmap->N != X->cmap->N) {
3439566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[0]));
3449566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[1]));
34553022affSStefano Zampini   }
34653022affSStefano Zampini   if (!pch2opus->wnsmat[0]) {
3479566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(X,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[0]));
3489566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(Y,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[1]));
34953022affSStefano Zampini   }
35053022affSStefano Zampini   if (pch2opus->wnsmat[2] && pch2opus->wnsmat[2]->cmap->N != X->cmap->N) {
3519566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[2]));
3529566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[3]));
35353022affSStefano Zampini   }
35453022affSStefano Zampini   if (!pch2opus->wnsmat[2]) {
3559566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(X,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[2]));
3569566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(Y,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[3]));
35753022affSStefano Zampini   }
3589566063dSJacob Faibussowitsch   PetscCall(MatCopy(X,pch2opus->wnsmat[0],SAME_NONZERO_PATTERN));
3599566063dSJacob Faibussowitsch   PetscCall(MatCopy(X,pch2opus->wnsmat[3],SAME_NONZERO_PATTERN));
36053022affSStefano Zampini   if (t) {
36153022affSStefano Zampini     for (i=0;i<pch2opus->hyperorder-1;i++) {
3629566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(A,pch2opus->wnsmat[0],MAT_REUSE_MATRIX,PETSC_DEFAULT,&pch2opus->wnsmat[1]));
3639566063dSJacob Faibussowitsch       PetscCall(PCApplyTransposeMat_H2OPUS(pc,pch2opus->wnsmat[1],pch2opus->wnsmat[2]));
3649566063dSJacob Faibussowitsch       PetscCall(MatAXPY(pch2opus->wnsmat[0],-1.,pch2opus->wnsmat[2],SAME_NONZERO_PATTERN));
3659566063dSJacob Faibussowitsch       PetscCall(MatAXPY(pch2opus->wnsmat[3],1.,pch2opus->wnsmat[0],SAME_NONZERO_PATTERN));
36653022affSStefano Zampini     }
3679566063dSJacob Faibussowitsch     PetscCall(PCApplyTransposeMat_H2OPUS(pc,pch2opus->wnsmat[3],Y));
36853022affSStefano Zampini   } else {
36953022affSStefano Zampini     for (i=0;i<pch2opus->hyperorder-1;i++) {
3709566063dSJacob Faibussowitsch       PetscCall(PCApplyMat_H2OPUS(pc,pch2opus->wnsmat[0],pch2opus->wnsmat[1]));
3719566063dSJacob Faibussowitsch       PetscCall(MatMatMult(A,pch2opus->wnsmat[1],MAT_REUSE_MATRIX,PETSC_DEFAULT,&pch2opus->wnsmat[2]));
3729566063dSJacob Faibussowitsch       PetscCall(MatAXPY(pch2opus->wnsmat[0],-1.,pch2opus->wnsmat[2],SAME_NONZERO_PATTERN));
3739566063dSJacob Faibussowitsch       PetscCall(MatAXPY(pch2opus->wnsmat[3],1.,pch2opus->wnsmat[0],SAME_NONZERO_PATTERN));
37453022affSStefano Zampini     }
3759566063dSJacob Faibussowitsch     PetscCall(PCApplyMat_H2OPUS(pc,pch2opus->wnsmat[3],Y));
37653022affSStefano Zampini   }
37753022affSStefano Zampini   PetscFunctionReturn(0);
37853022affSStefano Zampini }
37953022affSStefano Zampini 
38053022affSStefano Zampini static PetscErrorCode MatMatMultNumeric_Hyper(Mat M, Mat X, Mat Y,void *ctx)
38153022affSStefano Zampini {
38253022affSStefano Zampini   PetscFunctionBegin;
3839566063dSJacob Faibussowitsch   PetscCall(MatMatMultKernel_Hyper(M,X,Y,PETSC_FALSE));
38453022affSStefano Zampini   PetscFunctionReturn(0);
38553022affSStefano Zampini }
38653022affSStefano Zampini 
38753022affSStefano Zampini /* Basic Newton-Schultz sampler: (2 * I - M * A)*M */
38853022affSStefano Zampini static PetscErrorCode MatMultKernel_NS(Mat M, Vec x, Vec y, PetscBool t)
38953022affSStefano Zampini {
39053022affSStefano Zampini   PC         pc;
39153022affSStefano Zampini   Mat        A;
39253022affSStefano Zampini   PC_H2OPUS *pch2opus;
39353022affSStefano Zampini 
39453022affSStefano Zampini   PetscFunctionBegin;
3959566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M,(void*)&pc));
39653022affSStefano Zampini   pch2opus = (PC_H2OPUS*)pc->data;
39753022affSStefano Zampini   A = pch2opus->A;
3989566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pch2opus->M,pch2opus->wns[0] ? NULL : &pch2opus->wns[0],pch2opus->wns[1] ? NULL : &pch2opus->wns[1]));
3999566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[0],pch2opus->boundtocpu));
4009566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[1],pch2opus->boundtocpu));
40153022affSStefano Zampini   if (t) {
4029566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose_H2OPUS(pc,x,y));
4039566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,y,pch2opus->wns[1]));
4049566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose_H2OPUS(pc,pch2opus->wns[1],pch2opus->wns[0]));
4059566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(y,-1.,2.,pch2opus->wns[0]));
40653022affSStefano Zampini   } else {
4079566063dSJacob Faibussowitsch     PetscCall(PCApply_H2OPUS(pc,x,y));
4089566063dSJacob Faibussowitsch     PetscCall(MatMult(A,y,pch2opus->wns[0]));
4099566063dSJacob Faibussowitsch     PetscCall(PCApply_H2OPUS(pc,pch2opus->wns[0],pch2opus->wns[1]));
4109566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(y,-1.,2.,pch2opus->wns[1]));
41153022affSStefano Zampini   }
41253022affSStefano Zampini   PetscFunctionReturn(0);
41353022affSStefano Zampini }
41453022affSStefano Zampini 
41553022affSStefano Zampini static PetscErrorCode MatMult_NS(Mat M, Vec x, Vec y)
41653022affSStefano Zampini {
41753022affSStefano Zampini   PetscFunctionBegin;
4189566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_NS(M,x,y,PETSC_FALSE));
41953022affSStefano Zampini   PetscFunctionReturn(0);
42053022affSStefano Zampini }
42153022affSStefano Zampini 
42253022affSStefano Zampini static PetscErrorCode MatMultTranspose_NS(Mat M, Vec x, Vec y)
42353022affSStefano Zampini {
42453022affSStefano Zampini   PetscFunctionBegin;
4259566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_NS(M,x,y,PETSC_TRUE));
42653022affSStefano Zampini   PetscFunctionReturn(0);
42753022affSStefano Zampini }
42853022affSStefano Zampini 
42953022affSStefano Zampini /* Basic Newton-Schultz sampler: (2 * I - M * A)*M, MatMat version */
43053022affSStefano Zampini static PetscErrorCode MatMatMultKernel_NS(Mat M, Mat X, Mat Y, PetscBool t)
43153022affSStefano Zampini {
43253022affSStefano Zampini   PC         pc;
43353022affSStefano Zampini   Mat        A;
43453022affSStefano Zampini   PC_H2OPUS *pch2opus;
43553022affSStefano Zampini 
43653022affSStefano Zampini   PetscFunctionBegin;
4379566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M,(void*)&pc));
43853022affSStefano Zampini   pch2opus = (PC_H2OPUS*)pc->data;
43953022affSStefano Zampini   A = pch2opus->A;
44053022affSStefano Zampini   if (pch2opus->wnsmat[0] && pch2opus->wnsmat[0]->cmap->N != X->cmap->N) {
4419566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[0]));
4429566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[1]));
44353022affSStefano Zampini   }
44453022affSStefano Zampini   if (!pch2opus->wnsmat[0]) {
4459566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(X,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[0]));
4469566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(Y,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[1]));
44753022affSStefano Zampini   }
44853022affSStefano Zampini   if (t) {
4499566063dSJacob Faibussowitsch     PetscCall(PCApplyTransposeMat_H2OPUS(pc,X,Y));
4509566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A,Y,MAT_REUSE_MATRIX,PETSC_DEFAULT,&pch2opus->wnsmat[1]));
4519566063dSJacob Faibussowitsch     PetscCall(PCApplyTransposeMat_H2OPUS(pc,pch2opus->wnsmat[1],pch2opus->wnsmat[0]));
4529566063dSJacob Faibussowitsch     PetscCall(MatScale(Y,2.));
4539566063dSJacob Faibussowitsch     PetscCall(MatAXPY(Y,-1.,pch2opus->wnsmat[0],SAME_NONZERO_PATTERN));
45453022affSStefano Zampini   } else {
4559566063dSJacob Faibussowitsch     PetscCall(PCApplyMat_H2OPUS(pc,X,Y));
4569566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A,Y,MAT_REUSE_MATRIX,PETSC_DEFAULT,&pch2opus->wnsmat[0]));
4579566063dSJacob Faibussowitsch     PetscCall(PCApplyMat_H2OPUS(pc,pch2opus->wnsmat[0],pch2opus->wnsmat[1]));
4589566063dSJacob Faibussowitsch     PetscCall(MatScale(Y,2.));
4599566063dSJacob Faibussowitsch     PetscCall(MatAXPY(Y,-1.,pch2opus->wnsmat[1],SAME_NONZERO_PATTERN));
46053022affSStefano Zampini   }
46153022affSStefano Zampini   PetscFunctionReturn(0);
46253022affSStefano Zampini }
46353022affSStefano Zampini 
46453022affSStefano Zampini static PetscErrorCode MatMatMultNumeric_NS(Mat M, Mat X, Mat Y, void *ctx)
46553022affSStefano Zampini {
46653022affSStefano Zampini   PetscFunctionBegin;
4679566063dSJacob Faibussowitsch   PetscCall(MatMatMultKernel_NS(M,X,Y,PETSC_FALSE));
46853022affSStefano Zampini   PetscFunctionReturn(0);
46953022affSStefano Zampini }
47053022affSStefano Zampini 
47153022affSStefano Zampini static PetscErrorCode PCH2OpusSetUpSampler_Private(PC pc)
47253022affSStefano Zampini {
47353022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
47453022affSStefano Zampini   Mat        A        = pc->useAmat ? pc->mat : pc->pmat;
47553022affSStefano Zampini 
47653022affSStefano Zampini   PetscFunctionBegin;
47753022affSStefano Zampini   if (!pch2opus->S) {
47853022affSStefano Zampini     PetscInt M,N,m,n;
47953022affSStefano Zampini 
4809566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A,&M,&N));
4819566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A,&m,&n));
4829566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)A),m,n,M,N,pc,&pch2opus->S));
4839566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(pch2opus->S,A,A));
48453022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU)
4859566063dSJacob Faibussowitsch     PetscCall(MatShellSetVecType(pch2opus->S,VECCUDA));
48653022affSStefano Zampini #endif
48753022affSStefano Zampini   }
48853022affSStefano Zampini   if (pch2opus->hyperorder >= 2) {
4899566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->S,MATOP_MULT,(void (*)(void))MatMult_Hyper));
4909566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_Hyper));
4919566063dSJacob Faibussowitsch     PetscCall(MatShellSetMatProductOperation(pch2opus->S,MATPRODUCT_AB,NULL,MatMatMultNumeric_Hyper,NULL,MATDENSE,MATDENSE));
4929566063dSJacob Faibussowitsch     PetscCall(MatShellSetMatProductOperation(pch2opus->S,MATPRODUCT_AB,NULL,MatMatMultNumeric_Hyper,NULL,MATDENSECUDA,MATDENSECUDA));
49353022affSStefano Zampini   } else {
4949566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->S,MATOP_MULT,(void (*)(void))MatMult_NS));
4959566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_NS));
4969566063dSJacob Faibussowitsch     PetscCall(MatShellSetMatProductOperation(pch2opus->S,MATPRODUCT_AB,NULL,MatMatMultNumeric_NS,NULL,MATDENSE,MATDENSE));
4979566063dSJacob Faibussowitsch     PetscCall(MatShellSetMatProductOperation(pch2opus->S,MATPRODUCT_AB,NULL,MatMatMultNumeric_NS,NULL,MATDENSECUDA,MATDENSECUDA));
49853022affSStefano Zampini   }
4999566063dSJacob Faibussowitsch   PetscCall(MatPropagateSymmetryOptions(A,pch2opus->S));
5009566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(pch2opus->S,pch2opus->boundtocpu));
50153022affSStefano Zampini   /* XXX */
5029566063dSJacob Faibussowitsch   PetscCall(MatSetOption(pch2opus->S,MAT_SYMMETRIC,PETSC_TRUE));
50353022affSStefano Zampini   PetscFunctionReturn(0);
50453022affSStefano Zampini }
50553022affSStefano Zampini 
50653022affSStefano Zampini static PetscErrorCode PCSetUp_H2OPUS(PC pc)
50753022affSStefano Zampini {
50853022affSStefano Zampini   PC_H2OPUS *pch2opus  = (PC_H2OPUS*)pc->data;
50953022affSStefano Zampini   Mat        A         = pc->useAmat ? pc->mat : pc->pmat;
51053022affSStefano Zampini   NormType   norm      = pch2opus->normtype;
51153022affSStefano Zampini   PetscReal  initerr   = 0.0,err;
51253022affSStefano Zampini   PetscBool  ish2opus;
51353022affSStefano Zampini 
51453022affSStefano Zampini   PetscFunctionBegin;
51553022affSStefano Zampini   if (!pch2opus->T) {
51653022affSStefano Zampini     PetscInt    M,N,m,n;
51753022affSStefano Zampini     const char *prefix;
51853022affSStefano Zampini 
5199566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc,&prefix));
5209566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A,&M,&N));
5219566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A,&m,&n));
5229566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc->pmat),m,n,M,N,pc,&pch2opus->T));
5239566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(pch2opus->T,A,A));
5249566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->T,MATOP_MULT,(void (*)(void))MatMult_MAmI));
5259566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->T,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_MAmI));
5269566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->T,MATOP_NORM,(void (*)(void))MatNorm_H2OPUS));
52753022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU)
5289566063dSJacob Faibussowitsch     PetscCall(MatShellSetVecType(pch2opus->T,VECCUDA));
52953022affSStefano Zampini #endif
5309566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)pch2opus->T));
5319566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(pch2opus->T,prefix));
5329566063dSJacob Faibussowitsch     PetscCall(MatAppendOptionsPrefix(pch2opus->T,"pc_h2opus_est_"));
53353022affSStefano Zampini   }
5349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->A));
5359566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATH2OPUS,&ish2opus));
53653022affSStefano Zampini   if (ish2opus) {
5379566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
53853022affSStefano Zampini     pch2opus->A = A;
53953022affSStefano Zampini   } else {
54053022affSStefano Zampini     const char *prefix;
5419566063dSJacob Faibussowitsch     PetscCall(MatCreateH2OpusFromMat(A,pch2opus->sdim,pch2opus->coords,PETSC_FALSE,pch2opus->eta,pch2opus->leafsize,pch2opus->max_rank,pch2opus->bs,pch2opus->mrtol,&pch2opus->A));
54253022affSStefano Zampini     /* XXX */
5439566063dSJacob Faibussowitsch     PetscCall(MatSetOption(pch2opus->A,MAT_SYMMETRIC,PETSC_TRUE));
5449566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc,&prefix));
5459566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(pch2opus->A,prefix));
5469566063dSJacob Faibussowitsch     PetscCall(MatAppendOptionsPrefix(pch2opus->A,"pc_h2opus_init_"));
5479566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(pch2opus->A));
5489566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(pch2opus->A,MAT_FINAL_ASSEMBLY));
5499566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(pch2opus->A,MAT_FINAL_ASSEMBLY));
55053022affSStefano Zampini     /* XXX */
5519566063dSJacob Faibussowitsch     PetscCall(MatSetOption(pch2opus->A,MAT_SYMMETRIC,PETSC_TRUE));
552*8db51263SStefano Zampini 
553*8db51263SStefano Zampini     /* always perform construction on the GPU unless forcecpu is true */
554*8db51263SStefano Zampini     PetscCall(MatBindToCPU(pch2opus->A,pch2opus->forcecpu));
55553022affSStefano Zampini   }
55653022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU)
557*8db51263SStefano Zampini   pch2opus->boundtocpu = pch2opus->forcecpu ? PETSC_TRUE : pch2opus->A->boundtocpu;
55853022affSStefano Zampini #endif
5599566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(pch2opus->T,pch2opus->boundtocpu));
56053022affSStefano Zampini   if (pch2opus->M) { /* see if we can reuse M as initial guess */
56153022affSStefano Zampini     PetscReal reuse;
56253022affSStefano Zampini 
5639566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(pch2opus->M,pch2opus->boundtocpu));
5649566063dSJacob Faibussowitsch     PetscCall(MatNorm(pch2opus->T,norm,&reuse));
5659566063dSJacob Faibussowitsch     if (reuse >= 1.0) PetscCall(MatDestroy(&pch2opus->M));
56653022affSStefano Zampini   }
56753022affSStefano Zampini   if (!pch2opus->M) {
56853022affSStefano Zampini     const char *prefix;
5699566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(pch2opus->A,MAT_COPY_VALUES,&pch2opus->M));
5709566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc,&prefix));
5719566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(pch2opus->M,prefix));
5729566063dSJacob Faibussowitsch     PetscCall(MatAppendOptionsPrefix(pch2opus->M,"pc_h2opus_inv_"));
5739566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(pch2opus->M));
5749566063dSJacob Faibussowitsch     PetscCall(PCH2OpusSetUpInit(pc));
5759566063dSJacob Faibussowitsch     PetscCall(MatScale(pch2opus->M,pch2opus->s0));
57653022affSStefano Zampini   }
57753022affSStefano Zampini   /* A and M have the same h2 matrix structure, save on reordering routines */
5789566063dSJacob Faibussowitsch   PetscCall(MatH2OpusSetNativeMult(pch2opus->A,PETSC_TRUE));
5799566063dSJacob Faibussowitsch   PetscCall(MatH2OpusSetNativeMult(pch2opus->M,PETSC_TRUE));
580*8db51263SStefano Zampini   if (norm == NORM_1 || norm == NORM_2 || norm == NORM_INFINITY) PetscCall(MatNorm(pch2opus->T,norm,&initerr));
58153022affSStefano Zampini   if (PetscIsInfOrNanReal(initerr)) pc->failedreason = PC_SETUP_ERROR;
58253022affSStefano Zampini   err   = initerr;
58353022affSStefano Zampini   if (pch2opus->monitor) {
58463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pc),"%" PetscInt_FMT ": ||M*A - I|| NORM%s abs %g rel %g\n",0,NormTypes[norm],(double)err,(double)(err/initerr)));
58553022affSStefano Zampini   }
58653022affSStefano Zampini   if (initerr > pch2opus->atol && !pc->failedreason) {
58753022affSStefano Zampini     PetscInt i;
58853022affSStefano Zampini 
5899566063dSJacob Faibussowitsch     PetscCall(PCH2OpusSetUpSampler_Private(pc));
59053022affSStefano Zampini     for (i = 0; i < pch2opus->maxits; i++) {
59153022affSStefano Zampini       Mat         M;
59253022affSStefano Zampini       const char* prefix;
59353022affSStefano Zampini 
5949566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(pch2opus->M,MAT_SHARE_NONZERO_PATTERN,&M));
5959566063dSJacob Faibussowitsch       PetscCall(MatGetOptionsPrefix(pch2opus->M,&prefix));
5969566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(M,prefix));
5979566063dSJacob Faibussowitsch       PetscCall(MatH2OpusSetSamplingMat(M,pch2opus->S,PETSC_DECIDE,PETSC_DECIDE));
5989566063dSJacob Faibussowitsch       PetscCall(MatSetFromOptions(M));
5999566063dSJacob Faibussowitsch       PetscCall(MatH2OpusSetNativeMult(M,PETSC_TRUE));
6009566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
6019566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
6029566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pch2opus->M));
60353022affSStefano Zampini       pch2opus->M = M;
604*8db51263SStefano Zampini       if (norm == NORM_1 || norm == NORM_2 || norm == NORM_INFINITY) PetscCall(MatNorm(pch2opus->T,norm,&err));
60553022affSStefano Zampini       if (pch2opus->monitor) {
60663a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pc),"%" PetscInt_FMT ": ||M*A - I|| NORM%s abs %g rel %g\n",i+1,NormTypes[norm],(double)err,(double)(err/initerr)));
60753022affSStefano Zampini       }
60853022affSStefano Zampini       if (PetscIsInfOrNanReal(err)) pc->failedreason = PC_SETUP_ERROR;
60953022affSStefano Zampini       if (err < pch2opus->atol || err < pch2opus->rtol*initerr || pc->failedreason) break;
61053022affSStefano Zampini     }
61153022affSStefano Zampini   }
61253022affSStefano Zampini   /* cleanup setup workspace */
6139566063dSJacob Faibussowitsch   PetscCall(MatH2OpusSetNativeMult(pch2opus->A,PETSC_FALSE));
6149566063dSJacob Faibussowitsch   PetscCall(MatH2OpusSetNativeMult(pch2opus->M,PETSC_FALSE));
6159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[0]));
6169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[1]));
6179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[2]));
6189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[3]));
6199566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[0]));
6209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[1]));
6219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[2]));
6229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[3]));
62353022affSStefano Zampini   PetscFunctionReturn(0);
62453022affSStefano Zampini }
62553022affSStefano Zampini 
62653022affSStefano Zampini static PetscErrorCode PCView_H2OPUS(PC pc, PetscViewer viewer)
62753022affSStefano Zampini {
62853022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
62953022affSStefano Zampini   PetscBool  isascii;
63053022affSStefano Zampini 
63153022affSStefano Zampini   PetscFunctionBegin;
6329566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
63353022affSStefano Zampini   if (isascii) {
63453022affSStefano Zampini     if (pch2opus->A && pch2opus->A != pc->mat && pch2opus->A != pc->pmat) {
6359566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Initial approximation matrix\n"));
6369566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
6379566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
6389566063dSJacob Faibussowitsch       PetscCall(MatView(pch2opus->A,viewer));
6399566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
6409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
64153022affSStefano Zampini     }
64253022affSStefano Zampini     if (pch2opus->M) {
6439566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Inner matrix constructed\n"));
6449566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
6459566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
6469566063dSJacob Faibussowitsch       PetscCall(MatView(pch2opus->M,viewer));
6479566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
6489566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
64953022affSStefano Zampini     }
6509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Initial scaling: %g\n",pch2opus->s0));
65153022affSStefano Zampini   }
65253022affSStefano Zampini   PetscFunctionReturn(0);
65353022affSStefano Zampini }
65453022affSStefano Zampini 
65553022affSStefano Zampini PETSC_EXTERN PetscErrorCode PCCreate_H2OPUS(PC pc)
65653022affSStefano Zampini {
65753022affSStefano Zampini   PC_H2OPUS *pch2opus;
65853022affSStefano Zampini 
65953022affSStefano Zampini   PetscFunctionBegin;
6609566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&pch2opus));
66153022affSStefano Zampini   pc->data = (void*)pch2opus;
66253022affSStefano Zampini 
66353022affSStefano Zampini   pch2opus->atol       = 1.e-2;
66453022affSStefano Zampini   pch2opus->rtol       = 1.e-6;
66553022affSStefano Zampini   pch2opus->maxits     = 50;
66653022affSStefano Zampini   pch2opus->hyperorder = 1; /* defaults to basic NewtonSchultz */
66753022affSStefano Zampini   pch2opus->normtype   = NORM_2;
66853022affSStefano Zampini 
66953022affSStefano Zampini   /* these are needed when we are sampling the pmat */
67053022affSStefano Zampini   pch2opus->eta        = PETSC_DECIDE;
67153022affSStefano Zampini   pch2opus->leafsize   = PETSC_DECIDE;
67253022affSStefano Zampini   pch2opus->max_rank   = PETSC_DECIDE;
67353022affSStefano Zampini   pch2opus->bs         = PETSC_DECIDE;
67453022affSStefano Zampini   pch2opus->mrtol      = PETSC_DECIDE;
67553022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU)
67653022affSStefano Zampini   pch2opus->boundtocpu = PETSC_FALSE;
67753022affSStefano Zampini #else
67853022affSStefano Zampini   pch2opus->boundtocpu = PETSC_TRUE;
67953022affSStefano Zampini #endif
68053022affSStefano Zampini   pc->ops->destroy        = PCDestroy_H2OPUS;
68153022affSStefano Zampini   pc->ops->setup          = PCSetUp_H2OPUS;
68253022affSStefano Zampini   pc->ops->apply          = PCApply_H2OPUS;
68353022affSStefano Zampini   pc->ops->matapply       = PCApplyMat_H2OPUS;
68453022affSStefano Zampini   pc->ops->applytranspose = PCApplyTranspose_H2OPUS;
68553022affSStefano Zampini   pc->ops->reset          = PCReset_H2OPUS;
68653022affSStefano Zampini   pc->ops->setfromoptions = PCSetFromOptions_H2OPUS;
68753022affSStefano Zampini   pc->ops->view           = PCView_H2OPUS;
68853022affSStefano Zampini 
6899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_H2OPUS));
69053022affSStefano Zampini   PetscFunctionReturn(0);
69153022affSStefano Zampini }
692