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 */ 6d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBDDCNullSpaceCorrPreSolve(KSP ksp, Vec y, Vec x, void *ctx) 7d71ae5a4SJacob Faibussowitsch { 892cccca0SStefano Zampini NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 992cccca0SStefano Zampini Mat K; 10bd5e1169SStefano Zampini 11bd5e1169SStefano Zampini PetscFunctionBegin; 129566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(corr_ctx->evapply, ksp, 0, 0, 0)); 139566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(corr_ctx->basis_mat, y, corr_ctx->sw[0])); 147ebab0bbSStefano Zampini if (corr_ctx->symm) { 159566063dSJacob Faibussowitsch PetscCall(MatMult(corr_ctx->inv_smat, corr_ctx->sw[0], corr_ctx->sw[1])); 167ebab0bbSStefano Zampini } else { 179566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(corr_ctx->inv_smat, corr_ctx->sw[0], corr_ctx->sw[1])); 187ebab0bbSStefano Zampini } 199566063dSJacob Faibussowitsch PetscCall(VecScale(corr_ctx->sw[1], -1.0)); 209566063dSJacob Faibussowitsch PetscCall(MatMult(corr_ctx->basis_mat, corr_ctx->sw[1], corr_ctx->fw[0])); 219566063dSJacob Faibussowitsch PetscCall(VecScale(corr_ctx->sw[1], -1.0)); 229566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(ksp, &K, NULL)); 239566063dSJacob Faibussowitsch PetscCall(MatMultAdd(K, corr_ctx->fw[0], y, y)); 249566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(corr_ctx->evapply, ksp, 0, 0, 0)); 253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26bd5e1169SStefano Zampini } 27bd5e1169SStefano Zampini 2892cccca0SStefano Zampini /* E^t + small */ 29d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBDDCNullSpaceCorrPostSolve(KSP ksp, Vec y, Vec x, void *ctx) 30d71ae5a4SJacob Faibussowitsch { 3192cccca0SStefano Zampini NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 3292cccca0SStefano Zampini Mat K; 33674ae819SStefano Zampini 34674ae819SStefano Zampini PetscFunctionBegin; 359566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(corr_ctx->evapply, ksp, 0, 0, 0)); 369566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(ksp, &K, NULL)); 377ebab0bbSStefano Zampini if (corr_ctx->symm) { 389566063dSJacob Faibussowitsch PetscCall(MatMult(K, x, corr_ctx->fw[0])); 397ebab0bbSStefano Zampini } else { 409566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(K, x, corr_ctx->fw[0])); 417ebab0bbSStefano Zampini } 429566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(corr_ctx->basis_mat, corr_ctx->fw[0], corr_ctx->sw[0])); 439566063dSJacob Faibussowitsch PetscCall(VecScale(corr_ctx->sw[0], -1.0)); 449566063dSJacob Faibussowitsch PetscCall(MatMult(corr_ctx->inv_smat, corr_ctx->sw[0], corr_ctx->sw[2])); 459566063dSJacob Faibussowitsch PetscCall(MatMultAdd(corr_ctx->basis_mat, corr_ctx->sw[2], x, corr_ctx->fw[0])); 469566063dSJacob Faibussowitsch PetscCall(VecScale(corr_ctx->fw[0], corr_ctx->scale)); 4792cccca0SStefano Zampini /* Sum contributions from approximate solver and projected system */ 489566063dSJacob Faibussowitsch PetscCall(MatMultAdd(corr_ctx->basis_mat, corr_ctx->sw[1], corr_ctx->fw[0], x)); 499566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(corr_ctx->evapply, ksp, 0, 0, 0)); 503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51674ae819SStefano Zampini } 52674ae819SStefano Zampini 53d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBDDCNullSpaceCorrDestroy(void *ctx) 54d71ae5a4SJacob Faibussowitsch { 5592cccca0SStefano Zampini NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 56674ae819SStefano Zampini 57674ae819SStefano Zampini PetscFunctionBegin; 589566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(3, &corr_ctx->sw)); 599566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(1, &corr_ctx->fw)); 609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&corr_ctx->basis_mat)); 619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&corr_ctx->inv_smat)); 629566063dSJacob Faibussowitsch PetscCall(PetscFree(corr_ctx)); 633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64674ae819SStefano Zampini } 65674ae819SStefano Zampini 66d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling) 67d71ae5a4SJacob Faibussowitsch { 68674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 69c7017625SStefano Zampini MatNullSpace NullSpace = NULL; 70185763b3SStefano Zampini KSP local_ksp; 71674ae819SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 7292cccca0SStefano Zampini Mat local_mat, local_pmat, dmat, Kbasis_mat; 7392cccca0SStefano Zampini Vec v; 7492cccca0SStefano Zampini PetscContainer c; 7592cccca0SStefano Zampini PetscInt basis_size; 7692cccca0SStefano Zampini IS zerorows; 7780fdaca0SStefano Zampini PetscBool iscusp; 78674ae819SStefano Zampini 79674ae819SStefano Zampini PetscFunctionBegin; 808ead10e4SStefano Zampini if (isdir) local_ksp = pcbddc->ksp_D; /* Dirichlet solver */ 818ead10e4SStefano Zampini else local_ksp = pcbddc->ksp_R; /* Neumann solver */ 829566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(local_ksp, &local_mat, &local_pmat)); 839566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(local_pmat, &NullSpace)); 842958b453SStefano Zampini if (!NullSpace) { 8548a46eb9SPierre 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")); 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8735529e7bSStefano Zampini } 889566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)NullSpace, "_PBDDC_Null_dmat", (PetscObject *)&dmat)); 8928b400f6SJacob Faibussowitsch PetscCheck(dmat, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing dense matrix"); 909566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_BDDC_ApproxSetUp[pcbddc->current_level], pc, 0, 0, 0)); 912958b453SStefano Zampini 929566063dSJacob Faibussowitsch PetscCall(PetscNew(&shell_ctx)); 93bd5e1169SStefano Zampini shell_ctx->scale = 1.0; 949566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dmat)); 9592cccca0SStefano Zampini shell_ctx->basis_mat = dmat; 969566063dSJacob Faibussowitsch PetscCall(MatGetSize(dmat, NULL, &basis_size)); 978ead10e4SStefano Zampini shell_ctx->evapply = PC_BDDC_ApproxApply[pcbddc->current_level]; 9892cccca0SStefano Zampini 99b94d7dedSBarry Smith PetscCall(MatIsSymmetric(local_mat, 0.0, &shell_ctx->symm)); 1007ebab0bbSStefano Zampini 10192cccca0SStefano Zampini /* explicit construct (Phi^T K Phi)^-1 */ 1029566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)local_mat, MATSEQAIJCUSPARSE, &iscusp)); 10348a46eb9SPierre Jolivet if (iscusp) PetscCall(MatConvert(shell_ctx->basis_mat, MATSEQDENSECUDA, MAT_INPLACE_MATRIX, &shell_ctx->basis_mat)); 1049566063dSJacob Faibussowitsch PetscCall(MatMatMult(local_mat, shell_ctx->basis_mat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Kbasis_mat)); 1059566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(Kbasis_mat, shell_ctx->basis_mat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &shell_ctx->inv_smat)); 1069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Kbasis_mat)); 1079566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(shell_ctx->inv_smat, PETSC_TRUE)); 1089566063dSJacob Faibussowitsch PetscCall(MatFindZeroRows(shell_ctx->inv_smat, &zerorows)); 10992cccca0SStefano Zampini if (zerorows) { /* linearly dependent basis */ 11092cccca0SStefano Zampini const PetscInt *idxs; 11192cccca0SStefano Zampini PetscInt i, nz; 11292cccca0SStefano Zampini 1139566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(zerorows, &nz)); 1149566063dSJacob Faibussowitsch PetscCall(ISGetIndices(zerorows, &idxs)); 11548a46eb9SPierre Jolivet for (i = 0; i < nz; i++) PetscCall(MatSetValue(shell_ctx->inv_smat, idxs[i], idxs[i], 1.0, INSERT_VALUES)); 1169566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(zerorows, &idxs)); 1179566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(shell_ctx->inv_smat, MAT_FINAL_ASSEMBLY)); 1189566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(shell_ctx->inv_smat, MAT_FINAL_ASSEMBLY)); 11992cccca0SStefano Zampini } 1209566063dSJacob Faibussowitsch PetscCall(MatLUFactor(shell_ctx->inv_smat, NULL, NULL, NULL)); 1219566063dSJacob Faibussowitsch PetscCall(MatSeqDenseInvertFactors_Private(shell_ctx->inv_smat)); 12292cccca0SStefano Zampini if (zerorows) { /* linearly dependent basis */ 12392cccca0SStefano Zampini const PetscInt *idxs; 12492cccca0SStefano Zampini PetscInt i, nz; 12592cccca0SStefano Zampini 1269566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(zerorows, &nz)); 1279566063dSJacob Faibussowitsch PetscCall(ISGetIndices(zerorows, &idxs)); 12848a46eb9SPierre Jolivet for (i = 0; i < nz; i++) PetscCall(MatSetValue(shell_ctx->inv_smat, idxs[i], idxs[i], 0.0, INSERT_VALUES)); 1299566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(zerorows, &idxs)); 1309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(shell_ctx->inv_smat, MAT_FINAL_ASSEMBLY)); 1319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(shell_ctx->inv_smat, MAT_FINAL_ASSEMBLY)); 13292cccca0SStefano Zampini } 1339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&zerorows)); 134674ae819SStefano Zampini 135674ae819SStefano Zampini /* Create work vectors in shell context */ 1369566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(shell_ctx->inv_smat, &v, NULL)); 1379566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(local_ksp, 1, &shell_ctx->fw, 0, NULL)); 1389566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(v, 3, &shell_ctx->sw)); 1399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 140674ae819SStefano Zampini 14192cccca0SStefano Zampini /* add special pre/post solve to KSP (see [1], eq. 48) */ 1429566063dSJacob Faibussowitsch PetscCall(KSPSetPreSolve(local_ksp, PCBDDCNullSpaceCorrPreSolve, shell_ctx)); 1439566063dSJacob Faibussowitsch PetscCall(KSPSetPostSolve(local_ksp, PCBDDCNullSpaceCorrPostSolve, shell_ctx)); 1449566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)local_ksp), &c)); 1459566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(c, shell_ctx)); 1469566063dSJacob Faibussowitsch PetscCall(PetscContainerSetUserDestroy(c, PCBDDCNullSpaceCorrDestroy)); 1479566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)local_ksp, "_PCBDDC_Null_PrePost_ctx", (PetscObject)c)); 1489566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&c)); 149c7017625SStefano Zampini 150c7017625SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 15192cccca0SStefano Zampini if (needscaling || pcbddc->dbg_flag) { 152674ae819SStefano Zampini KSP check_ksp; 15392cccca0SStefano Zampini PC local_pc; 15492cccca0SStefano Zampini Vec work1, work2; 155bd5e1169SStefano Zampini const char *prefix; 15692cccca0SStefano Zampini PetscReal test_err, lambda_min, lambda_max; 15792cccca0SStefano Zampini PetscInt k, maxit; 158b94d7dedSBarry Smith PetscBool isspd, isset; 159c7017625SStefano Zampini 1609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell_ctx->fw[0], &work1)); 1619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell_ctx->fw[0], &work2)); 1629566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &check_ksp)); 163*3821be0aSBarry Smith PetscCall(KSPSetNestLevel(check_ksp, pc->kspnestlevel)); 164b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(local_mat, &isset, &isspd)); 165b94d7dedSBarry Smith if (isset && isspd) PetscCall(KSPSetType(check_ksp, KSPCG)); 1669566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)check_ksp, (PetscObject)local_ksp, 0)); 1679566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(local_ksp, &prefix)); 1689566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(check_ksp, prefix)); 1699566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(check_ksp, "approximate_scale_")); 1709566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(check_ksp, PETSC_FALSE)); 1719566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(check_ksp, local_mat, local_pmat)); 1729566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(check_ksp, PETSC_TRUE)); 1739566063dSJacob Faibussowitsch PetscCall(KSPSetPreSolve(check_ksp, PCBDDCNullSpaceCorrPreSolve, shell_ctx)); 1749566063dSJacob Faibussowitsch PetscCall(KSPSetPostSolve(check_ksp, PCBDDCNullSpaceCorrPostSolve, shell_ctx)); 1759566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(check_ksp, PETSC_SMALL, PETSC_SMALL, PETSC_DEFAULT, PETSC_DEFAULT)); 1769566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(check_ksp)); 17782b5ce2aSStefano 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 */ 1789566063dSJacob Faibussowitsch PetscCall(KSPSetUp(check_ksp)); 1799566063dSJacob Faibussowitsch PetscCall(KSPGetPC(local_ksp, &local_pc)); 1809566063dSJacob Faibussowitsch PetscCall(KSPSetPC(check_ksp, local_pc)); 1819566063dSJacob Faibussowitsch PetscCall(KSPGetTolerances(check_ksp, NULL, NULL, NULL, &maxit)); 1829566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(check_ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, PetscMin(10, maxit))); 1839566063dSJacob Faibussowitsch PetscCall(VecSetRandom(work2, NULL)); 1849566063dSJacob Faibussowitsch PetscCall(MatMult(local_mat, work2, work1)); 1859566063dSJacob Faibussowitsch PetscCall(KSPSolve(check_ksp, work1, work1)); 1869566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(check_ksp, pc, work1)); 1879566063dSJacob Faibussowitsch PetscCall(VecAXPY(work1, -1., work2)); 1889566063dSJacob Faibussowitsch PetscCall(VecNorm(work1, NORM_INFINITY, &test_err)); 1899566063dSJacob Faibussowitsch PetscCall(KSPComputeExtremeSingularValues(check_ksp, &lambda_max, &lambda_min)); 1909566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(check_ksp, &k)); 191c7017625SStefano Zampini if (pcbddc->dbg_flag) { 192c7017625SStefano Zampini if (isdir) { 19363a3b9bcSJacob 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)); 194c7017625SStefano Zampini } else { 19563a3b9bcSJacob 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)); 196c7017625SStefano Zampini } 197c7017625SStefano Zampini } 19892cccca0SStefano Zampini if (needscaling) shell_ctx->scale = 1.0 / lambda_max; 19992cccca0SStefano Zampini 20092cccca0SStefano Zampini if (needscaling && pcbddc->dbg_flag) { /* test for scaling factor */ 20192cccca0SStefano Zampini PC new_pc; 20292cccca0SStefano Zampini 2039566063dSJacob Faibussowitsch PetscCall(VecSetRandom(work2, NULL)); 2049566063dSJacob Faibussowitsch PetscCall(MatMult(local_mat, work2, work1)); 2059566063dSJacob Faibussowitsch PetscCall(PCCreate(PetscObjectComm((PetscObject)check_ksp), &new_pc)); 2069566063dSJacob Faibussowitsch PetscCall(PCSetType(new_pc, PCKSP)); 2079566063dSJacob Faibussowitsch PetscCall(PCSetOperators(new_pc, local_mat, local_pmat)); 2089566063dSJacob Faibussowitsch PetscCall(PCKSPSetKSP(new_pc, local_ksp)); 2099566063dSJacob Faibussowitsch PetscCall(KSPSetPC(check_ksp, new_pc)); 2109566063dSJacob Faibussowitsch PetscCall(PCDestroy(&new_pc)); 2119566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(check_ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, maxit)); 2129566063dSJacob Faibussowitsch PetscCall(KSPSetPreSolve(check_ksp, NULL, NULL)); 2139566063dSJacob Faibussowitsch PetscCall(KSPSetPostSolve(check_ksp, NULL, NULL)); 2149566063dSJacob Faibussowitsch PetscCall(KSPSolve(check_ksp, work1, work1)); 2159566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(check_ksp, pc, work1)); 2169566063dSJacob Faibussowitsch PetscCall(VecAXPY(work1, -1., work2)); 2179566063dSJacob Faibussowitsch PetscCall(VecNorm(work1, NORM_INFINITY, &test_err)); 2189566063dSJacob Faibussowitsch PetscCall(KSPComputeExtremeSingularValues(check_ksp, &lambda_max, &lambda_min)); 2199566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(check_ksp, &k)); 22092cccca0SStefano Zampini if (pcbddc->dbg_flag) { 22192cccca0SStefano Zampini if (isdir) { 22263a3b9bcSJacob 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)); 22392cccca0SStefano Zampini } else { 22463a3b9bcSJacob 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)); 22592cccca0SStefano Zampini } 22692cccca0SStefano Zampini } 22792cccca0SStefano Zampini } 2289566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&check_ksp)); 2299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work1)); 2309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work2)); 231c7017625SStefano Zampini } 2329566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_BDDC_ApproxSetUp[pcbddc->current_level], pc, 0, 0, 0)); 233c7017625SStefano Zampini 23492cccca0SStefano Zampini if (pcbddc->dbg_flag) { 235c7017625SStefano Zampini Vec work1, work2, work3; 23692cccca0SStefano Zampini PetscReal test_err; 237674ae819SStefano Zampini 23892cccca0SStefano Zampini /* check nullspace basis is solved exactly */ 2399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell_ctx->fw[0], &work1)); 2409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell_ctx->fw[0], &work2)); 2419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell_ctx->fw[0], &work3)); 2429566063dSJacob Faibussowitsch PetscCall(VecSetRandom(shell_ctx->sw[0], NULL)); 2439566063dSJacob Faibussowitsch PetscCall(MatMult(shell_ctx->basis_mat, shell_ctx->sw[0], work1)); 2449566063dSJacob Faibussowitsch PetscCall(VecCopy(work1, work2)); 2459566063dSJacob Faibussowitsch PetscCall(MatMult(local_mat, work1, work3)); 2469566063dSJacob Faibussowitsch PetscCall(KSPSolve(local_ksp, work3, work1)); 2479566063dSJacob Faibussowitsch PetscCall(VecAXPY(work1, -1., work2)); 2489566063dSJacob Faibussowitsch PetscCall(VecNorm(work1, NORM_INFINITY, &test_err)); 249185763b3SStefano Zampini if (isdir) { 25063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer, "Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n", PetscGlobalRank, (double)test_err)); 251674ae819SStefano Zampini } else { 25263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer, "Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n", PetscGlobalRank, (double)test_err)); 253674ae819SStefano Zampini } 2549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work1)); 2559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work2)); 2569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work3)); 25792cccca0SStefano Zampini } 2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 259674ae819SStefano Zampini } 260