153022affSStefano Zampini #include <petsc/private/pcimpl.h> 253022affSStefano Zampini #include <petsc/private/matimpl.h> 38ba4effdSStefano Zampini #include <petscdm.h> 453022affSStefano Zampini #include <h2opusconf.h> 553022affSStefano Zampini 653022affSStefano Zampini /* Use GPU only if H2OPUS is configured for GPU */ 753022affSStefano Zampini #if defined(PETSC_HAVE_CUDA) && defined(H2OPUS_USE_GPU) 853022affSStefano Zampini #define PETSC_H2OPUS_USE_GPU 953022affSStefano Zampini #endif 1053022affSStefano Zampini 1153022affSStefano Zampini typedef struct { 1253022affSStefano Zampini Mat A; 1353022affSStefano Zampini Mat M; 1453022affSStefano Zampini PetscScalar s0; 1553022affSStefano Zampini 1653022affSStefano Zampini /* sampler for Newton-Schultz */ 1753022affSStefano Zampini Mat S; 1853022affSStefano Zampini PetscInt hyperorder; 1953022affSStefano Zampini Vec wns[4]; 2053022affSStefano Zampini Mat wnsmat[4]; 2153022affSStefano Zampini 2253022affSStefano Zampini /* convergence testing */ 2353022affSStefano Zampini Mat T; 2453022affSStefano Zampini Vec w; 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 458db51263SStefano Zampini /* CPU/GPU */ 468db51263SStefano Zampini PetscBool forcecpu; 4753022affSStefano Zampini PetscBool boundtocpu; 4853022affSStefano Zampini } PC_H2OPUS; 4953022affSStefano Zampini 5053022affSStefano Zampini PETSC_EXTERN PetscErrorCode MatNorm_H2OPUS(Mat, NormType, PetscReal *); 5153022affSStefano Zampini 528ba4effdSStefano Zampini static PetscErrorCode PCH2OpusInferCoordinates_Private(PC pc) 538ba4effdSStefano Zampini { 548ba4effdSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 558ba4effdSStefano Zampini DM dm; 568ba4effdSStefano Zampini PetscBool isdmda; 578ba4effdSStefano Zampini 588ba4effdSStefano Zampini PetscFunctionBegin; 598ba4effdSStefano Zampini if (pch2opus->sdim) PetscFunctionReturn(PETSC_SUCCESS); 608ba4effdSStefano Zampini PetscCall(PCGetDM(pc, &dm)); 618ba4effdSStefano Zampini if (!dm) PetscCall(MatGetDM(pc->useAmat ? pc->mat : pc->pmat, &dm)); 628ba4effdSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isdmda)); 638ba4effdSStefano Zampini if (isdmda) { 648ba4effdSStefano Zampini Vec c; 658ba4effdSStefano Zampini const PetscScalar *coords; 668ba4effdSStefano Zampini PetscInt n, sdim; 678ba4effdSStefano Zampini 688ba4effdSStefano Zampini PetscCall(DMGetCoordinates(dm, &c)); 698ba4effdSStefano Zampini PetscCall(DMGetDimension(dm, &sdim)); 708ba4effdSStefano Zampini PetscCall(VecGetLocalSize(c, &n)); 718ba4effdSStefano Zampini PetscCall(VecGetArrayRead(c, &coords)); 728ba4effdSStefano Zampini PetscCall(PCSetCoordinates(pc, sdim, n / sdim, (PetscScalar *)coords)); 738ba4effdSStefano Zampini PetscCall(VecRestoreArrayRead(c, &coords)); 748ba4effdSStefano Zampini } 758ba4effdSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 768ba4effdSStefano Zampini } 778ba4effdSStefano Zampini 78d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_H2OPUS(PC pc) 79d71ae5a4SJacob Faibussowitsch { 8053022affSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 8153022affSStefano Zampini 8253022affSStefano Zampini PetscFunctionBegin; 8353022affSStefano Zampini pch2opus->sdim = 0; 8453022affSStefano Zampini pch2opus->nlocc = 0; 859566063dSJacob Faibussowitsch PetscCall(PetscFree(pch2opus->coords)); 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8753022affSStefano Zampini } 8853022affSStefano Zampini 89d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_H2OPUS(PC pc, PetscInt sdim, PetscInt nlocc, PetscReal *coords) 90d71ae5a4SJacob Faibussowitsch { 9153022affSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 9253022affSStefano Zampini PetscBool reset = PETSC_TRUE; 9353022affSStefano Zampini 9453022affSStefano Zampini PetscFunctionBegin; 9553022affSStefano Zampini if (pch2opus->sdim && sdim == pch2opus->sdim && nlocc == pch2opus->nlocc) { 969566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(pch2opus->coords, coords, sdim * nlocc, &reset)); 9753022affSStefano Zampini reset = (PetscBool)!reset; 9853022affSStefano Zampini } 991c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &reset, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)pc))); 10053022affSStefano Zampini if (reset) { 1019566063dSJacob Faibussowitsch PetscCall(PCReset_H2OPUS(pc)); 1029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(sdim * nlocc, &pch2opus->coords)); 1039566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(pch2opus->coords, coords, sdim * nlocc)); 10453022affSStefano Zampini pch2opus->sdim = sdim; 10553022affSStefano Zampini pch2opus->nlocc = nlocc; 10653022affSStefano Zampini } 1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10853022affSStefano Zampini } 10953022affSStefano Zampini 110d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_H2OPUS(PC pc) 111d71ae5a4SJacob Faibussowitsch { 1128ba4effdSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 1138ba4effdSStefano Zampini 11453022affSStefano Zampini PetscFunctionBegin; 1159566063dSJacob Faibussowitsch PetscCall(PCReset_H2OPUS(pc)); 1168ba4effdSStefano Zampini PetscCall(MatDestroy(&pch2opus->A)); 1178ba4effdSStefano Zampini PetscCall(MatDestroy(&pch2opus->M)); 1188ba4effdSStefano Zampini PetscCall(MatDestroy(&pch2opus->T)); 1198ba4effdSStefano Zampini PetscCall(VecDestroy(&pch2opus->w)); 1208ba4effdSStefano Zampini PetscCall(MatDestroy(&pch2opus->S)); 1218ba4effdSStefano Zampini PetscCall(VecDestroy(&pch2opus->wns[0])); 1228ba4effdSStefano Zampini PetscCall(VecDestroy(&pch2opus->wns[1])); 1238ba4effdSStefano Zampini PetscCall(VecDestroy(&pch2opus->wns[2])); 1248ba4effdSStefano Zampini PetscCall(VecDestroy(&pch2opus->wns[3])); 1258ba4effdSStefano Zampini PetscCall(MatDestroy(&pch2opus->wnsmat[0])); 1268ba4effdSStefano Zampini PetscCall(MatDestroy(&pch2opus->wnsmat[1])); 1278ba4effdSStefano Zampini PetscCall(MatDestroy(&pch2opus->wnsmat[2])); 1288ba4effdSStefano Zampini PetscCall(MatDestroy(&pch2opus->wnsmat[3])); 1299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13253022affSStefano Zampini } 13353022affSStefano Zampini 134d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_H2OPUS(PC pc, PetscOptionItems *PetscOptionsObject) 135d71ae5a4SJacob Faibussowitsch { 13653022affSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 13753022affSStefano Zampini 13853022affSStefano Zampini PetscFunctionBegin; 139d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "H2OPUS options"); 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_h2opus_maxits", "Maximum number of iterations for Newton-Schultz", NULL, pch2opus->maxits, &pch2opus->maxits, NULL)); 1419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_h2opus_monitor", "Monitor Newton-Schultz convergence", NULL, pch2opus->monitor, &pch2opus->monitor, NULL)); 1429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_h2opus_atol", "Absolute tolerance", NULL, pch2opus->atol, &pch2opus->atol, NULL)); 1439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_h2opus_rtol", "Relative tolerance", NULL, pch2opus->rtol, &pch2opus->rtol, NULL)); 1449566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_h2opus_norm_type", "Norm type for convergence monitoring", NULL, NormTypes, (PetscEnum)pch2opus->normtype, (PetscEnum *)&pch2opus->normtype, NULL)); 1459566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_h2opus_hyperorder", "Hyper power order of sampling", NULL, pch2opus->hyperorder, &pch2opus->hyperorder, NULL)); 1469566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_h2opus_leafsize", "Leaf size when constructed from kernel", NULL, pch2opus->leafsize, &pch2opus->leafsize, NULL)); 1479566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_h2opus_eta", "Admissibility condition tolerance", NULL, pch2opus->eta, &pch2opus->eta, NULL)); 1489566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_h2opus_maxrank", "Maximum rank when constructed from matvecs", NULL, pch2opus->max_rank, &pch2opus->max_rank, NULL)); 1499566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_h2opus_samples", "Number of samples to be taken concurrently when constructing from matvecs", NULL, pch2opus->bs, &pch2opus->bs, NULL)); 1509566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_h2opus_mrtol", "Relative tolerance for construction from sampling", NULL, pch2opus->mrtol, &pch2opus->mrtol, NULL)); 1518db51263SStefano Zampini PetscCall(PetscOptionsBool("-pc_h2opus_forcecpu", "Force construction on CPU", NULL, pch2opus->forcecpu, &pch2opus->forcecpu, NULL)); 152d0609cedSBarry Smith PetscOptionsHeadEnd(); 1533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15453022affSStefano Zampini } 15553022affSStefano Zampini 15653022affSStefano Zampini typedef struct { 15753022affSStefano Zampini Mat A; 15853022affSStefano Zampini Mat M; 15953022affSStefano Zampini Vec w; 16053022affSStefano Zampini } AAtCtx; 16153022affSStefano Zampini 162d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_AAt(Mat A, Vec x, Vec y) 163d71ae5a4SJacob Faibussowitsch { 16453022affSStefano Zampini AAtCtx *aat; 16553022affSStefano Zampini 16653022affSStefano Zampini PetscFunctionBegin; 1679566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, (void *)&aat)); 1689566063dSJacob Faibussowitsch /* PetscCall(MatMultTranspose(aat->M,x,aat->w)); */ 1699566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(aat->A, x, aat->w)); 1709566063dSJacob Faibussowitsch PetscCall(MatMult(aat->A, aat->w, y)); 1713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17253022affSStefano Zampini } 17353022affSStefano Zampini 174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCH2OpusSetUpInit(PC pc) 175d71ae5a4SJacob Faibussowitsch { 17653022affSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 17753022affSStefano Zampini Mat A = pc->useAmat ? pc->mat : pc->pmat, AAt; 17853022affSStefano Zampini PetscInt M, m; 17953022affSStefano Zampini VecType vtype; 18053022affSStefano Zampini PetscReal n; 18153022affSStefano Zampini AAtCtx aat; 18253022affSStefano Zampini 18353022affSStefano Zampini PetscFunctionBegin; 18453022affSStefano Zampini aat.A = A; 18553022affSStefano Zampini aat.M = pch2opus->M; /* unused so far */ 1869566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &aat.w)); 1879566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, NULL)); 1889566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 1899566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)A), m, m, M, M, &aat, &AAt)); 1909566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(AAt, pch2opus->boundtocpu)); 1919566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(AAt, MATOP_MULT, (void (*)(void))MatMult_AAt)); 1929566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(AAt, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_AAt)); 1939566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(AAt, MATOP_NORM, (void (*)(void))MatNorm_H2OPUS)); 1949566063dSJacob Faibussowitsch PetscCall(MatGetVecType(A, &vtype)); 1959566063dSJacob Faibussowitsch PetscCall(MatShellSetVecType(AAt, vtype)); 1969566063dSJacob Faibussowitsch PetscCall(MatNorm(AAt, NORM_1, &n)); 19753022affSStefano Zampini pch2opus->s0 = 1. / n; 1989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aat.w)); 1999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AAt)); 2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20153022affSStefano Zampini } 20253022affSStefano Zampini 203d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyKernel_H2OPUS(PC pc, Vec x, Vec y, PetscBool t) 204d71ae5a4SJacob Faibussowitsch { 20553022affSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 20653022affSStefano Zampini 20753022affSStefano Zampini PetscFunctionBegin; 2089566063dSJacob Faibussowitsch if (t) PetscCall(MatMultTranspose(pch2opus->M, x, y)); 2099566063dSJacob Faibussowitsch else PetscCall(MatMult(pch2opus->M, x, y)); 2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21153022affSStefano Zampini } 21253022affSStefano Zampini 213d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyMatKernel_H2OPUS(PC pc, Mat X, Mat Y, PetscBool t) 214d71ae5a4SJacob Faibussowitsch { 21553022affSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 21653022affSStefano Zampini 21753022affSStefano Zampini PetscFunctionBegin; 2189566063dSJacob Faibussowitsch if (t) PetscCall(MatTransposeMatMult(pch2opus->M, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 2199566063dSJacob Faibussowitsch else PetscCall(MatMatMult(pch2opus->M, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 2203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22153022affSStefano Zampini } 22253022affSStefano Zampini 223d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyMat_H2OPUS(PC pc, Mat X, Mat Y) 224d71ae5a4SJacob Faibussowitsch { 22553022affSStefano Zampini PetscFunctionBegin; 2269566063dSJacob Faibussowitsch PetscCall(PCApplyMatKernel_H2OPUS(pc, X, Y, PETSC_FALSE)); 2273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22853022affSStefano Zampini } 22953022affSStefano Zampini 230d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTransposeMat_H2OPUS(PC pc, Mat X, Mat Y) 231d71ae5a4SJacob Faibussowitsch { 23253022affSStefano Zampini PetscFunctionBegin; 2339566063dSJacob Faibussowitsch PetscCall(PCApplyMatKernel_H2OPUS(pc, X, Y, PETSC_TRUE)); 2343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23553022affSStefano Zampini } 23653022affSStefano Zampini 237d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_H2OPUS(PC pc, Vec x, Vec y) 238d71ae5a4SJacob Faibussowitsch { 23953022affSStefano Zampini PetscFunctionBegin; 2409566063dSJacob Faibussowitsch PetscCall(PCApplyKernel_H2OPUS(pc, x, y, PETSC_FALSE)); 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24253022affSStefano Zampini } 24353022affSStefano Zampini 244d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_H2OPUS(PC pc, Vec x, Vec y) 245d71ae5a4SJacob Faibussowitsch { 24653022affSStefano Zampini PetscFunctionBegin; 2479566063dSJacob Faibussowitsch PetscCall(PCApplyKernel_H2OPUS(pc, x, y, PETSC_TRUE)); 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24953022affSStefano Zampini } 25053022affSStefano Zampini 25153022affSStefano Zampini /* used to test the norm of (M^-1 A - I) */ 252d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultKernel_MAmI(Mat M, Vec x, Vec y, PetscBool t) 253d71ae5a4SJacob Faibussowitsch { 25453022affSStefano Zampini PC pc; 25553022affSStefano Zampini Mat A; 25653022affSStefano Zampini PC_H2OPUS *pch2opus; 25753022affSStefano Zampini PetscBool sideleft = PETSC_TRUE; 25853022affSStefano Zampini 25953022affSStefano Zampini PetscFunctionBegin; 2609566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(M, (void *)&pc)); 26153022affSStefano Zampini pch2opus = (PC_H2OPUS *)pc->data; 2629566063dSJacob Faibussowitsch if (!pch2opus->w) PetscCall(MatCreateVecs(pch2opus->M, &pch2opus->w, NULL)); 26353022affSStefano Zampini A = pch2opus->A; 2649566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(pch2opus->w, pch2opus->boundtocpu)); 26553022affSStefano Zampini if (t) { 26653022affSStefano Zampini if (sideleft) { 2679566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose_H2OPUS(pc, x, pch2opus->w)); 2689566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, pch2opus->w, y)); 26953022affSStefano Zampini } else { 2709566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, pch2opus->w)); 2719566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose_H2OPUS(pc, pch2opus->w, y)); 27253022affSStefano Zampini } 27353022affSStefano Zampini } else { 27453022affSStefano Zampini if (sideleft) { 2759566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, pch2opus->w)); 2769566063dSJacob Faibussowitsch PetscCall(PCApply_H2OPUS(pc, pch2opus->w, y)); 27753022affSStefano Zampini } else { 2789566063dSJacob Faibussowitsch PetscCall(PCApply_H2OPUS(pc, x, pch2opus->w)); 2799566063dSJacob Faibussowitsch PetscCall(MatMult(A, pch2opus->w, y)); 28053022affSStefano Zampini } 28153022affSStefano Zampini } 2828db51263SStefano Zampini PetscCall(VecAXPY(y, -1.0, x)); 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28453022affSStefano Zampini } 28553022affSStefano Zampini 286d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_MAmI(Mat A, Vec x, Vec y) 287d71ae5a4SJacob Faibussowitsch { 28853022affSStefano Zampini PetscFunctionBegin; 2899566063dSJacob Faibussowitsch PetscCall(MatMultKernel_MAmI(A, x, y, PETSC_FALSE)); 2903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29153022affSStefano Zampini } 29253022affSStefano Zampini 293d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_MAmI(Mat A, Vec x, Vec y) 294d71ae5a4SJacob Faibussowitsch { 29553022affSStefano Zampini PetscFunctionBegin; 2969566063dSJacob Faibussowitsch PetscCall(MatMultKernel_MAmI(A, x, y, PETSC_TRUE)); 2973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29853022affSStefano Zampini } 29953022affSStefano Zampini 30053022affSStefano Zampini /* HyperPower kernel: 30153022affSStefano Zampini Y = R = x 30253022affSStefano Zampini for i = 1 . . . l - 1 do 30353022affSStefano Zampini R = (I - A * Xk) * R 30453022affSStefano Zampini Y = Y + R 30553022affSStefano Zampini Y = Xk * Y 30653022affSStefano Zampini */ 307d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultKernel_Hyper(Mat M, Vec x, Vec y, PetscBool t) 308d71ae5a4SJacob Faibussowitsch { 30953022affSStefano Zampini PC pc; 31053022affSStefano Zampini Mat A; 31153022affSStefano Zampini PC_H2OPUS *pch2opus; 31253022affSStefano Zampini 31353022affSStefano Zampini PetscFunctionBegin; 3149566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(M, (void *)&pc)); 31553022affSStefano Zampini pch2opus = (PC_H2OPUS *)pc->data; 31653022affSStefano Zampini A = pch2opus->A; 3179566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pch2opus->M, pch2opus->wns[0] ? NULL : &pch2opus->wns[0], pch2opus->wns[1] ? NULL : &pch2opus->wns[1])); 3189566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pch2opus->M, pch2opus->wns[2] ? NULL : &pch2opus->wns[2], pch2opus->wns[3] ? NULL : &pch2opus->wns[3])); 3199566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(pch2opus->wns[0], pch2opus->boundtocpu)); 3209566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(pch2opus->wns[1], pch2opus->boundtocpu)); 3219566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(pch2opus->wns[2], pch2opus->boundtocpu)); 3229566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(pch2opus->wns[3], pch2opus->boundtocpu)); 3239566063dSJacob Faibussowitsch PetscCall(VecCopy(x, pch2opus->wns[0])); 3249566063dSJacob Faibussowitsch PetscCall(VecCopy(x, pch2opus->wns[3])); 32553022affSStefano Zampini if (t) { 3265f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < pch2opus->hyperorder - 1; i++) { 3279566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, pch2opus->wns[0], pch2opus->wns[1])); 3289566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose_H2OPUS(pc, pch2opus->wns[1], pch2opus->wns[2])); 3299566063dSJacob Faibussowitsch PetscCall(VecAXPY(pch2opus->wns[0], -1., pch2opus->wns[2])); 3309566063dSJacob Faibussowitsch PetscCall(VecAXPY(pch2opus->wns[3], 1., pch2opus->wns[0])); 33153022affSStefano Zampini } 3329566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose_H2OPUS(pc, pch2opus->wns[3], y)); 33353022affSStefano Zampini } else { 3345f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < pch2opus->hyperorder - 1; i++) { 3359566063dSJacob Faibussowitsch PetscCall(PCApply_H2OPUS(pc, pch2opus->wns[0], pch2opus->wns[1])); 3369566063dSJacob Faibussowitsch PetscCall(MatMult(A, pch2opus->wns[1], pch2opus->wns[2])); 3379566063dSJacob Faibussowitsch PetscCall(VecAXPY(pch2opus->wns[0], -1., pch2opus->wns[2])); 3389566063dSJacob Faibussowitsch PetscCall(VecAXPY(pch2opus->wns[3], 1., pch2opus->wns[0])); 33953022affSStefano Zampini } 3409566063dSJacob Faibussowitsch PetscCall(PCApply_H2OPUS(pc, pch2opus->wns[3], y)); 34153022affSStefano Zampini } 3423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34353022affSStefano Zampini } 34453022affSStefano Zampini 345d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_Hyper(Mat M, Vec x, Vec y) 346d71ae5a4SJacob Faibussowitsch { 34753022affSStefano Zampini PetscFunctionBegin; 3489566063dSJacob Faibussowitsch PetscCall(MatMultKernel_Hyper(M, x, y, PETSC_FALSE)); 3493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35053022affSStefano Zampini } 35153022affSStefano Zampini 352d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Hyper(Mat M, Vec x, Vec y) 353d71ae5a4SJacob Faibussowitsch { 35453022affSStefano Zampini PetscFunctionBegin; 3559566063dSJacob Faibussowitsch PetscCall(MatMultKernel_Hyper(M, x, y, PETSC_TRUE)); 3563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35753022affSStefano Zampini } 35853022affSStefano Zampini 35953022affSStefano Zampini /* Hyper power kernel, MatMat version */ 360d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatMultKernel_Hyper(Mat M, Mat X, Mat Y, PetscBool t) 361d71ae5a4SJacob Faibussowitsch { 36253022affSStefano Zampini PC pc; 36353022affSStefano Zampini Mat A; 36453022affSStefano Zampini PC_H2OPUS *pch2opus; 36553022affSStefano Zampini PetscInt i; 36653022affSStefano Zampini 36753022affSStefano Zampini PetscFunctionBegin; 3689566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(M, (void *)&pc)); 36953022affSStefano Zampini pch2opus = (PC_H2OPUS *)pc->data; 37053022affSStefano Zampini A = pch2opus->A; 37153022affSStefano Zampini if (pch2opus->wnsmat[0] && pch2opus->wnsmat[0]->cmap->N != X->cmap->N) { 3729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->wnsmat[0])); 3739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->wnsmat[1])); 37453022affSStefano Zampini } 37553022affSStefano Zampini if (!pch2opus->wnsmat[0]) { 3769566063dSJacob Faibussowitsch PetscCall(MatDuplicate(X, MAT_SHARE_NONZERO_PATTERN, &pch2opus->wnsmat[0])); 3779566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Y, MAT_SHARE_NONZERO_PATTERN, &pch2opus->wnsmat[1])); 37853022affSStefano Zampini } 37953022affSStefano Zampini if (pch2opus->wnsmat[2] && pch2opus->wnsmat[2]->cmap->N != X->cmap->N) { 3809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->wnsmat[2])); 3819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->wnsmat[3])); 38253022affSStefano Zampini } 38353022affSStefano Zampini if (!pch2opus->wnsmat[2]) { 3849566063dSJacob Faibussowitsch PetscCall(MatDuplicate(X, MAT_SHARE_NONZERO_PATTERN, &pch2opus->wnsmat[2])); 3859566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Y, MAT_SHARE_NONZERO_PATTERN, &pch2opus->wnsmat[3])); 38653022affSStefano Zampini } 3879566063dSJacob Faibussowitsch PetscCall(MatCopy(X, pch2opus->wnsmat[0], SAME_NONZERO_PATTERN)); 3889566063dSJacob Faibussowitsch PetscCall(MatCopy(X, pch2opus->wnsmat[3], SAME_NONZERO_PATTERN)); 38953022affSStefano Zampini if (t) { 39053022affSStefano Zampini for (i = 0; i < pch2opus->hyperorder - 1; i++) { 3919566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(A, pch2opus->wnsmat[0], MAT_REUSE_MATRIX, PETSC_DEFAULT, &pch2opus->wnsmat[1])); 3929566063dSJacob Faibussowitsch PetscCall(PCApplyTransposeMat_H2OPUS(pc, pch2opus->wnsmat[1], pch2opus->wnsmat[2])); 3939566063dSJacob Faibussowitsch PetscCall(MatAXPY(pch2opus->wnsmat[0], -1., pch2opus->wnsmat[2], SAME_NONZERO_PATTERN)); 3949566063dSJacob Faibussowitsch PetscCall(MatAXPY(pch2opus->wnsmat[3], 1., pch2opus->wnsmat[0], SAME_NONZERO_PATTERN)); 39553022affSStefano Zampini } 3969566063dSJacob Faibussowitsch PetscCall(PCApplyTransposeMat_H2OPUS(pc, pch2opus->wnsmat[3], Y)); 39753022affSStefano Zampini } else { 39853022affSStefano Zampini for (i = 0; i < pch2opus->hyperorder - 1; i++) { 3999566063dSJacob Faibussowitsch PetscCall(PCApplyMat_H2OPUS(pc, pch2opus->wnsmat[0], pch2opus->wnsmat[1])); 4009566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, pch2opus->wnsmat[1], MAT_REUSE_MATRIX, PETSC_DEFAULT, &pch2opus->wnsmat[2])); 4019566063dSJacob Faibussowitsch PetscCall(MatAXPY(pch2opus->wnsmat[0], -1., pch2opus->wnsmat[2], SAME_NONZERO_PATTERN)); 4029566063dSJacob Faibussowitsch PetscCall(MatAXPY(pch2opus->wnsmat[3], 1., pch2opus->wnsmat[0], SAME_NONZERO_PATTERN)); 40353022affSStefano Zampini } 4049566063dSJacob Faibussowitsch PetscCall(PCApplyMat_H2OPUS(pc, pch2opus->wnsmat[3], Y)); 40553022affSStefano Zampini } 4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40753022affSStefano Zampini } 40853022affSStefano Zampini 409d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatMultNumeric_Hyper(Mat M, Mat X, Mat Y, void *ctx) 410d71ae5a4SJacob Faibussowitsch { 41153022affSStefano Zampini PetscFunctionBegin; 4129566063dSJacob Faibussowitsch PetscCall(MatMatMultKernel_Hyper(M, X, Y, PETSC_FALSE)); 4133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41453022affSStefano Zampini } 41553022affSStefano Zampini 41653022affSStefano Zampini /* Basic Newton-Schultz sampler: (2 * I - M * A)*M */ 417d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultKernel_NS(Mat M, Vec x, Vec y, PetscBool t) 418d71ae5a4SJacob Faibussowitsch { 41953022affSStefano Zampini PC pc; 42053022affSStefano Zampini Mat A; 42153022affSStefano Zampini PC_H2OPUS *pch2opus; 42253022affSStefano Zampini 42353022affSStefano Zampini PetscFunctionBegin; 4249566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(M, (void *)&pc)); 42553022affSStefano Zampini pch2opus = (PC_H2OPUS *)pc->data; 42653022affSStefano Zampini A = pch2opus->A; 4279566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pch2opus->M, pch2opus->wns[0] ? NULL : &pch2opus->wns[0], pch2opus->wns[1] ? NULL : &pch2opus->wns[1])); 4289566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(pch2opus->wns[0], pch2opus->boundtocpu)); 4299566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(pch2opus->wns[1], pch2opus->boundtocpu)); 43053022affSStefano Zampini if (t) { 4319566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose_H2OPUS(pc, x, y)); 4329566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, y, pch2opus->wns[1])); 4339566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose_H2OPUS(pc, pch2opus->wns[1], pch2opus->wns[0])); 4349566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y, -1., 2., pch2opus->wns[0])); 43553022affSStefano Zampini } else { 4369566063dSJacob Faibussowitsch PetscCall(PCApply_H2OPUS(pc, x, y)); 4379566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, pch2opus->wns[0])); 4389566063dSJacob Faibussowitsch PetscCall(PCApply_H2OPUS(pc, pch2opus->wns[0], pch2opus->wns[1])); 4399566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y, -1., 2., pch2opus->wns[1])); 44053022affSStefano Zampini } 4413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44253022affSStefano Zampini } 44353022affSStefano Zampini 444d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_NS(Mat M, Vec x, Vec y) 445d71ae5a4SJacob Faibussowitsch { 44653022affSStefano Zampini PetscFunctionBegin; 4479566063dSJacob Faibussowitsch PetscCall(MatMultKernel_NS(M, x, y, PETSC_FALSE)); 4483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44953022affSStefano Zampini } 45053022affSStefano Zampini 451d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_NS(Mat M, Vec x, Vec y) 452d71ae5a4SJacob Faibussowitsch { 45353022affSStefano Zampini PetscFunctionBegin; 4549566063dSJacob Faibussowitsch PetscCall(MatMultKernel_NS(M, x, y, PETSC_TRUE)); 4553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45653022affSStefano Zampini } 45753022affSStefano Zampini 45853022affSStefano Zampini /* Basic Newton-Schultz sampler: (2 * I - M * A)*M, MatMat version */ 459d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatMultKernel_NS(Mat M, Mat X, Mat Y, PetscBool t) 460d71ae5a4SJacob Faibussowitsch { 46153022affSStefano Zampini PC pc; 46253022affSStefano Zampini Mat A; 46353022affSStefano Zampini PC_H2OPUS *pch2opus; 46453022affSStefano Zampini 46553022affSStefano Zampini PetscFunctionBegin; 4669566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(M, (void *)&pc)); 46753022affSStefano Zampini pch2opus = (PC_H2OPUS *)pc->data; 46853022affSStefano Zampini A = pch2opus->A; 46953022affSStefano Zampini if (pch2opus->wnsmat[0] && pch2opus->wnsmat[0]->cmap->N != X->cmap->N) { 4709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->wnsmat[0])); 4719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->wnsmat[1])); 47253022affSStefano Zampini } 47353022affSStefano Zampini if (!pch2opus->wnsmat[0]) { 4749566063dSJacob Faibussowitsch PetscCall(MatDuplicate(X, MAT_SHARE_NONZERO_PATTERN, &pch2opus->wnsmat[0])); 4759566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Y, MAT_SHARE_NONZERO_PATTERN, &pch2opus->wnsmat[1])); 47653022affSStefano Zampini } 47753022affSStefano Zampini if (t) { 4789566063dSJacob Faibussowitsch PetscCall(PCApplyTransposeMat_H2OPUS(pc, X, Y)); 4799566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(A, Y, MAT_REUSE_MATRIX, PETSC_DEFAULT, &pch2opus->wnsmat[1])); 4809566063dSJacob Faibussowitsch PetscCall(PCApplyTransposeMat_H2OPUS(pc, pch2opus->wnsmat[1], pch2opus->wnsmat[0])); 4819566063dSJacob Faibussowitsch PetscCall(MatScale(Y, 2.)); 4829566063dSJacob Faibussowitsch PetscCall(MatAXPY(Y, -1., pch2opus->wnsmat[0], SAME_NONZERO_PATTERN)); 48353022affSStefano Zampini } else { 4849566063dSJacob Faibussowitsch PetscCall(PCApplyMat_H2OPUS(pc, X, Y)); 4859566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, Y, MAT_REUSE_MATRIX, PETSC_DEFAULT, &pch2opus->wnsmat[0])); 4869566063dSJacob Faibussowitsch PetscCall(PCApplyMat_H2OPUS(pc, pch2opus->wnsmat[0], pch2opus->wnsmat[1])); 4879566063dSJacob Faibussowitsch PetscCall(MatScale(Y, 2.)); 4889566063dSJacob Faibussowitsch PetscCall(MatAXPY(Y, -1., pch2opus->wnsmat[1], SAME_NONZERO_PATTERN)); 48953022affSStefano Zampini } 4903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49153022affSStefano Zampini } 49253022affSStefano Zampini 493d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatMultNumeric_NS(Mat M, Mat X, Mat Y, void *ctx) 494d71ae5a4SJacob Faibussowitsch { 49553022affSStefano Zampini PetscFunctionBegin; 4969566063dSJacob Faibussowitsch PetscCall(MatMatMultKernel_NS(M, X, Y, PETSC_FALSE)); 4973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49853022affSStefano Zampini } 49953022affSStefano Zampini 500d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCH2OpusSetUpSampler_Private(PC pc) 501d71ae5a4SJacob Faibussowitsch { 50253022affSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 50353022affSStefano Zampini Mat A = pc->useAmat ? pc->mat : pc->pmat; 50453022affSStefano Zampini 50553022affSStefano Zampini PetscFunctionBegin; 50653022affSStefano Zampini if (!pch2opus->S) { 50753022affSStefano Zampini PetscInt M, N, m, n; 50853022affSStefano Zampini 5099566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 5109566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 5119566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)A), m, n, M, N, pc, &pch2opus->S)); 5129566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizesFromMats(pch2opus->S, A, A)); 51353022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU) 5149566063dSJacob Faibussowitsch PetscCall(MatShellSetVecType(pch2opus->S, VECCUDA)); 51553022affSStefano Zampini #endif 51653022affSStefano Zampini } 51753022affSStefano Zampini if (pch2opus->hyperorder >= 2) { 5189566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(pch2opus->S, MATOP_MULT, (void (*)(void))MatMult_Hyper)); 5199566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(pch2opus->S, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Hyper)); 5209566063dSJacob Faibussowitsch PetscCall(MatShellSetMatProductOperation(pch2opus->S, MATPRODUCT_AB, NULL, MatMatMultNumeric_Hyper, NULL, MATDENSE, MATDENSE)); 5219566063dSJacob Faibussowitsch PetscCall(MatShellSetMatProductOperation(pch2opus->S, MATPRODUCT_AB, NULL, MatMatMultNumeric_Hyper, NULL, MATDENSECUDA, MATDENSECUDA)); 52253022affSStefano Zampini } else { 5239566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(pch2opus->S, MATOP_MULT, (void (*)(void))MatMult_NS)); 5249566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(pch2opus->S, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_NS)); 5259566063dSJacob Faibussowitsch PetscCall(MatShellSetMatProductOperation(pch2opus->S, MATPRODUCT_AB, NULL, MatMatMultNumeric_NS, NULL, MATDENSE, MATDENSE)); 5269566063dSJacob Faibussowitsch PetscCall(MatShellSetMatProductOperation(pch2opus->S, MATPRODUCT_AB, NULL, MatMatMultNumeric_NS, NULL, MATDENSECUDA, MATDENSECUDA)); 52753022affSStefano Zampini } 5289566063dSJacob Faibussowitsch PetscCall(MatPropagateSymmetryOptions(A, pch2opus->S)); 5299566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(pch2opus->S, pch2opus->boundtocpu)); 53053022affSStefano Zampini /* XXX */ 5319566063dSJacob Faibussowitsch PetscCall(MatSetOption(pch2opus->S, MAT_SYMMETRIC, PETSC_TRUE)); 5323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53353022affSStefano Zampini } 53453022affSStefano Zampini 535d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_H2OPUS(PC pc) 536d71ae5a4SJacob Faibussowitsch { 53753022affSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 53853022affSStefano Zampini Mat A = pc->useAmat ? pc->mat : pc->pmat; 53953022affSStefano Zampini NormType norm = pch2opus->normtype; 54053022affSStefano Zampini PetscReal initerr = 0.0, err; 54153022affSStefano Zampini PetscBool ish2opus; 54253022affSStefano Zampini 54353022affSStefano Zampini PetscFunctionBegin; 54453022affSStefano Zampini if (!pch2opus->T) { 54553022affSStefano Zampini PetscInt M, N, m, n; 54653022affSStefano Zampini const char *prefix; 54753022affSStefano Zampini 5489566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 5499566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 5509566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 5519566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc->pmat), m, n, M, N, pc, &pch2opus->T)); 5529566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizesFromMats(pch2opus->T, A, A)); 5539566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(pch2opus->T, MATOP_MULT, (void (*)(void))MatMult_MAmI)); 5549566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(pch2opus->T, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_MAmI)); 5559566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(pch2opus->T, MATOP_NORM, (void (*)(void))MatNorm_H2OPUS)); 55653022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU) 5579566063dSJacob Faibussowitsch PetscCall(MatShellSetVecType(pch2opus->T, VECCUDA)); 55853022affSStefano Zampini #endif 5599566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(pch2opus->T, prefix)); 5609566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(pch2opus->T, "pc_h2opus_est_")); 56153022affSStefano Zampini } 5629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->A)); 5639566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus)); 56453022affSStefano Zampini if (ish2opus) { 5659566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 56653022affSStefano Zampini pch2opus->A = A; 56753022affSStefano Zampini } else { 56853022affSStefano Zampini const char *prefix; 5698ba4effdSStefano Zampini 5708ba4effdSStefano Zampini /* See if we can infer coordinates from the DM */ 5718ba4effdSStefano Zampini if (!pch2opus->sdim) PetscCall(PCH2OpusInferCoordinates_Private(pc)); 5729566063dSJacob Faibussowitsch PetscCall(MatCreateH2OpusFromMat(A, pch2opus->sdim, pch2opus->coords, PETSC_FALSE, pch2opus->eta, pch2opus->leafsize, pch2opus->max_rank, pch2opus->bs, pch2opus->mrtol, &pch2opus->A)); 57353022affSStefano Zampini /* XXX */ 5749566063dSJacob Faibussowitsch PetscCall(MatSetOption(pch2opus->A, MAT_SYMMETRIC, PETSC_TRUE)); 5759566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 5769566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(pch2opus->A, prefix)); 5779566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(pch2opus->A, "pc_h2opus_init_")); 5789566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pch2opus->A)); 5799566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pch2opus->A, MAT_FINAL_ASSEMBLY)); 5809566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pch2opus->A, MAT_FINAL_ASSEMBLY)); 58153022affSStefano Zampini /* XXX */ 5829566063dSJacob Faibussowitsch PetscCall(MatSetOption(pch2opus->A, MAT_SYMMETRIC, PETSC_TRUE)); 5838db51263SStefano Zampini 5848db51263SStefano Zampini /* always perform construction on the GPU unless forcecpu is true */ 5858db51263SStefano Zampini PetscCall(MatBindToCPU(pch2opus->A, pch2opus->forcecpu)); 58653022affSStefano Zampini } 58753022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU) 5888db51263SStefano Zampini pch2opus->boundtocpu = pch2opus->forcecpu ? PETSC_TRUE : pch2opus->A->boundtocpu; 58953022affSStefano Zampini #endif 5909566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(pch2opus->T, pch2opus->boundtocpu)); 59153022affSStefano Zampini if (pch2opus->M) { /* see if we can reuse M as initial guess */ 59253022affSStefano Zampini PetscReal reuse; 59353022affSStefano Zampini 5949566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(pch2opus->M, pch2opus->boundtocpu)); 5959566063dSJacob Faibussowitsch PetscCall(MatNorm(pch2opus->T, norm, &reuse)); 5969566063dSJacob Faibussowitsch if (reuse >= 1.0) PetscCall(MatDestroy(&pch2opus->M)); 59753022affSStefano Zampini } 59853022affSStefano Zampini if (!pch2opus->M) { 59953022affSStefano Zampini const char *prefix; 6009566063dSJacob Faibussowitsch PetscCall(MatDuplicate(pch2opus->A, MAT_COPY_VALUES, &pch2opus->M)); 6019566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 6029566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(pch2opus->M, prefix)); 6039566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(pch2opus->M, "pc_h2opus_inv_")); 6049566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pch2opus->M)); 6059566063dSJacob Faibussowitsch PetscCall(PCH2OpusSetUpInit(pc)); 6069566063dSJacob Faibussowitsch PetscCall(MatScale(pch2opus->M, pch2opus->s0)); 60753022affSStefano Zampini } 60853022affSStefano Zampini /* A and M have the same h2 matrix structure, save on reordering routines */ 6099566063dSJacob Faibussowitsch PetscCall(MatH2OpusSetNativeMult(pch2opus->A, PETSC_TRUE)); 6109566063dSJacob Faibussowitsch PetscCall(MatH2OpusSetNativeMult(pch2opus->M, PETSC_TRUE)); 6118db51263SStefano Zampini if (norm == NORM_1 || norm == NORM_2 || norm == NORM_INFINITY) PetscCall(MatNorm(pch2opus->T, norm, &initerr)); 61253022affSStefano Zampini if (PetscIsInfOrNanReal(initerr)) pc->failedreason = PC_SETUP_ERROR; 61353022affSStefano Zampini err = initerr; 61448a46eb9SPierre Jolivet if (pch2opus->monitor) 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))); 61553022affSStefano Zampini if (initerr > pch2opus->atol && !pc->failedreason) { 61653022affSStefano Zampini PetscInt i; 61753022affSStefano Zampini 6189566063dSJacob Faibussowitsch PetscCall(PCH2OpusSetUpSampler_Private(pc)); 61953022affSStefano Zampini for (i = 0; i < pch2opus->maxits; i++) { 62053022affSStefano Zampini Mat M; 62153022affSStefano Zampini const char *prefix; 62253022affSStefano Zampini 6239566063dSJacob Faibussowitsch PetscCall(MatDuplicate(pch2opus->M, MAT_SHARE_NONZERO_PATTERN, &M)); 6249566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(pch2opus->M, &prefix)); 6259566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(M, prefix)); 6269566063dSJacob Faibussowitsch PetscCall(MatH2OpusSetSamplingMat(M, pch2opus->S, PETSC_DECIDE, PETSC_DECIDE)); 6279566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(M)); 6289566063dSJacob Faibussowitsch PetscCall(MatH2OpusSetNativeMult(M, PETSC_TRUE)); 6299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY)); 6309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY)); 6319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->M)); 63253022affSStefano Zampini pch2opus->M = M; 6338db51263SStefano Zampini if (norm == NORM_1 || norm == NORM_2 || norm == NORM_INFINITY) PetscCall(MatNorm(pch2opus->T, norm, &err)); 63448a46eb9SPierre Jolivet if (pch2opus->monitor) 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))); 63553022affSStefano Zampini if (PetscIsInfOrNanReal(err)) pc->failedreason = PC_SETUP_ERROR; 63653022affSStefano Zampini if (err < pch2opus->atol || err < pch2opus->rtol * initerr || pc->failedreason) break; 63753022affSStefano Zampini } 63853022affSStefano Zampini } 63953022affSStefano Zampini /* cleanup setup workspace */ 6409566063dSJacob Faibussowitsch PetscCall(MatH2OpusSetNativeMult(pch2opus->A, PETSC_FALSE)); 6419566063dSJacob Faibussowitsch PetscCall(MatH2OpusSetNativeMult(pch2opus->M, PETSC_FALSE)); 6429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pch2opus->wns[0])); 6439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pch2opus->wns[1])); 6449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pch2opus->wns[2])); 6459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pch2opus->wns[3])); 6469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->wnsmat[0])); 6479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->wnsmat[1])); 6489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->wnsmat[2])); 6499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pch2opus->wnsmat[3])); 6503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65153022affSStefano Zampini } 65253022affSStefano Zampini 653d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_H2OPUS(PC pc, PetscViewer viewer) 654d71ae5a4SJacob Faibussowitsch { 65553022affSStefano Zampini PC_H2OPUS *pch2opus = (PC_H2OPUS *)pc->data; 65653022affSStefano Zampini PetscBool isascii; 65753022affSStefano Zampini 65853022affSStefano Zampini PetscFunctionBegin; 6599566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 66053022affSStefano Zampini if (isascii) { 66153022affSStefano Zampini if (pch2opus->A && pch2opus->A != pc->mat && pch2opus->A != pc->pmat) { 6629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Initial approximation matrix\n")); 6639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 6649566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL)); 6659566063dSJacob Faibussowitsch PetscCall(MatView(pch2opus->A, viewer)); 6669566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 6679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 66853022affSStefano Zampini } 66953022affSStefano Zampini if (pch2opus->M) { 6709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Inner matrix constructed\n")); 6719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 6729566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL)); 6739566063dSJacob Faibussowitsch PetscCall(MatView(pch2opus->M, viewer)); 6749566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 6759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 67653022affSStefano Zampini } 6779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Initial scaling: %g\n", pch2opus->s0)); 67853022affSStefano Zampini } 6793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68053022affSStefano Zampini } 68153022affSStefano Zampini 682f1580f4eSBarry Smith /*MC 683f1580f4eSBarry Smith PCH2OPUS = "h2opus" - A preconditioner type for, `MATH2OPUS`, hierarchical matrices using the H2Opus package. 684f1580f4eSBarry Smith 685f1580f4eSBarry Smith Options Database Keys: 686f1580f4eSBarry Smith + -pc_type h2opus - pc type to "h2opus" during a call to `PCSetFromOptions()` 687f1580f4eSBarry Smith . -pc_h2opus_maxits - maximum number of iterations for Newton-Schultz 688f1580f4eSBarry Smith . -pc_h2opus_monitor - monitor Newton-Schultz convergence 689f1580f4eSBarry Smith . -pc_h2opus_atol - absolute tolerance 690f1580f4eSBarry Smith . -pc_h2opus_rtol - relative tolerance 691f1580f4eSBarry Smith . -pc_h2opus_norm_type - normtype 692f1580f4eSBarry Smith . -pc_h2opus_hyperorder - Hyper power order of sampling 693f1580f4eSBarry Smith . -pc_h2opus_leafsize - leaf size when constructed from kernel 694f1580f4eSBarry Smith . -pc_h2opus_eta - admissibility condition tolerance 695f1580f4eSBarry Smith . -pc_h2opus_maxrank - maximum rank when constructed from matvecs 696f1580f4eSBarry Smith . -pc_h2opus_samples - number of samples to be taken concurrently when constructing from matvecs 697f1580f4eSBarry Smith . -pc_h2opus_mrtol - relative tolerance for construction from sampling 698f1580f4eSBarry Smith - -pc_h2opus_forcecpu - force construction of preconditioner on CPU 699f1580f4eSBarry Smith 700f1580f4eSBarry Smith Level: intermediate 701f1580f4eSBarry Smith 702*562efe2eSBarry Smith .seealso: [](ch_ksp), `MATH2OPUS`, `MATHTOOL`, `MATDENSE`, `MatCreateH2OpusFromKernel()`, `MatCreateH2OpusFromMat()` 703f1580f4eSBarry Smith M*/ 704f1580f4eSBarry Smith 705d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_H2OPUS(PC pc) 706d71ae5a4SJacob Faibussowitsch { 70753022affSStefano Zampini PC_H2OPUS *pch2opus; 70853022affSStefano Zampini 70953022affSStefano Zampini PetscFunctionBegin; 7104dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pch2opus)); 71153022affSStefano Zampini pc->data = (void *)pch2opus; 71253022affSStefano Zampini 71353022affSStefano Zampini pch2opus->atol = 1.e-2; 71453022affSStefano Zampini pch2opus->rtol = 1.e-6; 71553022affSStefano Zampini pch2opus->maxits = 50; 716f1580f4eSBarry Smith pch2opus->hyperorder = 1; /* defaults to basic Newton-Schultz */ 71753022affSStefano Zampini pch2opus->normtype = NORM_2; 71853022affSStefano Zampini 71953022affSStefano Zampini /* these are needed when we are sampling the pmat */ 72053022affSStefano Zampini pch2opus->eta = PETSC_DECIDE; 72153022affSStefano Zampini pch2opus->leafsize = PETSC_DECIDE; 72253022affSStefano Zampini pch2opus->max_rank = PETSC_DECIDE; 72353022affSStefano Zampini pch2opus->bs = PETSC_DECIDE; 72453022affSStefano Zampini pch2opus->mrtol = PETSC_DECIDE; 72553022affSStefano Zampini #if defined(PETSC_H2OPUS_USE_GPU) 72653022affSStefano Zampini pch2opus->boundtocpu = PETSC_FALSE; 72753022affSStefano Zampini #else 72853022affSStefano Zampini pch2opus->boundtocpu = PETSC_TRUE; 72953022affSStefano Zampini #endif 73053022affSStefano Zampini pc->ops->destroy = PCDestroy_H2OPUS; 73153022affSStefano Zampini pc->ops->setup = PCSetUp_H2OPUS; 73253022affSStefano Zampini pc->ops->apply = PCApply_H2OPUS; 73353022affSStefano Zampini pc->ops->matapply = PCApplyMat_H2OPUS; 73453022affSStefano Zampini pc->ops->applytranspose = PCApplyTranspose_H2OPUS; 73553022affSStefano Zampini pc->ops->reset = PCReset_H2OPUS; 73653022affSStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_H2OPUS; 73753022affSStefano Zampini pc->ops->view = PCView_H2OPUS; 73853022affSStefano Zampini 7399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_H2OPUS)); 7403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74153022affSStefano Zampini } 742