xref: /petsc/src/ksp/pc/impls/h2opus/pch2opus.c (revision 1c2dc1cbabe5212f80cf309c4ca5a26f0cadc660)
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   PetscBool testMA;
2553022affSStefano Zampini 
2653022affSStefano Zampini   /* Support for PCSetCoordinates */
2753022affSStefano Zampini   PetscInt  sdim;
2853022affSStefano Zampini   PetscInt  nlocc;
2953022affSStefano Zampini   PetscReal *coords;
3053022affSStefano Zampini 
3153022affSStefano Zampini   /* Newton-Schultz customization */
3253022affSStefano Zampini   PetscInt  maxits;
3353022affSStefano Zampini   PetscReal rtol,atol;
3453022affSStefano Zampini   PetscBool monitor;
3553022affSStefano Zampini   PetscBool useapproximatenorms;
3653022affSStefano Zampini   NormType  normtype;
3753022affSStefano Zampini 
3853022affSStefano Zampini   /* Used when pmat != MATH2OPUS */
3953022affSStefano Zampini   PetscReal eta;
4053022affSStefano Zampini   PetscInt  leafsize;
4153022affSStefano Zampini   PetscInt  max_rank;
4253022affSStefano Zampini   PetscInt  bs;
4353022affSStefano Zampini   PetscReal mrtol;
4453022affSStefano Zampini 
4553022affSStefano Zampini   PetscBool boundtocpu;
4653022affSStefano Zampini } PC_H2OPUS;
4753022affSStefano Zampini 
4853022affSStefano Zampini PETSC_EXTERN PetscErrorCode MatNorm_H2OPUS(Mat,NormType,PetscReal*);
4953022affSStefano Zampini 
5053022affSStefano Zampini static PetscErrorCode PCReset_H2OPUS(PC pc)
5153022affSStefano Zampini {
5253022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
5353022affSStefano Zampini 
5453022affSStefano Zampini   PetscFunctionBegin;
5553022affSStefano Zampini   pch2opus->sdim  = 0;
5653022affSStefano Zampini   pch2opus->nlocc = 0;
579566063dSJacob Faibussowitsch   PetscCall(PetscFree(pch2opus->coords));
589566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->A));
599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->M));
609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->T));
619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->w));
629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->S));
639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[0]));
649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[1]));
659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[2]));
669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[3]));
679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[0]));
689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[1]));
699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[2]));
709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[3]));
7153022affSStefano Zampini   PetscFunctionReturn(0);
7253022affSStefano Zampini }
7353022affSStefano Zampini 
7453022affSStefano Zampini static PetscErrorCode PCSetCoordinates_H2OPUS(PC pc, PetscInt sdim, PetscInt nlocc, PetscReal *coords)
7553022affSStefano Zampini {
7653022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
7753022affSStefano Zampini   PetscBool  reset    = PETSC_TRUE;
7853022affSStefano Zampini 
7953022affSStefano Zampini   PetscFunctionBegin;
8053022affSStefano Zampini   if (pch2opus->sdim && sdim == pch2opus->sdim && nlocc == pch2opus->nlocc) {
819566063dSJacob Faibussowitsch     PetscCall(PetscArraycmp(pch2opus->coords,coords,sdim*nlocc,&reset));
8253022affSStefano Zampini     reset = (PetscBool)!reset;
8353022affSStefano Zampini   }
84*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(MPI_IN_PLACE,&reset,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc)));
8553022affSStefano Zampini   if (reset) {
869566063dSJacob Faibussowitsch     PetscCall(PCReset_H2OPUS(pc));
879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(sdim*nlocc,&pch2opus->coords));
889566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pch2opus->coords,coords,sdim*nlocc));
8953022affSStefano Zampini     pch2opus->sdim  = sdim;
9053022affSStefano Zampini     pch2opus->nlocc = nlocc;
9153022affSStefano Zampini   }
9253022affSStefano Zampini   PetscFunctionReturn(0);
9353022affSStefano Zampini }
9453022affSStefano Zampini 
9553022affSStefano Zampini static PetscErrorCode PCDestroy_H2OPUS(PC pc)
9653022affSStefano Zampini {
9753022affSStefano Zampini   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(PCReset_H2OPUS(pc));
999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL));
1009566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
10153022affSStefano Zampini   PetscFunctionReturn(0);
10253022affSStefano Zampini }
10353022affSStefano Zampini 
10453022affSStefano Zampini static PetscErrorCode PCSetFromOptions_H2OPUS(PetscOptionItems *PetscOptionsObject,PC pc)
10553022affSStefano Zampini {
10653022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
10753022affSStefano Zampini 
10853022affSStefano Zampini   PetscFunctionBegin;
1099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"H2OPUS options"));
1109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_h2opus_maxits","Maximum number of iterations for Newton-Schultz",NULL,pch2opus->maxits,&pch2opus->maxits,NULL));
1119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_h2opus_monitor","Monitor Newton-Schultz convergence",NULL,pch2opus->monitor,&pch2opus->monitor,NULL));
1129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_h2opus_atol","Absolute tolerance",NULL,pch2opus->atol,&pch2opus->atol,NULL));
1139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_h2opus_rtol","Relative tolerance",NULL,pch2opus->rtol,&pch2opus->rtol,NULL));
1149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_h2opus_norm_type","Norm type for convergence monitoring",NULL,NormTypes,(PetscEnum)pch2opus->normtype,(PetscEnum*)&pch2opus->normtype,NULL));
1159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_h2opus_hyperorder","Hyper power order of sampling",NULL,pch2opus->hyperorder,&pch2opus->hyperorder,NULL));
1169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_h2opus_leafsize","Leaf size when constructed from kernel",NULL,pch2opus->leafsize,&pch2opus->leafsize,NULL));
1179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_h2opus_eta","Admissibility condition tolerance",NULL,pch2opus->eta,&pch2opus->eta,NULL));
1189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_h2opus_maxrank","Maximum rank when constructed from matvecs",NULL,pch2opus->max_rank,&pch2opus->max_rank,NULL));
1199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_h2opus_samples","Number of samples to be taken concurrently when constructing from matvecs",NULL,pch2opus->bs,&pch2opus->bs,NULL));
1209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_h2opus_mrtol","Relative tolerance for construction from sampling",NULL,pch2opus->mrtol,&pch2opus->mrtol,NULL));
1219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
12253022affSStefano Zampini   PetscFunctionReturn(0);
12353022affSStefano Zampini }
12453022affSStefano Zampini 
12553022affSStefano Zampini typedef struct {
12653022affSStefano Zampini   Mat A;
12753022affSStefano Zampini   Mat M;
12853022affSStefano Zampini   Vec w;
12953022affSStefano Zampini } AAtCtx;
13053022affSStefano Zampini 
13153022affSStefano Zampini static PetscErrorCode MatMult_AAt(Mat A, Vec x, Vec y)
13253022affSStefano Zampini {
13353022affSStefano Zampini   AAtCtx *aat;
13453022affSStefano Zampini 
13553022affSStefano Zampini   PetscFunctionBegin;
1369566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A,(void*)&aat));
1379566063dSJacob Faibussowitsch   /* PetscCall(MatMultTranspose(aat->M,x,aat->w)); */
1389566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(aat->A,x,aat->w));
1399566063dSJacob Faibussowitsch   PetscCall(MatMult(aat->A,aat->w,y));
14053022affSStefano Zampini   PetscFunctionReturn(0);
14153022affSStefano Zampini }
14253022affSStefano Zampini 
14353022affSStefano Zampini static PetscErrorCode PCH2OpusSetUpInit(PC pc)
14453022affSStefano Zampini {
14553022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
14653022affSStefano Zampini   Mat        A        = pc->useAmat ? pc->mat : pc->pmat, AAt;
14753022affSStefano Zampini   PetscInt   M,m;
14853022affSStefano Zampini   VecType    vtype;
14953022affSStefano Zampini   PetscReal  n;
15053022affSStefano Zampini   AAtCtx     aat;
15153022affSStefano Zampini 
15253022affSStefano Zampini   PetscFunctionBegin;
15353022affSStefano Zampini   aat.A = A;
15453022affSStefano Zampini   aat.M = pch2opus->M; /* unused so far */
1559566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,NULL,&aat.w));
1569566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&M,NULL));
1579566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,NULL));
1589566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)A),m,m,M,M,&aat,&AAt));
1599566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(AAt,pch2opus->boundtocpu));
1609566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(AAt,MATOP_MULT,(void (*)(void))MatMult_AAt));
1619566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(AAt,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMult_AAt));
1629566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(AAt,MATOP_NORM,(void (*)(void))MatNorm_H2OPUS));
1639566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A,&vtype));
1649566063dSJacob Faibussowitsch   PetscCall(MatShellSetVecType(AAt,vtype));
1659566063dSJacob Faibussowitsch   PetscCall(MatNorm(AAt,NORM_1,&n));
16653022affSStefano Zampini   pch2opus->s0 = 1./n;
1679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aat.w));
1689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&AAt));
16953022affSStefano Zampini   PetscFunctionReturn(0);
17053022affSStefano Zampini }
17153022affSStefano Zampini 
17253022affSStefano Zampini static PetscErrorCode PCApplyKernel_H2OPUS(PC pc, Vec x, Vec y, PetscBool t)
17353022affSStefano Zampini {
17453022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
17553022affSStefano Zampini 
17653022affSStefano Zampini   PetscFunctionBegin;
1779566063dSJacob Faibussowitsch   if (t) PetscCall(MatMultTranspose(pch2opus->M,x,y));
1789566063dSJacob Faibussowitsch   else  PetscCall(MatMult(pch2opus->M,x,y));
17953022affSStefano Zampini   PetscFunctionReturn(0);
18053022affSStefano Zampini }
18153022affSStefano Zampini 
18253022affSStefano Zampini static PetscErrorCode PCApplyMatKernel_H2OPUS(PC pc, Mat X, Mat Y, PetscBool t)
18353022affSStefano Zampini {
18453022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
18553022affSStefano Zampini 
18653022affSStefano Zampini   PetscFunctionBegin;
1879566063dSJacob Faibussowitsch   if (t) PetscCall(MatTransposeMatMult(pch2opus->M,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y));
1889566063dSJacob Faibussowitsch   else   PetscCall(MatMatMult(pch2opus->M,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y));
18953022affSStefano Zampini   PetscFunctionReturn(0);
19053022affSStefano Zampini }
19153022affSStefano Zampini 
19253022affSStefano Zampini static PetscErrorCode PCApplyMat_H2OPUS(PC pc, Mat X, Mat Y)
19353022affSStefano Zampini {
19453022affSStefano Zampini   PetscFunctionBegin;
1959566063dSJacob Faibussowitsch   PetscCall(PCApplyMatKernel_H2OPUS(pc,X,Y,PETSC_FALSE));
19653022affSStefano Zampini   PetscFunctionReturn(0);
19753022affSStefano Zampini }
19853022affSStefano Zampini 
19953022affSStefano Zampini static PetscErrorCode PCApplyTransposeMat_H2OPUS(PC pc, Mat X, Mat Y)
20053022affSStefano Zampini {
20153022affSStefano Zampini   PetscFunctionBegin;
2029566063dSJacob Faibussowitsch   PetscCall(PCApplyMatKernel_H2OPUS(pc,X,Y,PETSC_TRUE));
20353022affSStefano Zampini   PetscFunctionReturn(0);
20453022affSStefano Zampini }
20553022affSStefano Zampini 
20653022affSStefano Zampini static PetscErrorCode PCApply_H2OPUS(PC pc, Vec x, Vec y)
20753022affSStefano Zampini {
20853022affSStefano Zampini   PetscFunctionBegin;
2099566063dSJacob Faibussowitsch   PetscCall(PCApplyKernel_H2OPUS(pc,x,y,PETSC_FALSE));
21053022affSStefano Zampini   PetscFunctionReturn(0);
21153022affSStefano Zampini }
21253022affSStefano Zampini 
21353022affSStefano Zampini static PetscErrorCode PCApplyTranspose_H2OPUS(PC pc, Vec x, Vec y)
21453022affSStefano Zampini {
21553022affSStefano Zampini   PetscFunctionBegin;
2169566063dSJacob Faibussowitsch   PetscCall(PCApplyKernel_H2OPUS(pc,x,y,PETSC_TRUE));
21753022affSStefano Zampini   PetscFunctionReturn(0);
21853022affSStefano Zampini }
21953022affSStefano Zampini 
22053022affSStefano Zampini /* used to test the norm of (M^-1 A - I) */
22153022affSStefano Zampini static PetscErrorCode MatMultKernel_MAmI(Mat M, Vec x, Vec y, PetscBool t)
22253022affSStefano Zampini {
22353022affSStefano Zampini   PC         pc;
22453022affSStefano Zampini   Mat        A;
22553022affSStefano Zampini   PC_H2OPUS *pch2opus;
22653022affSStefano Zampini   PetscBool  sideleft = PETSC_TRUE;
22753022affSStefano Zampini 
22853022affSStefano Zampini   PetscFunctionBegin;
2299566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M,(void*)&pc));
23053022affSStefano Zampini   pch2opus = (PC_H2OPUS*)pc->data;
2319566063dSJacob Faibussowitsch   if (!pch2opus->w) PetscCall(MatCreateVecs(pch2opus->M,&pch2opus->w,NULL));
23253022affSStefano Zampini   A = pch2opus->A;
2339566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->w,pch2opus->boundtocpu));
23453022affSStefano Zampini   if (t) {
23553022affSStefano Zampini     if (sideleft) {
2369566063dSJacob Faibussowitsch       PetscCall(PCApplyTranspose_H2OPUS(pc,x,pch2opus->w));
2379566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A,pch2opus->w,y));
23853022affSStefano Zampini     } else {
2399566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A,x,pch2opus->w));
2409566063dSJacob Faibussowitsch       PetscCall(PCApplyTranspose_H2OPUS(pc,pch2opus->w,y));
24153022affSStefano Zampini     }
24253022affSStefano Zampini   } else {
24353022affSStefano Zampini     if (sideleft) {
2449566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,pch2opus->w));
2459566063dSJacob Faibussowitsch       PetscCall(PCApply_H2OPUS(pc,pch2opus->w,y));
24653022affSStefano Zampini     } else {
2479566063dSJacob Faibussowitsch       PetscCall(PCApply_H2OPUS(pc,x,pch2opus->w));
2489566063dSJacob Faibussowitsch       PetscCall(MatMult(A,pch2opus->w,y));
24953022affSStefano Zampini     }
25053022affSStefano Zampini   }
2519566063dSJacob Faibussowitsch   if (!pch2opus->testMA) PetscCall(VecAXPY(y,-1.0,x));
25253022affSStefano Zampini   PetscFunctionReturn(0);
25353022affSStefano Zampini }
25453022affSStefano Zampini 
25553022affSStefano Zampini static PetscErrorCode MatMult_MAmI(Mat A, Vec x, Vec y)
25653022affSStefano Zampini {
25753022affSStefano Zampini   PetscFunctionBegin;
2589566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_MAmI(A,x,y,PETSC_FALSE));
25953022affSStefano Zampini   PetscFunctionReturn(0);
26053022affSStefano Zampini }
26153022affSStefano Zampini 
26253022affSStefano Zampini static PetscErrorCode MatMultTranspose_MAmI(Mat A, Vec x, Vec y)
26353022affSStefano Zampini {
26453022affSStefano Zampini   PetscFunctionBegin;
2659566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_MAmI(A,x,y,PETSC_TRUE));
26653022affSStefano Zampini   PetscFunctionReturn(0);
26753022affSStefano Zampini }
26853022affSStefano Zampini 
26953022affSStefano Zampini /* HyperPower kernel:
27053022affSStefano Zampini Y = R = x
27153022affSStefano Zampini for i = 1 . . . l - 1 do
27253022affSStefano Zampini   R = (I - A * Xk) * R
27353022affSStefano Zampini   Y = Y + R
27453022affSStefano Zampini Y = Xk * Y
27553022affSStefano Zampini */
27653022affSStefano Zampini static PetscErrorCode MatMultKernel_Hyper(Mat M, Vec x, Vec y, PetscBool t)
27753022affSStefano Zampini {
27853022affSStefano Zampini   PC         pc;
27953022affSStefano Zampini   Mat        A;
28053022affSStefano Zampini   PC_H2OPUS *pch2opus;
28153022affSStefano Zampini 
28253022affSStefano Zampini   PetscFunctionBegin;
2839566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M,(void*)&pc));
28453022affSStefano Zampini   pch2opus = (PC_H2OPUS*)pc->data;
28553022affSStefano Zampini   A = pch2opus->A;
2869566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pch2opus->M,pch2opus->wns[0] ? NULL : &pch2opus->wns[0],pch2opus->wns[1] ? NULL : &pch2opus->wns[1]));
2879566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pch2opus->M,pch2opus->wns[2] ? NULL : &pch2opus->wns[2],pch2opus->wns[3] ? NULL : &pch2opus->wns[3]));
2889566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[0],pch2opus->boundtocpu));
2899566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[1],pch2opus->boundtocpu));
2909566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[2],pch2opus->boundtocpu));
2919566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[3],pch2opus->boundtocpu));
2929566063dSJacob Faibussowitsch   PetscCall(VecCopy(x,pch2opus->wns[0]));
2939566063dSJacob Faibussowitsch   PetscCall(VecCopy(x,pch2opus->wns[3]));
29453022affSStefano Zampini   if (t) {
2955f80ce2aSJacob Faibussowitsch     for (PetscInt i=0;i<pch2opus->hyperorder-1;i++) {
2969566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A,pch2opus->wns[0],pch2opus->wns[1]));
2979566063dSJacob Faibussowitsch       PetscCall(PCApplyTranspose_H2OPUS(pc,pch2opus->wns[1],pch2opus->wns[2]));
2989566063dSJacob Faibussowitsch       PetscCall(VecAXPY(pch2opus->wns[0],-1.,pch2opus->wns[2]));
2999566063dSJacob Faibussowitsch       PetscCall(VecAXPY(pch2opus->wns[3], 1.,pch2opus->wns[0]));
30053022affSStefano Zampini     }
3019566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose_H2OPUS(pc,pch2opus->wns[3],y));
30253022affSStefano Zampini   } else {
3035f80ce2aSJacob Faibussowitsch     for (PetscInt i=0;i<pch2opus->hyperorder-1;i++) {
3049566063dSJacob Faibussowitsch       PetscCall(PCApply_H2OPUS(pc,pch2opus->wns[0],pch2opus->wns[1]));
3059566063dSJacob Faibussowitsch       PetscCall(MatMult(A,pch2opus->wns[1],pch2opus->wns[2]));
3069566063dSJacob Faibussowitsch       PetscCall(VecAXPY(pch2opus->wns[0],-1.,pch2opus->wns[2]));
3079566063dSJacob Faibussowitsch       PetscCall(VecAXPY(pch2opus->wns[3], 1.,pch2opus->wns[0]));
30853022affSStefano Zampini     }
3099566063dSJacob Faibussowitsch     PetscCall(PCApply_H2OPUS(pc,pch2opus->wns[3],y));
31053022affSStefano Zampini   }
31153022affSStefano Zampini   PetscFunctionReturn(0);
31253022affSStefano Zampini }
31353022affSStefano Zampini 
31453022affSStefano Zampini static PetscErrorCode MatMult_Hyper(Mat M, Vec x, Vec y)
31553022affSStefano Zampini {
31653022affSStefano Zampini   PetscFunctionBegin;
3179566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_Hyper(M,x,y,PETSC_FALSE));
31853022affSStefano Zampini   PetscFunctionReturn(0);
31953022affSStefano Zampini }
32053022affSStefano Zampini 
32153022affSStefano Zampini static PetscErrorCode MatMultTranspose_Hyper(Mat M, Vec x, Vec y)
32253022affSStefano Zampini {
32353022affSStefano Zampini   PetscFunctionBegin;
3249566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_Hyper(M,x,y,PETSC_TRUE));
32553022affSStefano Zampini   PetscFunctionReturn(0);
32653022affSStefano Zampini }
32753022affSStefano Zampini 
32853022affSStefano Zampini /* Hyper power kernel, MatMat version */
32953022affSStefano Zampini static PetscErrorCode MatMatMultKernel_Hyper(Mat M, Mat X, Mat Y, PetscBool t)
33053022affSStefano Zampini {
33153022affSStefano Zampini   PC         pc;
33253022affSStefano Zampini   Mat        A;
33353022affSStefano Zampini   PC_H2OPUS *pch2opus;
33453022affSStefano Zampini   PetscInt   i;
33553022affSStefano Zampini 
33653022affSStefano Zampini   PetscFunctionBegin;
3379566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M,(void*)&pc));
33853022affSStefano Zampini   pch2opus = (PC_H2OPUS*)pc->data;
33953022affSStefano Zampini   A = pch2opus->A;
34053022affSStefano Zampini   if (pch2opus->wnsmat[0] && pch2opus->wnsmat[0]->cmap->N != X->cmap->N) {
3419566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[0]));
3429566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[1]));
34353022affSStefano Zampini   }
34453022affSStefano Zampini   if (!pch2opus->wnsmat[0]) {
3459566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(X,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[0]));
3469566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(Y,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[1]));
34753022affSStefano Zampini   }
34853022affSStefano Zampini   if (pch2opus->wnsmat[2] && pch2opus->wnsmat[2]->cmap->N != X->cmap->N) {
3499566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[2]));
3509566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[3]));
35153022affSStefano Zampini   }
35253022affSStefano Zampini   if (!pch2opus->wnsmat[2]) {
3539566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(X,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[2]));
3549566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(Y,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[3]));
35553022affSStefano Zampini   }
3569566063dSJacob Faibussowitsch   PetscCall(MatCopy(X,pch2opus->wnsmat[0],SAME_NONZERO_PATTERN));
3579566063dSJacob Faibussowitsch   PetscCall(MatCopy(X,pch2opus->wnsmat[3],SAME_NONZERO_PATTERN));
35853022affSStefano Zampini   if (t) {
35953022affSStefano Zampini     for (i=0;i<pch2opus->hyperorder-1;i++) {
3609566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(A,pch2opus->wnsmat[0],MAT_REUSE_MATRIX,PETSC_DEFAULT,&pch2opus->wnsmat[1]));
3619566063dSJacob Faibussowitsch       PetscCall(PCApplyTransposeMat_H2OPUS(pc,pch2opus->wnsmat[1],pch2opus->wnsmat[2]));
3629566063dSJacob Faibussowitsch       PetscCall(MatAXPY(pch2opus->wnsmat[0],-1.,pch2opus->wnsmat[2],SAME_NONZERO_PATTERN));
3639566063dSJacob Faibussowitsch       PetscCall(MatAXPY(pch2opus->wnsmat[3],1.,pch2opus->wnsmat[0],SAME_NONZERO_PATTERN));
36453022affSStefano Zampini     }
3659566063dSJacob Faibussowitsch     PetscCall(PCApplyTransposeMat_H2OPUS(pc,pch2opus->wnsmat[3],Y));
36653022affSStefano Zampini   } else {
36753022affSStefano Zampini     for (i=0;i<pch2opus->hyperorder-1;i++) {
3689566063dSJacob Faibussowitsch       PetscCall(PCApplyMat_H2OPUS(pc,pch2opus->wnsmat[0],pch2opus->wnsmat[1]));
3699566063dSJacob Faibussowitsch       PetscCall(MatMatMult(A,pch2opus->wnsmat[1],MAT_REUSE_MATRIX,PETSC_DEFAULT,&pch2opus->wnsmat[2]));
3709566063dSJacob Faibussowitsch       PetscCall(MatAXPY(pch2opus->wnsmat[0],-1.,pch2opus->wnsmat[2],SAME_NONZERO_PATTERN));
3719566063dSJacob Faibussowitsch       PetscCall(MatAXPY(pch2opus->wnsmat[3],1.,pch2opus->wnsmat[0],SAME_NONZERO_PATTERN));
37253022affSStefano Zampini     }
3739566063dSJacob Faibussowitsch     PetscCall(PCApplyMat_H2OPUS(pc,pch2opus->wnsmat[3],Y));
37453022affSStefano Zampini   }
37553022affSStefano Zampini   PetscFunctionReturn(0);
37653022affSStefano Zampini }
37753022affSStefano Zampini 
37853022affSStefano Zampini static PetscErrorCode MatMatMultNumeric_Hyper(Mat M, Mat X, Mat Y,void *ctx)
37953022affSStefano Zampini {
38053022affSStefano Zampini   PetscFunctionBegin;
3819566063dSJacob Faibussowitsch   PetscCall(MatMatMultKernel_Hyper(M,X,Y,PETSC_FALSE));
38253022affSStefano Zampini   PetscFunctionReturn(0);
38353022affSStefano Zampini }
38453022affSStefano Zampini 
38553022affSStefano Zampini /* Basic Newton-Schultz sampler: (2 * I - M * A)*M */
38653022affSStefano Zampini static PetscErrorCode MatMultKernel_NS(Mat M, Vec x, Vec y, PetscBool t)
38753022affSStefano Zampini {
38853022affSStefano Zampini   PC         pc;
38953022affSStefano Zampini   Mat        A;
39053022affSStefano Zampini   PC_H2OPUS *pch2opus;
39153022affSStefano Zampini 
39253022affSStefano Zampini   PetscFunctionBegin;
3939566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M,(void*)&pc));
39453022affSStefano Zampini   pch2opus = (PC_H2OPUS*)pc->data;
39553022affSStefano Zampini   A = pch2opus->A;
3969566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pch2opus->M,pch2opus->wns[0] ? NULL : &pch2opus->wns[0],pch2opus->wns[1] ? NULL : &pch2opus->wns[1]));
3979566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[0],pch2opus->boundtocpu));
3989566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(pch2opus->wns[1],pch2opus->boundtocpu));
39953022affSStefano Zampini   if (t) {
4009566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose_H2OPUS(pc,x,y));
4019566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,y,pch2opus->wns[1]));
4029566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose_H2OPUS(pc,pch2opus->wns[1],pch2opus->wns[0]));
4039566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(y,-1.,2.,pch2opus->wns[0]));
40453022affSStefano Zampini   } else {
4059566063dSJacob Faibussowitsch     PetscCall(PCApply_H2OPUS(pc,x,y));
4069566063dSJacob Faibussowitsch     PetscCall(MatMult(A,y,pch2opus->wns[0]));
4079566063dSJacob Faibussowitsch     PetscCall(PCApply_H2OPUS(pc,pch2opus->wns[0],pch2opus->wns[1]));
4089566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(y,-1.,2.,pch2opus->wns[1]));
40953022affSStefano Zampini   }
41053022affSStefano Zampini   PetscFunctionReturn(0);
41153022affSStefano Zampini }
41253022affSStefano Zampini 
41353022affSStefano Zampini static PetscErrorCode MatMult_NS(Mat M, Vec x, Vec y)
41453022affSStefano Zampini {
41553022affSStefano Zampini   PetscFunctionBegin;
4169566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_NS(M,x,y,PETSC_FALSE));
41753022affSStefano Zampini   PetscFunctionReturn(0);
41853022affSStefano Zampini }
41953022affSStefano Zampini 
42053022affSStefano Zampini static PetscErrorCode MatMultTranspose_NS(Mat M, Vec x, Vec y)
42153022affSStefano Zampini {
42253022affSStefano Zampini   PetscFunctionBegin;
4239566063dSJacob Faibussowitsch   PetscCall(MatMultKernel_NS(M,x,y,PETSC_TRUE));
42453022affSStefano Zampini   PetscFunctionReturn(0);
42553022affSStefano Zampini }
42653022affSStefano Zampini 
42753022affSStefano Zampini /* Basic Newton-Schultz sampler: (2 * I - M * A)*M, MatMat version */
42853022affSStefano Zampini static PetscErrorCode MatMatMultKernel_NS(Mat M, Mat X, Mat Y, PetscBool t)
42953022affSStefano Zampini {
43053022affSStefano Zampini   PC         pc;
43153022affSStefano Zampini   Mat        A;
43253022affSStefano Zampini   PC_H2OPUS *pch2opus;
43353022affSStefano Zampini 
43453022affSStefano Zampini   PetscFunctionBegin;
4359566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M,(void*)&pc));
43653022affSStefano Zampini   pch2opus = (PC_H2OPUS*)pc->data;
43753022affSStefano Zampini   A = pch2opus->A;
43853022affSStefano Zampini   if (pch2opus->wnsmat[0] && pch2opus->wnsmat[0]->cmap->N != X->cmap->N) {
4399566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[0]));
4409566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pch2opus->wnsmat[1]));
44153022affSStefano Zampini   }
44253022affSStefano Zampini   if (!pch2opus->wnsmat[0]) {
4439566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(X,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[0]));
4449566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(Y,MAT_SHARE_NONZERO_PATTERN,&pch2opus->wnsmat[1]));
44553022affSStefano Zampini   }
44653022affSStefano Zampini   if (t) {
4479566063dSJacob Faibussowitsch     PetscCall(PCApplyTransposeMat_H2OPUS(pc,X,Y));
4489566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A,Y,MAT_REUSE_MATRIX,PETSC_DEFAULT,&pch2opus->wnsmat[1]));
4499566063dSJacob Faibussowitsch     PetscCall(PCApplyTransposeMat_H2OPUS(pc,pch2opus->wnsmat[1],pch2opus->wnsmat[0]));
4509566063dSJacob Faibussowitsch     PetscCall(MatScale(Y,2.));
4519566063dSJacob Faibussowitsch     PetscCall(MatAXPY(Y,-1.,pch2opus->wnsmat[0],SAME_NONZERO_PATTERN));
45253022affSStefano Zampini   } else {
4539566063dSJacob Faibussowitsch     PetscCall(PCApplyMat_H2OPUS(pc,X,Y));
4549566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A,Y,MAT_REUSE_MATRIX,PETSC_DEFAULT,&pch2opus->wnsmat[0]));
4559566063dSJacob Faibussowitsch     PetscCall(PCApplyMat_H2OPUS(pc,pch2opus->wnsmat[0],pch2opus->wnsmat[1]));
4569566063dSJacob Faibussowitsch     PetscCall(MatScale(Y,2.));
4579566063dSJacob Faibussowitsch     PetscCall(MatAXPY(Y,-1.,pch2opus->wnsmat[1],SAME_NONZERO_PATTERN));
45853022affSStefano Zampini   }
45953022affSStefano Zampini   PetscFunctionReturn(0);
46053022affSStefano Zampini }
46153022affSStefano Zampini 
46253022affSStefano Zampini static PetscErrorCode MatMatMultNumeric_NS(Mat M, Mat X, Mat Y, void *ctx)
46353022affSStefano Zampini {
46453022affSStefano Zampini   PetscFunctionBegin;
4659566063dSJacob Faibussowitsch   PetscCall(MatMatMultKernel_NS(M,X,Y,PETSC_FALSE));
46653022affSStefano Zampini   PetscFunctionReturn(0);
46753022affSStefano Zampini }
46853022affSStefano Zampini 
46953022affSStefano Zampini static PetscErrorCode PCH2OpusSetUpSampler_Private(PC pc)
47053022affSStefano Zampini {
47153022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
47253022affSStefano Zampini   Mat        A        = pc->useAmat ? pc->mat : pc->pmat;
47353022affSStefano Zampini 
47453022affSStefano Zampini   PetscFunctionBegin;
47553022affSStefano Zampini   if (!pch2opus->S) {
47653022affSStefano Zampini     PetscInt M,N,m,n;
47753022affSStefano Zampini 
4789566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A,&M,&N));
4799566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A,&m,&n));
4809566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)A),m,n,M,N,pc,&pch2opus->S));
4819566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(pch2opus->S,A,A));
48253022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU)
4839566063dSJacob Faibussowitsch     PetscCall(MatShellSetVecType(pch2opus->S,VECCUDA));
48453022affSStefano Zampini #endif
48553022affSStefano Zampini   }
48653022affSStefano Zampini   if (pch2opus->hyperorder >= 2) {
4879566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->S,MATOP_MULT,(void (*)(void))MatMult_Hyper));
4889566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_Hyper));
4899566063dSJacob Faibussowitsch     PetscCall(MatShellSetMatProductOperation(pch2opus->S,MATPRODUCT_AB,NULL,MatMatMultNumeric_Hyper,NULL,MATDENSE,MATDENSE));
4909566063dSJacob Faibussowitsch     PetscCall(MatShellSetMatProductOperation(pch2opus->S,MATPRODUCT_AB,NULL,MatMatMultNumeric_Hyper,NULL,MATDENSECUDA,MATDENSECUDA));
49153022affSStefano Zampini   } else {
4929566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->S,MATOP_MULT,(void (*)(void))MatMult_NS));
4939566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_NS));
4949566063dSJacob Faibussowitsch     PetscCall(MatShellSetMatProductOperation(pch2opus->S,MATPRODUCT_AB,NULL,MatMatMultNumeric_NS,NULL,MATDENSE,MATDENSE));
4959566063dSJacob Faibussowitsch     PetscCall(MatShellSetMatProductOperation(pch2opus->S,MATPRODUCT_AB,NULL,MatMatMultNumeric_NS,NULL,MATDENSECUDA,MATDENSECUDA));
49653022affSStefano Zampini   }
4979566063dSJacob Faibussowitsch   PetscCall(MatPropagateSymmetryOptions(A,pch2opus->S));
4989566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(pch2opus->S,pch2opus->boundtocpu));
49953022affSStefano Zampini   /* XXX */
5009566063dSJacob Faibussowitsch   PetscCall(MatSetOption(pch2opus->S,MAT_SYMMETRIC,PETSC_TRUE));
50153022affSStefano Zampini   PetscFunctionReturn(0);
50253022affSStefano Zampini }
50353022affSStefano Zampini 
50453022affSStefano Zampini static PetscErrorCode PCSetUp_H2OPUS(PC pc)
50553022affSStefano Zampini {
50653022affSStefano Zampini   PC_H2OPUS *pch2opus  = (PC_H2OPUS*)pc->data;
50753022affSStefano Zampini   Mat        A         = pc->useAmat ? pc->mat : pc->pmat;
50853022affSStefano Zampini   NormType   norm      = pch2opus->normtype;
50953022affSStefano Zampini   PetscReal  initerr   = 0.0,err;
51053022affSStefano Zampini   PetscReal  initerrMA = 0.0,errMA;
51153022affSStefano Zampini   PetscBool  ish2opus;
51253022affSStefano Zampini 
51353022affSStefano Zampini   PetscFunctionBegin;
51453022affSStefano Zampini   if (!pch2opus->T) {
51553022affSStefano Zampini     PetscInt    M,N,m,n;
51653022affSStefano Zampini     const char *prefix;
51753022affSStefano Zampini 
5189566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc,&prefix));
5199566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A,&M,&N));
5209566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A,&m,&n));
5219566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc->pmat),m,n,M,N,pc,&pch2opus->T));
5229566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(pch2opus->T,A,A));
5239566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->T,MATOP_MULT,(void (*)(void))MatMult_MAmI));
5249566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->T,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_MAmI));
5259566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(pch2opus->T,MATOP_NORM,(void (*)(void))MatNorm_H2OPUS));
52653022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU)
5279566063dSJacob Faibussowitsch     PetscCall(MatShellSetVecType(pch2opus->T,VECCUDA));
52853022affSStefano Zampini #endif
5299566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)pch2opus->T));
5309566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(pch2opus->T,prefix));
5319566063dSJacob Faibussowitsch     PetscCall(MatAppendOptionsPrefix(pch2opus->T,"pc_h2opus_est_"));
53253022affSStefano Zampini   }
5339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->A));
5349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATH2OPUS,&ish2opus));
53553022affSStefano Zampini   if (ish2opus) {
5369566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
53753022affSStefano Zampini     pch2opus->A = A;
53853022affSStefano Zampini   } else {
53953022affSStefano Zampini     const char *prefix;
5409566063dSJacob Faibussowitsch     PetscCall(MatCreateH2OpusFromMat(A,pch2opus->sdim,pch2opus->coords,PETSC_FALSE,pch2opus->eta,pch2opus->leafsize,pch2opus->max_rank,pch2opus->bs,pch2opus->mrtol,&pch2opus->A));
54153022affSStefano Zampini     /* XXX */
5429566063dSJacob Faibussowitsch     PetscCall(MatSetOption(pch2opus->A,MAT_SYMMETRIC,PETSC_TRUE));
5439566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc,&prefix));
5449566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(pch2opus->A,prefix));
5459566063dSJacob Faibussowitsch     PetscCall(MatAppendOptionsPrefix(pch2opus->A,"pc_h2opus_init_"));
5469566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(pch2opus->A));
5479566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(pch2opus->A,MAT_FINAL_ASSEMBLY));
5489566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(pch2opus->A,MAT_FINAL_ASSEMBLY));
54953022affSStefano Zampini     /* XXX */
5509566063dSJacob Faibussowitsch     PetscCall(MatSetOption(pch2opus->A,MAT_SYMMETRIC,PETSC_TRUE));
55153022affSStefano Zampini   }
55253022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU)
55353022affSStefano Zampini   pch2opus->boundtocpu = pch2opus->A->boundtocpu;
55453022affSStefano Zampini #endif
5559566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(pch2opus->T,pch2opus->boundtocpu));
55653022affSStefano Zampini   if (pch2opus->M) { /* see if we can reuse M as initial guess */
55753022affSStefano Zampini     PetscReal reuse;
55853022affSStefano Zampini 
5599566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(pch2opus->M,pch2opus->boundtocpu));
5609566063dSJacob Faibussowitsch     PetscCall(MatNorm(pch2opus->T,norm,&reuse));
5619566063dSJacob Faibussowitsch     if (reuse >= 1.0) PetscCall(MatDestroy(&pch2opus->M));
56253022affSStefano Zampini   }
56353022affSStefano Zampini   if (!pch2opus->M) {
56453022affSStefano Zampini     const char *prefix;
5659566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(pch2opus->A,MAT_COPY_VALUES,&pch2opus->M));
5669566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc,&prefix));
5679566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(pch2opus->M,prefix));
5689566063dSJacob Faibussowitsch     PetscCall(MatAppendOptionsPrefix(pch2opus->M,"pc_h2opus_inv_"));
5699566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(pch2opus->M));
5709566063dSJacob Faibussowitsch     PetscCall(PCH2OpusSetUpInit(pc));
5719566063dSJacob Faibussowitsch     PetscCall(MatScale(pch2opus->M,pch2opus->s0));
57253022affSStefano Zampini   }
57353022affSStefano Zampini   /* A and M have the same h2 matrix structure, save on reordering routines */
5749566063dSJacob Faibussowitsch   PetscCall(MatH2OpusSetNativeMult(pch2opus->A,PETSC_TRUE));
5759566063dSJacob Faibussowitsch   PetscCall(MatH2OpusSetNativeMult(pch2opus->M,PETSC_TRUE));
57653022affSStefano Zampini   if (norm == NORM_1 || norm == NORM_2 || norm == NORM_INFINITY) {
5779566063dSJacob Faibussowitsch     PetscCall(MatNorm(pch2opus->T,norm,&initerr));
57853022affSStefano Zampini     pch2opus->testMA = PETSC_TRUE;
5799566063dSJacob Faibussowitsch     PetscCall(MatNorm(pch2opus->T,norm,&initerrMA));
58053022affSStefano Zampini     pch2opus->testMA = PETSC_FALSE;
58153022affSStefano Zampini   }
58253022affSStefano Zampini   if (PetscIsInfOrNanReal(initerr)) pc->failedreason = PC_SETUP_ERROR;
58353022affSStefano Zampini   err   = initerr;
58453022affSStefano Zampini   errMA = initerrMA;
58553022affSStefano Zampini   if (pch2opus->monitor) {
5869566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pc),"%D: ||M*A - I|| NORM%s abs %g rel %g\n",0,NormTypes[norm],(double)err,(double)(err/initerr)));
5879566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pc),"%D: ||M*A||     NORM%s abs %g rel %g\n",0,NormTypes[norm],(double)errMA,(double)(errMA/initerrMA)));
58853022affSStefano Zampini   }
58953022affSStefano Zampini   if (initerr > pch2opus->atol && !pc->failedreason) {
59053022affSStefano Zampini     PetscInt i;
59153022affSStefano Zampini 
5929566063dSJacob Faibussowitsch     PetscCall(PCH2OpusSetUpSampler_Private(pc));
59353022affSStefano Zampini     for (i = 0; i < pch2opus->maxits; i++) {
59453022affSStefano Zampini       Mat         M;
59553022affSStefano Zampini       const char* prefix;
59653022affSStefano Zampini 
5979566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(pch2opus->M,MAT_SHARE_NONZERO_PATTERN,&M));
5989566063dSJacob Faibussowitsch       PetscCall(MatGetOptionsPrefix(pch2opus->M,&prefix));
5999566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(M,prefix));
6009566063dSJacob Faibussowitsch       PetscCall(MatH2OpusSetSamplingMat(M,pch2opus->S,PETSC_DECIDE,PETSC_DECIDE));
6019566063dSJacob Faibussowitsch       PetscCall(MatSetFromOptions(M));
6029566063dSJacob Faibussowitsch       PetscCall(MatH2OpusSetNativeMult(M,PETSC_TRUE));
6039566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
6049566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
60553022affSStefano Zampini 
6069566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pch2opus->M));
60753022affSStefano Zampini       pch2opus->M = M;
60853022affSStefano Zampini       if (norm == NORM_1 || norm == NORM_2 || norm == NORM_INFINITY) {
6099566063dSJacob Faibussowitsch         PetscCall(MatNorm(pch2opus->T,norm,&err));
61053022affSStefano Zampini         pch2opus->testMA = PETSC_TRUE;
6119566063dSJacob Faibussowitsch         PetscCall(MatNorm(pch2opus->T,norm,&errMA));
61253022affSStefano Zampini         pch2opus->testMA = PETSC_FALSE;
61353022affSStefano Zampini       }
61453022affSStefano Zampini       if (pch2opus->monitor) {
6159566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pc),"%D: ||M*A - I|| NORM%s abs %g rel %g\n",i+1,NormTypes[norm],(double)err,(double)(err/initerr)));
6169566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pc),"%D: ||M*A||     NORM%s abs %g rel %g\n",i+1,NormTypes[norm],(double)errMA,(double)(errMA/initerrMA)));
61753022affSStefano Zampini       }
61853022affSStefano Zampini       if (PetscIsInfOrNanReal(err)) pc->failedreason = PC_SETUP_ERROR;
61953022affSStefano Zampini       if (err < pch2opus->atol || err < pch2opus->rtol*initerr || pc->failedreason) break;
62053022affSStefano Zampini     }
62153022affSStefano Zampini   }
62253022affSStefano Zampini   /* cleanup setup workspace */
6239566063dSJacob Faibussowitsch   PetscCall(MatH2OpusSetNativeMult(pch2opus->A,PETSC_FALSE));
6249566063dSJacob Faibussowitsch   PetscCall(MatH2OpusSetNativeMult(pch2opus->M,PETSC_FALSE));
6259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[0]));
6269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[1]));
6279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[2]));
6289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pch2opus->wns[3]));
6299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[0]));
6309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[1]));
6319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[2]));
6329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pch2opus->wnsmat[3]));
63353022affSStefano Zampini   PetscFunctionReturn(0);
63453022affSStefano Zampini }
63553022affSStefano Zampini 
63653022affSStefano Zampini static PetscErrorCode PCView_H2OPUS(PC pc, PetscViewer viewer)
63753022affSStefano Zampini {
63853022affSStefano Zampini   PC_H2OPUS *pch2opus = (PC_H2OPUS*)pc->data;
63953022affSStefano Zampini   PetscBool  isascii;
64053022affSStefano Zampini 
64153022affSStefano Zampini   PetscFunctionBegin;
6429566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
64353022affSStefano Zampini   if (isascii) {
64453022affSStefano Zampini     if (pch2opus->A && pch2opus->A != pc->mat && pch2opus->A != pc->pmat) {
6459566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Initial approximation matrix\n"));
6469566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
6479566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
6489566063dSJacob Faibussowitsch       PetscCall(MatView(pch2opus->A,viewer));
6499566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
6509566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
65153022affSStefano Zampini     }
65253022affSStefano Zampini     if (pch2opus->M) {
6539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Inner matrix constructed\n"));
6549566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
6559566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
6569566063dSJacob Faibussowitsch       PetscCall(MatView(pch2opus->M,viewer));
6579566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
6589566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
65953022affSStefano Zampini     }
6609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Initial scaling: %g\n",pch2opus->s0));
66153022affSStefano Zampini   }
66253022affSStefano Zampini   PetscFunctionReturn(0);
66353022affSStefano Zampini }
66453022affSStefano Zampini 
66553022affSStefano Zampini PETSC_EXTERN PetscErrorCode PCCreate_H2OPUS(PC pc)
66653022affSStefano Zampini {
66753022affSStefano Zampini   PC_H2OPUS *pch2opus;
66853022affSStefano Zampini 
66953022affSStefano Zampini   PetscFunctionBegin;
6709566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&pch2opus));
67153022affSStefano Zampini   pc->data = (void*)pch2opus;
67253022affSStefano Zampini 
67353022affSStefano Zampini   pch2opus->atol       = 1.e-2;
67453022affSStefano Zampini   pch2opus->rtol       = 1.e-6;
67553022affSStefano Zampini   pch2opus->maxits     = 50;
67653022affSStefano Zampini   pch2opus->hyperorder = 1; /* defaults to basic NewtonSchultz */
67753022affSStefano Zampini   pch2opus->normtype   = NORM_2;
67853022affSStefano Zampini 
67953022affSStefano Zampini   /* these are needed when we are sampling the pmat */
68053022affSStefano Zampini   pch2opus->eta        = PETSC_DECIDE;
68153022affSStefano Zampini   pch2opus->leafsize   = PETSC_DECIDE;
68253022affSStefano Zampini   pch2opus->max_rank   = PETSC_DECIDE;
68353022affSStefano Zampini   pch2opus->bs         = PETSC_DECIDE;
68453022affSStefano Zampini   pch2opus->mrtol      = PETSC_DECIDE;
68553022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU)
68653022affSStefano Zampini   pch2opus->boundtocpu = PETSC_FALSE;
68753022affSStefano Zampini #else
68853022affSStefano Zampini   pch2opus->boundtocpu = PETSC_TRUE;
68953022affSStefano Zampini #endif
69053022affSStefano Zampini   pc->ops->destroy        = PCDestroy_H2OPUS;
69153022affSStefano Zampini   pc->ops->setup          = PCSetUp_H2OPUS;
69253022affSStefano Zampini   pc->ops->apply          = PCApply_H2OPUS;
69353022affSStefano Zampini   pc->ops->matapply       = PCApplyMat_H2OPUS;
69453022affSStefano Zampini   pc->ops->applytranspose = PCApplyTranspose_H2OPUS;
69553022affSStefano Zampini   pc->ops->reset          = PCReset_H2OPUS;
69653022affSStefano Zampini   pc->ops->setfromoptions = PCSetFromOptions_H2OPUS;
69753022affSStefano Zampini   pc->ops->view           = PCView_H2OPUS;
69853022affSStefano Zampini 
6999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_H2OPUS));
70053022affSStefano Zampini   PetscFunctionReturn(0);
70153022affSStefano Zampini }
702