xref: /petsc/src/ksp/pc/impls/bddc/bddcnullspace.c (revision 7ebab0bbd39729a4f6e2306e6cac2200f5e55332)
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   PetscErrorCode          ierr;
11bd5e1169SStefano Zampini 
12bd5e1169SStefano Zampini   PetscFunctionBegin;
138ead10e4SStefano Zampini   ierr = PetscLogEventBegin(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr);
1492cccca0SStefano Zampini   ierr = MatMultTranspose(corr_ctx->basis_mat,y,corr_ctx->sw[0]);CHKERRQ(ierr);
15*7ebab0bbSStefano Zampini   if (corr_ctx->symm) {
16*7ebab0bbSStefano Zampini     ierr = MatMult(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[1]);CHKERRQ(ierr);
17*7ebab0bbSStefano Zampini   } else {
1892cccca0SStefano Zampini     ierr = MatMultTranspose(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[1]);CHKERRQ(ierr);
19*7ebab0bbSStefano Zampini   }
2092cccca0SStefano Zampini   ierr = VecScale(corr_ctx->sw[1],-1.0);CHKERRQ(ierr);
2192cccca0SStefano Zampini   ierr = MatMult(corr_ctx->basis_mat,corr_ctx->sw[1],corr_ctx->fw[0]);CHKERRQ(ierr);
2292cccca0SStefano Zampini   ierr = VecScale(corr_ctx->sw[1],-1.0);CHKERRQ(ierr);
2392cccca0SStefano Zampini   ierr = KSPGetOperators(ksp,&K,NULL);CHKERRQ(ierr);
2492cccca0SStefano Zampini   ierr = MatMultAdd(K,corr_ctx->fw[0],y,y);CHKERRQ(ierr);
258ead10e4SStefano Zampini   ierr = PetscLogEventEnd(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr);
26bd5e1169SStefano Zampini   PetscFunctionReturn(0);
27bd5e1169SStefano Zampini }
28bd5e1169SStefano Zampini 
2992cccca0SStefano Zampini /* E^t + small */
3092cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrPostSolve(KSP ksp,Vec y,Vec x, void* ctx)
31674ae819SStefano Zampini {
3292cccca0SStefano Zampini   NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx;
33674ae819SStefano Zampini   PetscErrorCode          ierr;
3492cccca0SStefano Zampini   Mat                     K;
35674ae819SStefano Zampini 
36674ae819SStefano Zampini   PetscFunctionBegin;
378ead10e4SStefano Zampini   ierr = PetscLogEventBegin(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr);
3892cccca0SStefano Zampini   ierr = KSPGetOperators(ksp,&K,NULL);CHKERRQ(ierr);
39*7ebab0bbSStefano Zampini   if (corr_ctx->symm) {
40*7ebab0bbSStefano Zampini     ierr = MatMult(K,x,corr_ctx->fw[0]);CHKERRQ(ierr);
41*7ebab0bbSStefano Zampini   } else {
4292cccca0SStefano Zampini     ierr = MatMultTranspose(K,x,corr_ctx->fw[0]);CHKERRQ(ierr);
43*7ebab0bbSStefano Zampini   }
4492cccca0SStefano Zampini   ierr = MatMultTranspose(corr_ctx->basis_mat,corr_ctx->fw[0],corr_ctx->sw[0]);CHKERRQ(ierr);
4592cccca0SStefano Zampini   ierr = VecScale(corr_ctx->sw[0],-1.0);CHKERRQ(ierr);
4692cccca0SStefano Zampini   ierr = MatMult(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[2]);CHKERRQ(ierr);
4792cccca0SStefano Zampini   ierr = MatMultAdd(corr_ctx->basis_mat,corr_ctx->sw[2],x,corr_ctx->fw[0]);CHKERRQ(ierr);
4892cccca0SStefano Zampini   ierr = VecScale(corr_ctx->fw[0],corr_ctx->scale);CHKERRQ(ierr);
4992cccca0SStefano Zampini   /* Sum contributions from approximate solver and projected system */
5092cccca0SStefano Zampini   ierr = MatMultAdd(corr_ctx->basis_mat,corr_ctx->sw[1],corr_ctx->fw[0],x);CHKERRQ(ierr);
518ead10e4SStefano Zampini   ierr = PetscLogEventEnd(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr);
52674ae819SStefano Zampini   PetscFunctionReturn(0);
53674ae819SStefano Zampini }
54674ae819SStefano Zampini 
5592cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrDestroy(void * ctx)
56674ae819SStefano Zampini {
5792cccca0SStefano Zampini   NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx;
58674ae819SStefano Zampini   PetscErrorCode          ierr;
59674ae819SStefano Zampini 
60674ae819SStefano Zampini   PetscFunctionBegin;
6192cccca0SStefano Zampini   ierr = VecDestroyVecs(3,&corr_ctx->sw);CHKERRQ(ierr);
6292cccca0SStefano Zampini   ierr = VecDestroyVecs(1,&corr_ctx->fw);CHKERRQ(ierr);
6392cccca0SStefano Zampini   ierr = MatDestroy(&corr_ctx->basis_mat);CHKERRQ(ierr);
6492cccca0SStefano Zampini   ierr = MatDestroy(&corr_ctx->inv_smat);CHKERRQ(ierr);
6592cccca0SStefano Zampini   ierr = PetscFree(corr_ctx);CHKERRQ(ierr);
66674ae819SStefano Zampini   PetscFunctionReturn(0);
67674ae819SStefano Zampini }
68674ae819SStefano Zampini 
69c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling)
70674ae819SStefano Zampini {
71674ae819SStefano Zampini   PC_BDDC                  *pcbddc = (PC_BDDC*)pc->data;
72c7017625SStefano Zampini   MatNullSpace             NullSpace = NULL;
73185763b3SStefano Zampini   KSP                      local_ksp;
74674ae819SStefano Zampini   NullSpaceCorrection_ctx  shell_ctx;
7592cccca0SStefano Zampini   Mat                      local_mat,local_pmat,dmat,Kbasis_mat;
7692cccca0SStefano Zampini   Vec                      v;
7792cccca0SStefano Zampini   PetscContainer           c;
7892cccca0SStefano Zampini   PetscInt                 basis_size;
7992cccca0SStefano Zampini   IS                       zerorows;
80674ae819SStefano Zampini   PetscErrorCode           ierr;
81674ae819SStefano Zampini 
82674ae819SStefano Zampini   PetscFunctionBegin;
838ead10e4SStefano Zampini   if (isdir) local_ksp = pcbddc->ksp_D; /* Dirichlet solver */
848ead10e4SStefano Zampini   else local_ksp = pcbddc->ksp_R; /* Neumann solver */
8592cccca0SStefano Zampini   ierr = KSPGetOperators(local_ksp,&local_mat,&local_pmat);CHKERRQ(ierr);
8692cccca0SStefano Zampini   ierr = MatGetNearNullSpace(local_pmat,&NullSpace);CHKERRQ(ierr);
872958b453SStefano Zampini   if (!NullSpace) {
8835529e7bSStefano Zampini     if (pcbddc->dbg_flag) {
892958b453SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d doesn't have local (near) nullspace: no need for correction in %s solver \n",PetscGlobalRank,isdir ? "Dirichlet" : "Neumann");CHKERRQ(ierr);
9035529e7bSStefano Zampini     }
9135529e7bSStefano Zampini     PetscFunctionReturn(0);
9235529e7bSStefano Zampini   }
9392cccca0SStefano Zampini   ierr = PetscObjectQuery((PetscObject)NullSpace,"_PBDDC_Null_dmat",(PetscObject*)&dmat);CHKERRQ(ierr);
9492cccca0SStefano Zampini   if (!dmat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing dense matrix");CHKERRQ(ierr);
958ead10e4SStefano Zampini   ierr = PetscLogEventBegin(PC_BDDC_ApproxSetUp[pcbddc->current_level],pc,0,0,0);CHKERRQ(ierr);
962958b453SStefano Zampini 
97854ce69bSBarry Smith   ierr = PetscNew(&shell_ctx);CHKERRQ(ierr);
98bd5e1169SStefano Zampini   shell_ctx->scale = 1.0;
9992cccca0SStefano Zampini   ierr = PetscObjectReference((PetscObject)dmat);CHKERRQ(ierr);
10092cccca0SStefano Zampini   shell_ctx->basis_mat = dmat;
10192cccca0SStefano Zampini   ierr = MatGetSize(dmat,NULL,&basis_size);CHKERRQ(ierr);
1028ead10e4SStefano Zampini   shell_ctx->evapply = PC_BDDC_ApproxApply[pcbddc->current_level];
10392cccca0SStefano Zampini 
104*7ebab0bbSStefano Zampini   ierr = MatGetOption(local_mat,MAT_SYMMETRIC,&shell_ctx->symm);CHKERRQ(ierr);
105*7ebab0bbSStefano Zampini 
10692cccca0SStefano Zampini   /* explicit construct (Phi^T K Phi)^-1 */
10792cccca0SStefano Zampini   ierr = MatMatMult(local_mat,shell_ctx->basis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Kbasis_mat);CHKERRQ(ierr);
10892cccca0SStefano Zampini   ierr = MatTransposeMatMult(Kbasis_mat,shell_ctx->basis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->inv_smat);CHKERRQ(ierr);
10992cccca0SStefano Zampini   ierr = MatDestroy(&Kbasis_mat);CHKERRQ(ierr);
11092cccca0SStefano Zampini   ierr = MatFindZeroRows(shell_ctx->inv_smat,&zerorows);CHKERRQ(ierr);
11192cccca0SStefano Zampini   if (zerorows) { /* linearly dependent basis */
11292cccca0SStefano Zampini     const PetscInt *idxs;
11392cccca0SStefano Zampini     PetscInt       i,nz;
11492cccca0SStefano Zampini 
11592cccca0SStefano Zampini     ierr = ISGetLocalSize(zerorows,&nz);CHKERRQ(ierr);
11692cccca0SStefano Zampini     ierr = ISGetIndices(zerorows,&idxs);CHKERRQ(ierr);
11792cccca0SStefano Zampini     for (i=0;i<nz;i++) {
11892cccca0SStefano Zampini       ierr = MatSetValue(shell_ctx->inv_smat,idxs[i],idxs[i],1.0,INSERT_VALUES);CHKERRQ(ierr);
11992cccca0SStefano Zampini     }
12092cccca0SStefano Zampini     ierr = ISRestoreIndices(zerorows,&idxs);CHKERRQ(ierr);
12192cccca0SStefano Zampini     ierr = MatAssemblyBegin(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12292cccca0SStefano Zampini     ierr = MatAssemblyEnd(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12392cccca0SStefano Zampini   }
12492cccca0SStefano Zampini 
12592cccca0SStefano Zampini   ierr = MatLUFactor(shell_ctx->inv_smat,NULL,NULL,NULL);CHKERRQ(ierr);
12692cccca0SStefano Zampini   ierr = MatSeqDenseInvertFactors_Private(shell_ctx->inv_smat);CHKERRQ(ierr);
12792cccca0SStefano Zampini   if (zerorows) { /* linearly dependent basis */
12892cccca0SStefano Zampini     const PetscInt *idxs;
12992cccca0SStefano Zampini     PetscInt       i,nz;
13092cccca0SStefano Zampini 
13192cccca0SStefano Zampini     ierr = ISGetLocalSize(zerorows,&nz);CHKERRQ(ierr);
13292cccca0SStefano Zampini     ierr = ISGetIndices(zerorows,&idxs);CHKERRQ(ierr);
13392cccca0SStefano Zampini     for (i=0;i<nz;i++) {
13492cccca0SStefano Zampini       ierr = MatSetValue(shell_ctx->inv_smat,idxs[i],idxs[i],0.0,INSERT_VALUES);CHKERRQ(ierr);
13592cccca0SStefano Zampini     }
13692cccca0SStefano Zampini     ierr = ISRestoreIndices(zerorows,&idxs);CHKERRQ(ierr);
13792cccca0SStefano Zampini     ierr = MatAssemblyBegin(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13892cccca0SStefano Zampini     ierr = MatAssemblyEnd(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13992cccca0SStefano Zampini   }
14092cccca0SStefano Zampini   ierr = ISDestroy(&zerorows);CHKERRQ(ierr);
141674ae819SStefano Zampini 
142674ae819SStefano Zampini   /* Create work vectors in shell context */
14392cccca0SStefano Zampini   ierr = MatCreateVecs(shell_ctx->inv_smat,&v,NULL);CHKERRQ(ierr);
14492cccca0SStefano Zampini   ierr = KSPCreateVecs(local_ksp,1,&shell_ctx->fw,0,NULL);CHKERRQ(ierr);
14592cccca0SStefano Zampini   ierr = VecDuplicateVecs(v,3,&shell_ctx->sw);CHKERRQ(ierr);
14692cccca0SStefano Zampini   ierr = VecDestroy(&v);CHKERRQ(ierr);
147674ae819SStefano Zampini 
14892cccca0SStefano Zampini   /* add special pre/post solve to KSP (see [1], eq. 48) */
14992cccca0SStefano Zampini   ierr = KSPSetPreSolve(local_ksp,PCBDDCNullSpaceCorrPreSolve,shell_ctx);CHKERRQ(ierr);
15092cccca0SStefano Zampini   ierr = KSPSetPostSolve(local_ksp,PCBDDCNullSpaceCorrPostSolve,shell_ctx);CHKERRQ(ierr);
15192cccca0SStefano Zampini   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)local_ksp),&c);CHKERRQ(ierr);
15292cccca0SStefano Zampini   ierr = PetscContainerSetPointer(c,shell_ctx);CHKERRQ(ierr);
15392cccca0SStefano Zampini   ierr = PetscContainerSetUserDestroy(c,PCBDDCNullSpaceCorrDestroy);CHKERRQ(ierr);
15492cccca0SStefano Zampini   ierr = PetscObjectCompose((PetscObject)local_ksp,"_PCBDDC_Null_PrePost_ctx",(PetscObject)c);CHKERRQ(ierr);
15592cccca0SStefano Zampini   ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
156c7017625SStefano Zampini 
157c7017625SStefano Zampini   /* Create ksp object suitable for extreme eigenvalues' estimation */
15892cccca0SStefano Zampini   if (needscaling || pcbddc->dbg_flag) {
159674ae819SStefano Zampini     KSP         check_ksp;
16092cccca0SStefano Zampini     PC          local_pc;
16192cccca0SStefano Zampini     Vec         work1,work2;
162bd5e1169SStefano Zampini     const char* prefix;
16392cccca0SStefano Zampini     PetscReal   test_err,lambda_min,lambda_max;
16492cccca0SStefano Zampini     PetscInt    k,maxit;
165c7017625SStefano Zampini 
16692cccca0SStefano Zampini     ierr = VecDuplicate(shell_ctx->fw[0],&work1);CHKERRQ(ierr);
16792cccca0SStefano Zampini     ierr = VecDuplicate(shell_ctx->fw[0],&work2);CHKERRQ(ierr);
168c7017625SStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr);
16992cccca0SStefano Zampini     if (local_mat->spd) {
17092cccca0SStefano Zampini       ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr);
17192cccca0SStefano Zampini     }
17215579a77SStefano Zampini     ierr = PetscObjectIncrementTabLevel((PetscObject)check_ksp,(PetscObject)local_ksp,0);CHKERRQ(ierr);
173bd5e1169SStefano Zampini     ierr = KSPGetOptionsPrefix(local_ksp,&prefix);CHKERRQ(ierr);
174bd5e1169SStefano Zampini     ierr = KSPSetOptionsPrefix(check_ksp,prefix);CHKERRQ(ierr);
17592cccca0SStefano Zampini     ierr = KSPAppendOptionsPrefix(check_ksp,"approximate_scale_");CHKERRQ(ierr);
176bd5e1169SStefano Zampini     ierr = KSPSetErrorIfNotConverged(check_ksp,PETSC_FALSE);CHKERRQ(ierr);
177c1a3ebd0SStefano Zampini     ierr = KSPSetOperators(check_ksp,local_mat,local_pmat);CHKERRQ(ierr);
178c7017625SStefano Zampini     ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr);
17992cccca0SStefano Zampini     ierr = KSPSetPreSolve(check_ksp,PCBDDCNullSpaceCorrPreSolve,shell_ctx);CHKERRQ(ierr);
18092cccca0SStefano Zampini     ierr = KSPSetPostSolve(check_ksp,PCBDDCNullSpaceCorrPostSolve,shell_ctx);CHKERRQ(ierr);
18192cccca0SStefano Zampini     ierr = KSPSetTolerances(check_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
182bd5e1169SStefano Zampini     ierr = KSPSetFromOptions(check_ksp);CHKERRQ(ierr);
18392cccca0SStefano 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 */
184c7017625SStefano Zampini     ierr = KSPSetUp(check_ksp);CHKERRQ(ierr);
18592cccca0SStefano Zampini     ierr = KSPGetPC(local_ksp,&local_pc);CHKERRQ(ierr);
18692cccca0SStefano Zampini     ierr = KSPSetPC(check_ksp,local_pc);CHKERRQ(ierr);
18792cccca0SStefano Zampini     ierr = KSPGetTolerances(check_ksp,NULL,NULL,NULL,&maxit);CHKERRQ(ierr);
18892cccca0SStefano Zampini     ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PetscMin(10,maxit));CHKERRQ(ierr);
18992cccca0SStefano Zampini     ierr = VecSetRandom(work2,NULL);CHKERRQ(ierr);
19092cccca0SStefano Zampini     ierr = MatMult(local_mat,work2,work1);CHKERRQ(ierr);
19192cccca0SStefano Zampini     ierr = KSPSolve(check_ksp,work1,work1);CHKERRQ(ierr);
19292cccca0SStefano Zampini     ierr = KSPCheckSolve(check_ksp,pc,work1);CHKERRQ(ierr);
19392cccca0SStefano Zampini     ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr);
19492cccca0SStefano Zampini     ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr);
195c7017625SStefano Zampini     ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr);
196c7017625SStefano Zampini     ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr);
197c7017625SStefano Zampini     if (pcbddc->dbg_flag) {
198c7017625SStefano Zampini       if (isdir) {
19992cccca0SStefano Zampini         ierr = 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);CHKERRQ(ierr);
200c7017625SStefano Zampini       } else {
20192cccca0SStefano Zampini         ierr = 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);CHKERRQ(ierr);
202c7017625SStefano Zampini       }
203c7017625SStefano Zampini     }
20492cccca0SStefano Zampini     if (needscaling) shell_ctx->scale = 1.0/lambda_max;
20592cccca0SStefano Zampini 
20692cccca0SStefano Zampini     if (needscaling && pcbddc->dbg_flag) { /* test for scaling factor */
20792cccca0SStefano Zampini       PC new_pc;
20892cccca0SStefano Zampini 
20992cccca0SStefano Zampini       ierr = VecSetRandom(work2,NULL);CHKERRQ(ierr);
21092cccca0SStefano Zampini       ierr = MatMult(local_mat,work2,work1);CHKERRQ(ierr);
21192cccca0SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)check_ksp),&new_pc);CHKERRQ(ierr);
21292cccca0SStefano Zampini       ierr = PCSetType(new_pc,PCKSP);CHKERRQ(ierr);
21392cccca0SStefano Zampini       ierr = PCSetOperators(new_pc,local_mat,local_pmat);CHKERRQ(ierr);
21492cccca0SStefano Zampini       ierr = PCKSPSetKSP(new_pc,local_ksp);CHKERRQ(ierr);
21592cccca0SStefano Zampini       ierr = KSPSetPC(check_ksp,new_pc);CHKERRQ(ierr);
21692cccca0SStefano Zampini       ierr = PCDestroy(&new_pc);CHKERRQ(ierr);
21792cccca0SStefano Zampini       ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit);CHKERRQ(ierr);
21892cccca0SStefano Zampini       ierr = KSPSetPreSolve(check_ksp,NULL,NULL);CHKERRQ(ierr);
21992cccca0SStefano Zampini       ierr = KSPSetPostSolve(check_ksp,NULL,NULL);CHKERRQ(ierr);
22092cccca0SStefano Zampini       ierr = KSPSolve(check_ksp,work1,work1);CHKERRQ(ierr);
22192cccca0SStefano Zampini       ierr = KSPCheckSolve(check_ksp,pc,work1);CHKERRQ(ierr);
22292cccca0SStefano Zampini       ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr);
22392cccca0SStefano Zampini       ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr);
22492cccca0SStefano Zampini       ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr);
22592cccca0SStefano Zampini       ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr);
22692cccca0SStefano Zampini       if (pcbddc->dbg_flag) {
22792cccca0SStefano Zampini         if (isdir) {
22892cccca0SStefano Zampini           ierr = 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);CHKERRQ(ierr);
22992cccca0SStefano Zampini         } else {
23092cccca0SStefano Zampini           ierr = 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);CHKERRQ(ierr);
23192cccca0SStefano Zampini         }
23292cccca0SStefano Zampini       }
23392cccca0SStefano Zampini     }
234c7017625SStefano Zampini     ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr);
23592cccca0SStefano Zampini     ierr = VecDestroy(&work1);CHKERRQ(ierr);
23692cccca0SStefano Zampini     ierr = VecDestroy(&work2);CHKERRQ(ierr);
237c7017625SStefano Zampini   }
2388ead10e4SStefano Zampini   ierr = PetscLogEventEnd(PC_BDDC_ApproxSetUp[pcbddc->current_level],pc,0,0,0);CHKERRQ(ierr);
239c7017625SStefano Zampini 
24092cccca0SStefano Zampini   if (pcbddc->dbg_flag) {
241c7017625SStefano Zampini     Vec       work1,work2,work3;
24292cccca0SStefano Zampini     PetscReal test_err;
243674ae819SStefano Zampini 
24492cccca0SStefano Zampini     /* check nullspace basis is solved exactly */
24592cccca0SStefano Zampini     ierr = VecDuplicate(shell_ctx->fw[0],&work1);CHKERRQ(ierr);
24692cccca0SStefano Zampini     ierr = VecDuplicate(shell_ctx->fw[0],&work2);CHKERRQ(ierr);
24792cccca0SStefano Zampini     ierr = VecDuplicate(shell_ctx->fw[0],&work3);CHKERRQ(ierr);
24892cccca0SStefano Zampini     ierr = VecSetRandom(shell_ctx->sw[0],NULL);CHKERRQ(ierr);
24992cccca0SStefano Zampini     ierr = MatMult(shell_ctx->basis_mat,shell_ctx->sw[0],work1);CHKERRQ(ierr);
250674ae819SStefano Zampini     ierr = VecCopy(work1,work2);CHKERRQ(ierr);
251674ae819SStefano Zampini     ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr);
25292cccca0SStefano Zampini     ierr = KSPSolve(local_ksp,work3,work1);CHKERRQ(ierr);
253c7017625SStefano Zampini     ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr);
254674ae819SStefano Zampini     ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr);
255185763b3SStefano Zampini     if (isdir) {
256c7017625SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
257674ae819SStefano Zampini     } else {
258c7017625SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
259674ae819SStefano Zampini     }
260674ae819SStefano Zampini     ierr = VecDestroy(&work1);CHKERRQ(ierr);
261674ae819SStefano Zampini     ierr = VecDestroy(&work2);CHKERRQ(ierr);
262674ae819SStefano Zampini     ierr = VecDestroy(&work3);CHKERRQ(ierr);
26392cccca0SStefano Zampini   }
264674ae819SStefano Zampini   PetscFunctionReturn(0);
265674ae819SStefano Zampini }
266