15e5bbd0aSStefano Zampini #include <petsc/private/pcbddcimpl.h> 25e5bbd0aSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h> 392cccca0SStefano Zampini #include <../src/mat/impls/dense/seq/dense.h> 4674ae819SStefano Zampini 592cccca0SStefano Zampini /* E + small_solve */ 69371c9d4SSatish Balay static PetscErrorCode PCBDDCNullSpaceCorrPreSolve(KSP ksp, Vec y, Vec x, void *ctx) { 792cccca0SStefano Zampini NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 892cccca0SStefano Zampini Mat K; 9bd5e1169SStefano Zampini 10bd5e1169SStefano Zampini PetscFunctionBegin; 119566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(corr_ctx->evapply, ksp, 0, 0, 0)); 129566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(corr_ctx->basis_mat, y, corr_ctx->sw[0])); 137ebab0bbSStefano Zampini if (corr_ctx->symm) { 149566063dSJacob Faibussowitsch PetscCall(MatMult(corr_ctx->inv_smat, corr_ctx->sw[0], corr_ctx->sw[1])); 157ebab0bbSStefano Zampini } else { 169566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(corr_ctx->inv_smat, corr_ctx->sw[0], corr_ctx->sw[1])); 177ebab0bbSStefano Zampini } 189566063dSJacob Faibussowitsch PetscCall(VecScale(corr_ctx->sw[1], -1.0)); 199566063dSJacob Faibussowitsch PetscCall(MatMult(corr_ctx->basis_mat, corr_ctx->sw[1], corr_ctx->fw[0])); 209566063dSJacob Faibussowitsch PetscCall(VecScale(corr_ctx->sw[1], -1.0)); 219566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(ksp, &K, NULL)); 229566063dSJacob Faibussowitsch PetscCall(MatMultAdd(K, corr_ctx->fw[0], y, y)); 239566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(corr_ctx->evapply, ksp, 0, 0, 0)); 24bd5e1169SStefano Zampini PetscFunctionReturn(0); 25bd5e1169SStefano Zampini } 26bd5e1169SStefano Zampini 2792cccca0SStefano Zampini /* E^t + small */ 289371c9d4SSatish Balay static PetscErrorCode PCBDDCNullSpaceCorrPostSolve(KSP ksp, Vec y, Vec x, void *ctx) { 2992cccca0SStefano Zampini NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 3092cccca0SStefano Zampini Mat K; 31674ae819SStefano Zampini 32674ae819SStefano Zampini PetscFunctionBegin; 339566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(corr_ctx->evapply, ksp, 0, 0, 0)); 349566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(ksp, &K, NULL)); 357ebab0bbSStefano Zampini if (corr_ctx->symm) { 369566063dSJacob Faibussowitsch PetscCall(MatMult(K, x, corr_ctx->fw[0])); 377ebab0bbSStefano Zampini } else { 389566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(K, x, corr_ctx->fw[0])); 397ebab0bbSStefano Zampini } 409566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(corr_ctx->basis_mat, corr_ctx->fw[0], corr_ctx->sw[0])); 419566063dSJacob Faibussowitsch PetscCall(VecScale(corr_ctx->sw[0], -1.0)); 429566063dSJacob Faibussowitsch PetscCall(MatMult(corr_ctx->inv_smat, corr_ctx->sw[0], corr_ctx->sw[2])); 439566063dSJacob Faibussowitsch PetscCall(MatMultAdd(corr_ctx->basis_mat, corr_ctx->sw[2], x, corr_ctx->fw[0])); 449566063dSJacob Faibussowitsch PetscCall(VecScale(corr_ctx->fw[0], corr_ctx->scale)); 4592cccca0SStefano Zampini /* Sum contributions from approximate solver and projected system */ 469566063dSJacob Faibussowitsch PetscCall(MatMultAdd(corr_ctx->basis_mat, corr_ctx->sw[1], corr_ctx->fw[0], x)); 479566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(corr_ctx->evapply, ksp, 0, 0, 0)); 48674ae819SStefano Zampini PetscFunctionReturn(0); 49674ae819SStefano Zampini } 50674ae819SStefano Zampini 519371c9d4SSatish Balay static PetscErrorCode PCBDDCNullSpaceCorrDestroy(void *ctx) { 5292cccca0SStefano Zampini NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 53674ae819SStefano Zampini 54674ae819SStefano Zampini PetscFunctionBegin; 559566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(3, &corr_ctx->sw)); 569566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(1, &corr_ctx->fw)); 579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&corr_ctx->basis_mat)); 589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&corr_ctx->inv_smat)); 599566063dSJacob Faibussowitsch PetscCall(PetscFree(corr_ctx)); 60674ae819SStefano Zampini PetscFunctionReturn(0); 61674ae819SStefano Zampini } 62674ae819SStefano Zampini 639371c9d4SSatish Balay PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling) { 64674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 65c7017625SStefano Zampini MatNullSpace NullSpace = NULL; 66185763b3SStefano Zampini KSP local_ksp; 67674ae819SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 6892cccca0SStefano Zampini Mat local_mat, local_pmat, dmat, Kbasis_mat; 6992cccca0SStefano Zampini Vec v; 7092cccca0SStefano Zampini PetscContainer c; 7192cccca0SStefano Zampini PetscInt basis_size; 7292cccca0SStefano Zampini IS zerorows; 7380fdaca0SStefano Zampini PetscBool iscusp; 74674ae819SStefano Zampini 75674ae819SStefano Zampini PetscFunctionBegin; 768ead10e4SStefano Zampini if (isdir) local_ksp = pcbddc->ksp_D; /* Dirichlet solver */ 778ead10e4SStefano Zampini else local_ksp = pcbddc->ksp_R; /* Neumann solver */ 789566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(local_ksp, &local_mat, &local_pmat)); 799566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(local_pmat, &NullSpace)); 802958b453SStefano Zampini if (!NullSpace) { 81*48a46eb9SPierre Jolivet if (pcbddc->dbg_flag) PetscCall(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer, "Subdomain %04d doesn't have local (near) nullspace: no need for correction in %s solver \n", PetscGlobalRank, isdir ? "Dirichlet" : "Neumann")); 8235529e7bSStefano Zampini PetscFunctionReturn(0); 8335529e7bSStefano Zampini } 849566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)NullSpace, "_PBDDC_Null_dmat", (PetscObject *)&dmat)); 8528b400f6SJacob Faibussowitsch PetscCheck(dmat, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing dense matrix"); 869566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_BDDC_ApproxSetUp[pcbddc->current_level], pc, 0, 0, 0)); 872958b453SStefano Zampini 889566063dSJacob Faibussowitsch PetscCall(PetscNew(&shell_ctx)); 89bd5e1169SStefano Zampini shell_ctx->scale = 1.0; 909566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dmat)); 9192cccca0SStefano Zampini shell_ctx->basis_mat = dmat; 929566063dSJacob Faibussowitsch PetscCall(MatGetSize(dmat, NULL, &basis_size)); 938ead10e4SStefano Zampini shell_ctx->evapply = PC_BDDC_ApproxApply[pcbddc->current_level]; 9492cccca0SStefano Zampini 95b94d7dedSBarry Smith PetscCall(MatIsSymmetric(local_mat, 0.0, &shell_ctx->symm)); 967ebab0bbSStefano Zampini 9792cccca0SStefano Zampini /* explicit construct (Phi^T K Phi)^-1 */ 989566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)local_mat, MATSEQAIJCUSPARSE, &iscusp)); 99*48a46eb9SPierre Jolivet if (iscusp) PetscCall(MatConvert(shell_ctx->basis_mat, MATSEQDENSECUDA, MAT_INPLACE_MATRIX, &shell_ctx->basis_mat)); 1009566063dSJacob Faibussowitsch PetscCall(MatMatMult(local_mat, shell_ctx->basis_mat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Kbasis_mat)); 1019566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(Kbasis_mat, shell_ctx->basis_mat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &shell_ctx->inv_smat)); 1029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Kbasis_mat)); 1039566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(shell_ctx->inv_smat, PETSC_TRUE)); 1049566063dSJacob Faibussowitsch PetscCall(MatFindZeroRows(shell_ctx->inv_smat, &zerorows)); 10592cccca0SStefano Zampini if (zerorows) { /* linearly dependent basis */ 10692cccca0SStefano Zampini const PetscInt *idxs; 10792cccca0SStefano Zampini PetscInt i, nz; 10892cccca0SStefano Zampini 1099566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(zerorows, &nz)); 1109566063dSJacob Faibussowitsch PetscCall(ISGetIndices(zerorows, &idxs)); 111*48a46eb9SPierre Jolivet for (i = 0; i < nz; i++) PetscCall(MatSetValue(shell_ctx->inv_smat, idxs[i], idxs[i], 1.0, INSERT_VALUES)); 1129566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(zerorows, &idxs)); 1139566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(shell_ctx->inv_smat, MAT_FINAL_ASSEMBLY)); 1149566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(shell_ctx->inv_smat, MAT_FINAL_ASSEMBLY)); 11592cccca0SStefano Zampini } 1169566063dSJacob Faibussowitsch PetscCall(MatLUFactor(shell_ctx->inv_smat, NULL, NULL, NULL)); 1179566063dSJacob Faibussowitsch PetscCall(MatSeqDenseInvertFactors_Private(shell_ctx->inv_smat)); 11892cccca0SStefano Zampini if (zerorows) { /* linearly dependent basis */ 11992cccca0SStefano Zampini const PetscInt *idxs; 12092cccca0SStefano Zampini PetscInt i, nz; 12192cccca0SStefano Zampini 1229566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(zerorows, &nz)); 1239566063dSJacob Faibussowitsch PetscCall(ISGetIndices(zerorows, &idxs)); 124*48a46eb9SPierre Jolivet for (i = 0; i < nz; i++) PetscCall(MatSetValue(shell_ctx->inv_smat, idxs[i], idxs[i], 0.0, INSERT_VALUES)); 1259566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(zerorows, &idxs)); 1269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(shell_ctx->inv_smat, MAT_FINAL_ASSEMBLY)); 1279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(shell_ctx->inv_smat, MAT_FINAL_ASSEMBLY)); 12892cccca0SStefano Zampini } 1299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&zerorows)); 130674ae819SStefano Zampini 131674ae819SStefano Zampini /* Create work vectors in shell context */ 1329566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(shell_ctx->inv_smat, &v, NULL)); 1339566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(local_ksp, 1, &shell_ctx->fw, 0, NULL)); 1349566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(v, 3, &shell_ctx->sw)); 1359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 136674ae819SStefano Zampini 13792cccca0SStefano Zampini /* add special pre/post solve to KSP (see [1], eq. 48) */ 1389566063dSJacob Faibussowitsch PetscCall(KSPSetPreSolve(local_ksp, PCBDDCNullSpaceCorrPreSolve, shell_ctx)); 1399566063dSJacob Faibussowitsch PetscCall(KSPSetPostSolve(local_ksp, PCBDDCNullSpaceCorrPostSolve, shell_ctx)); 1409566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)local_ksp), &c)); 1419566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(c, shell_ctx)); 1429566063dSJacob Faibussowitsch PetscCall(PetscContainerSetUserDestroy(c, PCBDDCNullSpaceCorrDestroy)); 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)local_ksp, "_PCBDDC_Null_PrePost_ctx", (PetscObject)c)); 1449566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&c)); 145c7017625SStefano Zampini 146c7017625SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 14792cccca0SStefano Zampini if (needscaling || pcbddc->dbg_flag) { 148674ae819SStefano Zampini KSP check_ksp; 14992cccca0SStefano Zampini PC local_pc; 15092cccca0SStefano Zampini Vec work1, work2; 151bd5e1169SStefano Zampini const char *prefix; 15292cccca0SStefano Zampini PetscReal test_err, lambda_min, lambda_max; 15392cccca0SStefano Zampini PetscInt k, maxit; 154b94d7dedSBarry Smith PetscBool isspd, isset; 155c7017625SStefano Zampini 1569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell_ctx->fw[0], &work1)); 1579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell_ctx->fw[0], &work2)); 1589566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &check_ksp)); 159b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(local_mat, &isset, &isspd)); 160b94d7dedSBarry Smith if (isset && isspd) PetscCall(KSPSetType(check_ksp, KSPCG)); 1619566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)check_ksp, (PetscObject)local_ksp, 0)); 1629566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(local_ksp, &prefix)); 1639566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(check_ksp, prefix)); 1649566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(check_ksp, "approximate_scale_")); 1659566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(check_ksp, PETSC_FALSE)); 1669566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(check_ksp, local_mat, local_pmat)); 1679566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(check_ksp, PETSC_TRUE)); 1689566063dSJacob Faibussowitsch PetscCall(KSPSetPreSolve(check_ksp, PCBDDCNullSpaceCorrPreSolve, shell_ctx)); 1699566063dSJacob Faibussowitsch PetscCall(KSPSetPostSolve(check_ksp, PCBDDCNullSpaceCorrPostSolve, shell_ctx)); 1709566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(check_ksp, PETSC_SMALL, PETSC_SMALL, PETSC_DEFAULT, PETSC_DEFAULT)); 1719566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(check_ksp)); 17282b5ce2aSStefano Zampini /* setup with default maxit, then set maxit to min(10,any_set_from_command_line) (bug in computing eigenvalues when changing the number of iterations */ 1739566063dSJacob Faibussowitsch PetscCall(KSPSetUp(check_ksp)); 1749566063dSJacob Faibussowitsch PetscCall(KSPGetPC(local_ksp, &local_pc)); 1759566063dSJacob Faibussowitsch PetscCall(KSPSetPC(check_ksp, local_pc)); 1769566063dSJacob Faibussowitsch PetscCall(KSPGetTolerances(check_ksp, NULL, NULL, NULL, &maxit)); 1779566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(check_ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, PetscMin(10, maxit))); 1789566063dSJacob Faibussowitsch PetscCall(VecSetRandom(work2, NULL)); 1799566063dSJacob Faibussowitsch PetscCall(MatMult(local_mat, work2, work1)); 1809566063dSJacob Faibussowitsch PetscCall(KSPSolve(check_ksp, work1, work1)); 1819566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(check_ksp, pc, work1)); 1829566063dSJacob Faibussowitsch PetscCall(VecAXPY(work1, -1., work2)); 1839566063dSJacob Faibussowitsch PetscCall(VecNorm(work1, NORM_INFINITY, &test_err)); 1849566063dSJacob Faibussowitsch PetscCall(KSPComputeExtremeSingularValues(check_ksp, &lambda_max, &lambda_min)); 1859566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(check_ksp, &k)); 186c7017625SStefano Zampini if (pcbddc->dbg_flag) { 187c7017625SStefano Zampini if (isdir) { 18863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer, "Subdomain %04d infinity error for Dirichlet adapted solver (no scale) %1.14e (it %" PetscInt_FMT ", eigs %1.6e %1.6e)\n", PetscGlobalRank, (double)test_err, k, (double)lambda_min, (double)lambda_max)); 189c7017625SStefano Zampini } else { 19063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer, "Subdomain %04d infinity error for Neumann adapted solver (no scale) %1.14e (it %" PetscInt_FMT ", eigs %1.6e %1.6e)\n", PetscGlobalRank, (double)test_err, k, (double)lambda_min, (double)lambda_max)); 191c7017625SStefano Zampini } 192c7017625SStefano Zampini } 19392cccca0SStefano Zampini if (needscaling) shell_ctx->scale = 1.0 / lambda_max; 19492cccca0SStefano Zampini 19592cccca0SStefano Zampini if (needscaling && pcbddc->dbg_flag) { /* test for scaling factor */ 19692cccca0SStefano Zampini PC new_pc; 19792cccca0SStefano Zampini 1989566063dSJacob Faibussowitsch PetscCall(VecSetRandom(work2, NULL)); 1999566063dSJacob Faibussowitsch PetscCall(MatMult(local_mat, work2, work1)); 2009566063dSJacob Faibussowitsch PetscCall(PCCreate(PetscObjectComm((PetscObject)check_ksp), &new_pc)); 2019566063dSJacob Faibussowitsch PetscCall(PCSetType(new_pc, PCKSP)); 2029566063dSJacob Faibussowitsch PetscCall(PCSetOperators(new_pc, local_mat, local_pmat)); 2039566063dSJacob Faibussowitsch PetscCall(PCKSPSetKSP(new_pc, local_ksp)); 2049566063dSJacob Faibussowitsch PetscCall(KSPSetPC(check_ksp, new_pc)); 2059566063dSJacob Faibussowitsch PetscCall(PCDestroy(&new_pc)); 2069566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(check_ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, maxit)); 2079566063dSJacob Faibussowitsch PetscCall(KSPSetPreSolve(check_ksp, NULL, NULL)); 2089566063dSJacob Faibussowitsch PetscCall(KSPSetPostSolve(check_ksp, NULL, NULL)); 2099566063dSJacob Faibussowitsch PetscCall(KSPSolve(check_ksp, work1, work1)); 2109566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(check_ksp, pc, work1)); 2119566063dSJacob Faibussowitsch PetscCall(VecAXPY(work1, -1., work2)); 2129566063dSJacob Faibussowitsch PetscCall(VecNorm(work1, NORM_INFINITY, &test_err)); 2139566063dSJacob Faibussowitsch PetscCall(KSPComputeExtremeSingularValues(check_ksp, &lambda_max, &lambda_min)); 2149566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(check_ksp, &k)); 21592cccca0SStefano Zampini if (pcbddc->dbg_flag) { 21692cccca0SStefano Zampini if (isdir) { 21763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer, "Subdomain %04d infinity error for Dirichlet adapted solver (scale %g) %1.14e (it %" PetscInt_FMT ", eigs %1.6e %1.6e)\n", PetscGlobalRank, (double)PetscRealPart(shell_ctx->scale), (double)test_err, k, (double)lambda_min, (double)lambda_max)); 21892cccca0SStefano Zampini } else { 21963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer, "Subdomain %04d infinity error for Neumann adapted solver (scale %g) %1.14e (it %" PetscInt_FMT ", eigs %1.6e %1.6e)\n", PetscGlobalRank, (double)PetscRealPart(shell_ctx->scale), (double)test_err, k, (double)lambda_min, (double)lambda_max)); 22092cccca0SStefano Zampini } 22192cccca0SStefano Zampini } 22292cccca0SStefano Zampini } 2239566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&check_ksp)); 2249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work1)); 2259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work2)); 226c7017625SStefano Zampini } 2279566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_BDDC_ApproxSetUp[pcbddc->current_level], pc, 0, 0, 0)); 228c7017625SStefano Zampini 22992cccca0SStefano Zampini if (pcbddc->dbg_flag) { 230c7017625SStefano Zampini Vec work1, work2, work3; 23192cccca0SStefano Zampini PetscReal test_err; 232674ae819SStefano Zampini 23392cccca0SStefano Zampini /* check nullspace basis is solved exactly */ 2349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell_ctx->fw[0], &work1)); 2359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell_ctx->fw[0], &work2)); 2369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell_ctx->fw[0], &work3)); 2379566063dSJacob Faibussowitsch PetscCall(VecSetRandom(shell_ctx->sw[0], NULL)); 2389566063dSJacob Faibussowitsch PetscCall(MatMult(shell_ctx->basis_mat, shell_ctx->sw[0], work1)); 2399566063dSJacob Faibussowitsch PetscCall(VecCopy(work1, work2)); 2409566063dSJacob Faibussowitsch PetscCall(MatMult(local_mat, work1, work3)); 2419566063dSJacob Faibussowitsch PetscCall(KSPSolve(local_ksp, work3, work1)); 2429566063dSJacob Faibussowitsch PetscCall(VecAXPY(work1, -1., work2)); 2439566063dSJacob Faibussowitsch PetscCall(VecNorm(work1, NORM_INFINITY, &test_err)); 244185763b3SStefano Zampini if (isdir) { 24563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer, "Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n", PetscGlobalRank, (double)test_err)); 246674ae819SStefano Zampini } else { 24763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer, "Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n", PetscGlobalRank, (double)test_err)); 248674ae819SStefano Zampini } 2499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work1)); 2509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work2)); 2519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work3)); 25292cccca0SStefano Zampini } 253674ae819SStefano Zampini PetscFunctionReturn(0); 254674ae819SStefano Zampini } 255