xref: /petsc/src/ksp/pc/impls/bddc/bddcnullspace.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
392cccca0SStefano Zampini #include <../src/mat/impls/dense/seq/dense.h>
4674ae819SStefano Zampini 
592cccca0SStefano Zampini /* E + small_solve */
692cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrPreSolve(KSP ksp,Vec y,Vec x, void* ctx)
7bd5e1169SStefano Zampini {
892cccca0SStefano Zampini   NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx;
992cccca0SStefano Zampini   Mat                     K;
10bd5e1169SStefano Zampini 
11bd5e1169SStefano Zampini   PetscFunctionBegin;
125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(corr_ctx->evapply,ksp,0,0,0));
135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(corr_ctx->basis_mat,y,corr_ctx->sw[0]));
147ebab0bbSStefano Zampini   if (corr_ctx->symm) {
155f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[1]));
167ebab0bbSStefano Zampini   } else {
175f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTranspose(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[1]));
187ebab0bbSStefano Zampini   }
195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(corr_ctx->sw[1],-1.0));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(corr_ctx->basis_mat,corr_ctx->sw[1],corr_ctx->fw[0]));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(corr_ctx->sw[1],-1.0));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetOperators(ksp,&K,NULL));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(K,corr_ctx->fw[0],y,y));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(corr_ctx->evapply,ksp,0,0,0));
25bd5e1169SStefano Zampini   PetscFunctionReturn(0);
26bd5e1169SStefano Zampini }
27bd5e1169SStefano Zampini 
2892cccca0SStefano Zampini /* E^t + small */
2992cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrPostSolve(KSP ksp,Vec y,Vec x, void* ctx)
30674ae819SStefano Zampini {
3192cccca0SStefano Zampini   NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx;
3292cccca0SStefano Zampini   Mat                     K;
33674ae819SStefano Zampini 
34674ae819SStefano Zampini   PetscFunctionBegin;
355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(corr_ctx->evapply,ksp,0,0,0));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetOperators(ksp,&K,NULL));
377ebab0bbSStefano Zampini   if (corr_ctx->symm) {
385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(K,x,corr_ctx->fw[0]));
397ebab0bbSStefano Zampini   } else {
405f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTranspose(K,x,corr_ctx->fw[0]));
417ebab0bbSStefano Zampini   }
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(corr_ctx->basis_mat,corr_ctx->fw[0],corr_ctx->sw[0]));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(corr_ctx->sw[0],-1.0));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[2]));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(corr_ctx->basis_mat,corr_ctx->sw[2],x,corr_ctx->fw[0]));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(corr_ctx->fw[0],corr_ctx->scale));
4792cccca0SStefano Zampini   /* Sum contributions from approximate solver and projected system */
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(corr_ctx->basis_mat,corr_ctx->sw[1],corr_ctx->fw[0],x));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(corr_ctx->evapply,ksp,0,0,0));
50674ae819SStefano Zampini   PetscFunctionReturn(0);
51674ae819SStefano Zampini }
52674ae819SStefano Zampini 
5392cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrDestroy(void * ctx)
54674ae819SStefano Zampini {
5592cccca0SStefano Zampini   NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx;
56674ae819SStefano Zampini 
57674ae819SStefano Zampini   PetscFunctionBegin;
585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(3,&corr_ctx->sw));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(1,&corr_ctx->fw));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&corr_ctx->basis_mat));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&corr_ctx->inv_smat));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(corr_ctx));
63674ae819SStefano Zampini   PetscFunctionReturn(0);
64674ae819SStefano Zampini }
65674ae819SStefano Zampini 
66c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling)
67674ae819SStefano Zampini {
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 */
825f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetOperators(local_ksp,&local_mat,&local_pmat));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetNearNullSpace(local_pmat,&NullSpace));
842958b453SStefano Zampini   if (!NullSpace) {
8535529e7bSStefano Zampini     if (pcbddc->dbg_flag) {
865f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d doesn't have local (near) nullspace: no need for correction in %s solver \n",PetscGlobalRank,isdir ? "Dirichlet" : "Neumann"));
8735529e7bSStefano Zampini     }
8835529e7bSStefano Zampini     PetscFunctionReturn(0);
8935529e7bSStefano Zampini   }
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectQuery((PetscObject)NullSpace,"_PBDDC_Null_dmat",(PetscObject*)&dmat));
91*28b400f6SJacob Faibussowitsch   PetscCheck(dmat,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing dense matrix");
925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(PC_BDDC_ApproxSetUp[pcbddc->current_level],pc,0,0,0));
932958b453SStefano Zampini 
945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&shell_ctx));
95bd5e1169SStefano Zampini   shell_ctx->scale = 1.0;
965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)dmat));
9792cccca0SStefano Zampini   shell_ctx->basis_mat = dmat;
985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(dmat,NULL,&basis_size));
998ead10e4SStefano Zampini   shell_ctx->evapply = PC_BDDC_ApproxApply[pcbddc->current_level];
10092cccca0SStefano Zampini 
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOption(local_mat,MAT_SYMMETRIC,&shell_ctx->symm));
1027ebab0bbSStefano Zampini 
10392cccca0SStefano Zampini   /* explicit construct (Phi^T K Phi)^-1 */
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)local_mat,MATSEQAIJCUSPARSE,&iscusp));
10580fdaca0SStefano Zampini   if (iscusp) {
1065f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(shell_ctx->basis_mat,MATSEQDENSECUDA,MAT_INPLACE_MATRIX,&shell_ctx->basis_mat));
10780fdaca0SStefano Zampini   }
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(local_mat,shell_ctx->basis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Kbasis_mat));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransposeMatMult(Kbasis_mat,shell_ctx->basis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->inv_smat));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Kbasis_mat));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(MatBindToCPU(shell_ctx->inv_smat,PETSC_TRUE));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFindZeroRows(shell_ctx->inv_smat,&zerorows));
11392cccca0SStefano Zampini   if (zerorows) { /* linearly dependent basis */
11492cccca0SStefano Zampini     const PetscInt *idxs;
11592cccca0SStefano Zampini     PetscInt       i,nz;
11692cccca0SStefano Zampini 
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(zerorows,&nz));
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(zerorows,&idxs));
11992cccca0SStefano Zampini     for (i=0;i<nz;i++) {
1205f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValue(shell_ctx->inv_smat,idxs[i],idxs[i],1.0,INSERT_VALUES));
12192cccca0SStefano Zampini     }
1225f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(zerorows,&idxs));
1235f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY));
1245f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY));
12592cccca0SStefano Zampini   }
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactor(shell_ctx->inv_smat,NULL,NULL,NULL));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqDenseInvertFactors_Private(shell_ctx->inv_smat));
12892cccca0SStefano Zampini   if (zerorows) { /* linearly dependent basis */
12992cccca0SStefano Zampini     const PetscInt *idxs;
13092cccca0SStefano Zampini     PetscInt       i,nz;
13192cccca0SStefano Zampini 
1325f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(zerorows,&nz));
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(zerorows,&idxs));
13492cccca0SStefano Zampini     for (i=0;i<nz;i++) {
1355f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValue(shell_ctx->inv_smat,idxs[i],idxs[i],0.0,INSERT_VALUES));
13692cccca0SStefano Zampini     }
1375f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(zerorows,&idxs));
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY));
1395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY));
14092cccca0SStefano Zampini   }
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&zerorows));
142674ae819SStefano Zampini 
143674ae819SStefano Zampini   /* Create work vectors in shell context */
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(shell_ctx->inv_smat,&v,NULL));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPCreateVecs(local_ksp,1,&shell_ctx->fw,0,NULL));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(v,3,&shell_ctx->sw));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
148674ae819SStefano Zampini 
14992cccca0SStefano Zampini   /* add special pre/post solve to KSP (see [1], eq. 48) */
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetPreSolve(local_ksp,PCBDDCNullSpaceCorrPreSolve,shell_ctx));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetPostSolve(local_ksp,PCBDDCNullSpaceCorrPostSolve,shell_ctx));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscContainerCreate(PetscObjectComm((PetscObject)local_ksp),&c));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscContainerSetPointer(c,shell_ctx));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscContainerSetUserDestroy(c,PCBDDCNullSpaceCorrDestroy));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectCompose((PetscObject)local_ksp,"_PCBDDC_Null_PrePost_ctx",(PetscObject)c));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscContainerDestroy(&c));
157c7017625SStefano Zampini 
158c7017625SStefano Zampini   /* Create ksp object suitable for extreme eigenvalues' estimation */
15992cccca0SStefano Zampini   if (needscaling || pcbddc->dbg_flag) {
160674ae819SStefano Zampini     KSP         check_ksp;
16192cccca0SStefano Zampini     PC          local_pc;
16292cccca0SStefano Zampini     Vec         work1,work2;
163bd5e1169SStefano Zampini     const char* prefix;
16492cccca0SStefano Zampini     PetscReal   test_err,lambda_min,lambda_max;
16592cccca0SStefano Zampini     PetscInt    k,maxit;
166c7017625SStefano Zampini 
1675f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(shell_ctx->fw[0],&work1));
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(shell_ctx->fw[0],&work2));
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPCreate(PETSC_COMM_SELF,&check_ksp));
17092cccca0SStefano Zampini     if (local_mat->spd) {
1715f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSetType(check_ksp,KSPCG));
17292cccca0SStefano Zampini     }
1735f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectIncrementTabLevel((PetscObject)check_ksp,(PetscObject)local_ksp,0));
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetOptionsPrefix(local_ksp,&prefix));
1755f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetOptionsPrefix(check_ksp,prefix));
1765f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPAppendOptionsPrefix(check_ksp,"approximate_scale_"));
1775f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetErrorIfNotConverged(check_ksp,PETSC_FALSE));
1785f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetOperators(check_ksp,local_mat,local_pmat));
1795f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetComputeSingularValues(check_ksp,PETSC_TRUE));
1805f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetPreSolve(check_ksp,PCBDDCNullSpaceCorrPreSolve,shell_ctx));
1815f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetPostSolve(check_ksp,PCBDDCNullSpaceCorrPostSolve,shell_ctx));
1825f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetTolerances(check_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT));
1835f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetFromOptions(check_ksp));
18492cccca0SStefano Zampini     /* setup with default maxit, then set maxit to min(10,any_set_from_command_line) (bug in computing eigenvalues when chaning the number of iterations */
1855f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetUp(check_ksp));
1865f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(local_ksp,&local_pc));
1875f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetPC(check_ksp,local_pc));
1885f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetTolerances(check_ksp,NULL,NULL,NULL,&maxit));
1895f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PetscMin(10,maxit)));
1905f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(work2,NULL));
1915f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(local_mat,work2,work1));
1925f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSolve(check_ksp,work1,work1));
1935f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPCheckSolve(check_ksp,pc,work1));
1945f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(work1,-1.,work2));
1955f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(work1,NORM_INFINITY,&test_err));
1965f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min));
1975f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetIterationNumber(check_ksp,&k));
198c7017625SStefano Zampini     if (pcbddc->dbg_flag) {
199c7017625SStefano Zampini       if (isdir) {
2005f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet adapted solver (no scale) %1.14e (it %D, eigs %1.6e %1.6e)\n",PetscGlobalRank,test_err,k,lambda_min,lambda_max));
201c7017625SStefano Zampini       } else {
2025f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann adapted solver (no scale) %1.14e (it %D, eigs %1.6e %1.6e)\n",PetscGlobalRank,test_err,k,lambda_min,lambda_max));
203c7017625SStefano Zampini       }
204c7017625SStefano Zampini     }
20592cccca0SStefano Zampini     if (needscaling) shell_ctx->scale = 1.0/lambda_max;
20692cccca0SStefano Zampini 
20792cccca0SStefano Zampini     if (needscaling && pcbddc->dbg_flag) { /* test for scaling factor */
20892cccca0SStefano Zampini       PC new_pc;
20992cccca0SStefano Zampini 
2105f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(work2,NULL));
2115f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(local_mat,work2,work1));
2125f80ce2aSJacob Faibussowitsch       CHKERRQ(PCCreate(PetscObjectComm((PetscObject)check_ksp),&new_pc));
2135f80ce2aSJacob Faibussowitsch       CHKERRQ(PCSetType(new_pc,PCKSP));
2145f80ce2aSJacob Faibussowitsch       CHKERRQ(PCSetOperators(new_pc,local_mat,local_pmat));
2155f80ce2aSJacob Faibussowitsch       CHKERRQ(PCKSPSetKSP(new_pc,local_ksp));
2165f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSetPC(check_ksp,new_pc));
2175f80ce2aSJacob Faibussowitsch       CHKERRQ(PCDestroy(&new_pc));
2185f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit));
2195f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSetPreSolve(check_ksp,NULL,NULL));
2205f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSetPostSolve(check_ksp,NULL,NULL));
2215f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSolve(check_ksp,work1,work1));
2225f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPCheckSolve(check_ksp,pc,work1));
2235f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(work1,-1.,work2));
2245f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(work1,NORM_INFINITY,&test_err));
2255f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min));
2265f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPGetIterationNumber(check_ksp,&k));
22792cccca0SStefano Zampini       if (pcbddc->dbg_flag) {
22892cccca0SStefano Zampini         if (isdir) {
2295f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet adapted solver (scale %g) %1.14e (it %D, eigs %1.6e %1.6e)\n",PetscGlobalRank,(double)PetscRealPart(shell_ctx->scale),test_err,k,lambda_min,lambda_max));
23092cccca0SStefano Zampini         } else {
2315f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann adapted solver (scale %g) %1.14e (it %D, eigs %1.6e %1.6e)\n",PetscGlobalRank,(double)PetscRealPart(shell_ctx->scale),test_err,k,lambda_min,lambda_max));
23292cccca0SStefano Zampini         }
23392cccca0SStefano Zampini       }
23492cccca0SStefano Zampini     }
2355f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPDestroy(&check_ksp));
2365f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&work1));
2375f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&work2));
238c7017625SStefano Zampini   }
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(PC_BDDC_ApproxSetUp[pcbddc->current_level],pc,0,0,0));
240c7017625SStefano Zampini 
24192cccca0SStefano Zampini   if (pcbddc->dbg_flag) {
242c7017625SStefano Zampini     Vec       work1,work2,work3;
24392cccca0SStefano Zampini     PetscReal test_err;
244674ae819SStefano Zampini 
24592cccca0SStefano Zampini     /* check nullspace basis is solved exactly */
2465f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(shell_ctx->fw[0],&work1));
2475f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(shell_ctx->fw[0],&work2));
2485f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(shell_ctx->fw[0],&work3));
2495f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(shell_ctx->sw[0],NULL));
2505f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(shell_ctx->basis_mat,shell_ctx->sw[0],work1));
2515f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(work1,work2));
2525f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(local_mat,work1,work3));
2535f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSolve(local_ksp,work3,work1));
2545f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(work1,-1.,work2));
2555f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(work1,NORM_INFINITY,&test_err));
256185763b3SStefano Zampini     if (isdir) {
2575f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err));
258674ae819SStefano Zampini     } else {
2595f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err));
260674ae819SStefano Zampini     }
2615f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&work1));
2625f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&work2));
2635f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&work3));
26492cccca0SStefano Zampini   }
265674ae819SStefano Zampini   PetscFunctionReturn(0);
266674ae819SStefano Zampini }
267