xref: /petsc/src/ksp/pc/impls/bddc/bddcnullspace.c (revision 92cccca0902b191cd1bb1f99cb527a8d18161527)
1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3*92cccca0SStefano Zampini #include <../src/mat/impls/dense/seq/dense.h>
4674ae819SStefano Zampini 
5*92cccca0SStefano Zampini /* E + small_solve */
6*92cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrPreSolve(KSP ksp,Vec y,Vec x, void* ctx)
7bd5e1169SStefano Zampini {
8*92cccca0SStefano Zampini   NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx;
9*92cccca0SStefano Zampini   Mat                     K;
10bd5e1169SStefano Zampini   PetscErrorCode          ierr;
11bd5e1169SStefano Zampini 
12bd5e1169SStefano Zampini   PetscFunctionBegin;
13*92cccca0SStefano Zampini   ierr = MatMultTranspose(corr_ctx->basis_mat,y,corr_ctx->sw[0]);CHKERRQ(ierr);
14*92cccca0SStefano Zampini   ierr = MatMultTranspose(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[1]);CHKERRQ(ierr);
15*92cccca0SStefano Zampini   ierr = VecScale(corr_ctx->sw[1],-1.0);CHKERRQ(ierr);
16*92cccca0SStefano Zampini   ierr = MatMult(corr_ctx->basis_mat,corr_ctx->sw[1],corr_ctx->fw[0]);CHKERRQ(ierr);
17*92cccca0SStefano Zampini   ierr = VecScale(corr_ctx->sw[1],-1.0);CHKERRQ(ierr);
18*92cccca0SStefano Zampini   ierr = KSPGetOperators(ksp,&K,NULL);CHKERRQ(ierr);
19*92cccca0SStefano Zampini   ierr = MatMultAdd(K,corr_ctx->fw[0],y,y);CHKERRQ(ierr);
20bd5e1169SStefano Zampini   PetscFunctionReturn(0);
21bd5e1169SStefano Zampini }
22bd5e1169SStefano Zampini 
23*92cccca0SStefano Zampini /* E^t + small */
24*92cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrPostSolve(KSP ksp,Vec y,Vec x, void* ctx)
25674ae819SStefano Zampini {
26*92cccca0SStefano Zampini   NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx;
27674ae819SStefano Zampini   PetscErrorCode          ierr;
28*92cccca0SStefano Zampini   Mat                     K;
29674ae819SStefano Zampini 
30674ae819SStefano Zampini   PetscFunctionBegin;
31*92cccca0SStefano Zampini   ierr = KSPGetOperators(ksp,&K,NULL);CHKERRQ(ierr);
32*92cccca0SStefano Zampini   ierr = MatMultTranspose(K,x,corr_ctx->fw[0]);CHKERRQ(ierr);
33*92cccca0SStefano Zampini   ierr = MatMultTranspose(corr_ctx->basis_mat,corr_ctx->fw[0],corr_ctx->sw[0]);CHKERRQ(ierr);
34*92cccca0SStefano Zampini   ierr = VecScale(corr_ctx->sw[0],-1.0);CHKERRQ(ierr);
35*92cccca0SStefano Zampini   ierr = MatMult(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[2]);CHKERRQ(ierr);
36*92cccca0SStefano Zampini   ierr = MatMultAdd(corr_ctx->basis_mat,corr_ctx->sw[2],x,corr_ctx->fw[0]);CHKERRQ(ierr);
37*92cccca0SStefano Zampini   ierr = VecScale(corr_ctx->fw[0],corr_ctx->scale);CHKERRQ(ierr);
38*92cccca0SStefano Zampini   /* Sum contributions from approximate solver and projected system */
39*92cccca0SStefano Zampini   ierr = MatMultAdd(corr_ctx->basis_mat,corr_ctx->sw[1],corr_ctx->fw[0],x);CHKERRQ(ierr);
40674ae819SStefano Zampini   PetscFunctionReturn(0);
41674ae819SStefano Zampini }
42674ae819SStefano Zampini 
43*92cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrDestroy(void * ctx)
44674ae819SStefano Zampini {
45*92cccca0SStefano Zampini   NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx;
46674ae819SStefano Zampini   PetscErrorCode          ierr;
47674ae819SStefano Zampini 
48674ae819SStefano Zampini   PetscFunctionBegin;
49*92cccca0SStefano Zampini   ierr = VecDestroyVecs(3,&corr_ctx->sw);CHKERRQ(ierr);
50*92cccca0SStefano Zampini   ierr = VecDestroyVecs(1,&corr_ctx->fw);CHKERRQ(ierr);
51*92cccca0SStefano Zampini   ierr = MatDestroy(&corr_ctx->basis_mat);CHKERRQ(ierr);
52*92cccca0SStefano Zampini   ierr = MatDestroy(&corr_ctx->inv_smat);CHKERRQ(ierr);
53*92cccca0SStefano Zampini   ierr = PetscFree(corr_ctx);CHKERRQ(ierr);
54674ae819SStefano Zampini   PetscFunctionReturn(0);
55674ae819SStefano Zampini }
56674ae819SStefano Zampini 
57c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling)
58674ae819SStefano Zampini {
59674ae819SStefano Zampini   PC_BDDC                  *pcbddc = (PC_BDDC*)pc->data;
60c7017625SStefano Zampini   MatNullSpace             NullSpace = NULL;
61185763b3SStefano Zampini   KSP                      local_ksp;
62674ae819SStefano Zampini   NullSpaceCorrection_ctx  shell_ctx;
63*92cccca0SStefano Zampini   Mat                      local_mat,local_pmat,dmat,Kbasis_mat;
64*92cccca0SStefano Zampini   Vec                      v;
65*92cccca0SStefano Zampini   PetscContainer           c;
66*92cccca0SStefano Zampini   PetscInt                 basis_size;
67*92cccca0SStefano Zampini   IS                       zerorows;
68674ae819SStefano Zampini   PetscErrorCode           ierr;
69674ae819SStefano Zampini 
70674ae819SStefano Zampini   PetscFunctionBegin;
71*92cccca0SStefano Zampini   if (isdir) { /* Dirichlet solver */
72*92cccca0SStefano Zampini     local_ksp = pcbddc->ksp_D;
73*92cccca0SStefano Zampini   } else { /* Neumann solver */
74*92cccca0SStefano Zampini     local_ksp = pcbddc->ksp_R;
752958b453SStefano Zampini   }
76*92cccca0SStefano Zampini   ierr = KSPGetOperators(local_ksp,&local_mat,&local_pmat);CHKERRQ(ierr);
77*92cccca0SStefano Zampini   ierr = MatGetNearNullSpace(local_pmat,&NullSpace);CHKERRQ(ierr);
782958b453SStefano Zampini   if (!NullSpace) {
7935529e7bSStefano Zampini     if (pcbddc->dbg_flag) {
802958b453SStefano 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);
8135529e7bSStefano Zampini     }
8235529e7bSStefano Zampini     PetscFunctionReturn(0);
8335529e7bSStefano Zampini   }
84*92cccca0SStefano Zampini   ierr = PetscObjectQuery((PetscObject)NullSpace,"_PBDDC_Null_dmat",(PetscObject*)&dmat);CHKERRQ(ierr);
85*92cccca0SStefano Zampini   if (!dmat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing dense matrix");CHKERRQ(ierr);
862958b453SStefano Zampini 
87854ce69bSBarry Smith   ierr = PetscNew(&shell_ctx);CHKERRQ(ierr);
88bd5e1169SStefano Zampini   shell_ctx->scale = 1.0;
89*92cccca0SStefano Zampini   ierr = PetscObjectReference((PetscObject)dmat);CHKERRQ(ierr);
90*92cccca0SStefano Zampini   shell_ctx->basis_mat = dmat;
91*92cccca0SStefano Zampini   ierr = MatGetSize(dmat,NULL,&basis_size);CHKERRQ(ierr);
92*92cccca0SStefano Zampini 
93*92cccca0SStefano Zampini   /* explicit construct (Phi^T K Phi)^-1 */
94*92cccca0SStefano Zampini   ierr = MatMatMult(local_mat,shell_ctx->basis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Kbasis_mat);CHKERRQ(ierr);
95*92cccca0SStefano Zampini   ierr = MatTransposeMatMult(Kbasis_mat,shell_ctx->basis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->inv_smat);CHKERRQ(ierr);
96*92cccca0SStefano Zampini   ierr = MatDestroy(&Kbasis_mat);CHKERRQ(ierr);
97*92cccca0SStefano Zampini   ierr = MatFindZeroRows(shell_ctx->inv_smat,&zerorows);CHKERRQ(ierr);
98*92cccca0SStefano Zampini   if (zerorows) { /* linearly dependent basis */
99*92cccca0SStefano Zampini     const PetscInt *idxs;
100*92cccca0SStefano Zampini     PetscInt       i,nz;
101*92cccca0SStefano Zampini 
102*92cccca0SStefano Zampini     ierr = ISGetLocalSize(zerorows,&nz);CHKERRQ(ierr);
103*92cccca0SStefano Zampini     ierr = ISGetIndices(zerorows,&idxs);CHKERRQ(ierr);
104*92cccca0SStefano Zampini     for (i=0;i<nz;i++) {
105*92cccca0SStefano Zampini       ierr = MatSetValue(shell_ctx->inv_smat,idxs[i],idxs[i],1.0,INSERT_VALUES);CHKERRQ(ierr);
106*92cccca0SStefano Zampini     }
107*92cccca0SStefano Zampini     ierr = ISRestoreIndices(zerorows,&idxs);CHKERRQ(ierr);
108*92cccca0SStefano Zampini     ierr = MatAssemblyBegin(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
109*92cccca0SStefano Zampini     ierr = MatAssemblyEnd(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
110*92cccca0SStefano Zampini   }
111*92cccca0SStefano Zampini 
112*92cccca0SStefano Zampini   ierr = MatLUFactor(shell_ctx->inv_smat,NULL,NULL,NULL);CHKERRQ(ierr);
113*92cccca0SStefano Zampini   ierr = MatSeqDenseInvertFactors_Private(shell_ctx->inv_smat);CHKERRQ(ierr);
114*92cccca0SStefano Zampini   if (zerorows) { /* linearly dependent basis */
115*92cccca0SStefano Zampini     const PetscInt *idxs;
116*92cccca0SStefano Zampini     PetscInt       i,nz;
117*92cccca0SStefano Zampini 
118*92cccca0SStefano Zampini     ierr = ISGetLocalSize(zerorows,&nz);CHKERRQ(ierr);
119*92cccca0SStefano Zampini     ierr = ISGetIndices(zerorows,&idxs);CHKERRQ(ierr);
120*92cccca0SStefano Zampini     for (i=0;i<nz;i++) {
121*92cccca0SStefano Zampini       ierr = MatSetValue(shell_ctx->inv_smat,idxs[i],idxs[i],0.0,INSERT_VALUES);CHKERRQ(ierr);
122*92cccca0SStefano Zampini     }
123*92cccca0SStefano Zampini     ierr = ISRestoreIndices(zerorows,&idxs);CHKERRQ(ierr);
124*92cccca0SStefano Zampini     ierr = MatAssemblyBegin(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
125*92cccca0SStefano Zampini     ierr = MatAssemblyEnd(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
126*92cccca0SStefano Zampini   }
127*92cccca0SStefano Zampini   ierr = ISDestroy(&zerorows);CHKERRQ(ierr);
128674ae819SStefano Zampini 
129674ae819SStefano Zampini   /* Create work vectors in shell context */
130*92cccca0SStefano Zampini   ierr = MatCreateVecs(shell_ctx->inv_smat,&v,NULL);CHKERRQ(ierr);
131*92cccca0SStefano Zampini   ierr = KSPCreateVecs(local_ksp,1,&shell_ctx->fw,0,NULL);CHKERRQ(ierr);
132*92cccca0SStefano Zampini   ierr = VecDuplicateVecs(v,3,&shell_ctx->sw);CHKERRQ(ierr);
133*92cccca0SStefano Zampini   ierr = VecDestroy(&v);CHKERRQ(ierr);
134674ae819SStefano Zampini 
135*92cccca0SStefano Zampini   /* add special pre/post solve to KSP (see [1], eq. 48) */
136*92cccca0SStefano Zampini   ierr = KSPSetPreSolve(local_ksp,PCBDDCNullSpaceCorrPreSolve,shell_ctx);CHKERRQ(ierr);
137*92cccca0SStefano Zampini   ierr = KSPSetPostSolve(local_ksp,PCBDDCNullSpaceCorrPostSolve,shell_ctx);CHKERRQ(ierr);
138*92cccca0SStefano Zampini   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)local_ksp),&c);CHKERRQ(ierr);
139*92cccca0SStefano Zampini   ierr = PetscContainerSetPointer(c,shell_ctx);CHKERRQ(ierr);
140*92cccca0SStefano Zampini   ierr = PetscContainerSetUserDestroy(c,PCBDDCNullSpaceCorrDestroy);CHKERRQ(ierr);
141*92cccca0SStefano Zampini   ierr = PetscObjectCompose((PetscObject)local_ksp,"_PCBDDC_Null_PrePost_ctx",(PetscObject)c);CHKERRQ(ierr);
142*92cccca0SStefano Zampini   ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
143c7017625SStefano Zampini 
144c7017625SStefano Zampini   /* Create ksp object suitable for extreme eigenvalues' estimation */
145*92cccca0SStefano Zampini   if (needscaling || pcbddc->dbg_flag) {
146674ae819SStefano Zampini     KSP         check_ksp;
147*92cccca0SStefano Zampini     PC          local_pc;
148*92cccca0SStefano Zampini     Vec         work1,work2;
149bd5e1169SStefano Zampini     const char* prefix;
150*92cccca0SStefano Zampini     PetscReal   test_err,lambda_min,lambda_max;
151*92cccca0SStefano Zampini     PetscInt    k,maxit;
152c7017625SStefano Zampini 
153*92cccca0SStefano Zampini     ierr = VecDuplicate(shell_ctx->fw[0],&work1);CHKERRQ(ierr);
154*92cccca0SStefano Zampini     ierr = VecDuplicate(shell_ctx->fw[0],&work2);CHKERRQ(ierr);
155c7017625SStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr);
156*92cccca0SStefano Zampini     if (local_mat->spd) {
157*92cccca0SStefano Zampini       ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr);
158*92cccca0SStefano Zampini     }
15915579a77SStefano Zampini     ierr = PetscObjectIncrementTabLevel((PetscObject)check_ksp,(PetscObject)local_ksp,0);CHKERRQ(ierr);
160bd5e1169SStefano Zampini     ierr = KSPGetOptionsPrefix(local_ksp,&prefix);CHKERRQ(ierr);
161bd5e1169SStefano Zampini     ierr = KSPSetOptionsPrefix(check_ksp,prefix);CHKERRQ(ierr);
162*92cccca0SStefano Zampini     ierr = KSPAppendOptionsPrefix(check_ksp,"approximate_scale_");CHKERRQ(ierr);
163bd5e1169SStefano Zampini     ierr = KSPSetErrorIfNotConverged(check_ksp,PETSC_FALSE);CHKERRQ(ierr);
164c1a3ebd0SStefano Zampini     ierr = KSPSetOperators(check_ksp,local_mat,local_pmat);CHKERRQ(ierr);
165c7017625SStefano Zampini     ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr);
166*92cccca0SStefano Zampini     ierr = KSPSetPreSolve(check_ksp,PCBDDCNullSpaceCorrPreSolve,shell_ctx);CHKERRQ(ierr);
167*92cccca0SStefano Zampini     ierr = KSPSetPostSolve(check_ksp,PCBDDCNullSpaceCorrPostSolve,shell_ctx);CHKERRQ(ierr);
168*92cccca0SStefano Zampini     ierr = KSPSetTolerances(check_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
169bd5e1169SStefano Zampini     ierr = KSPSetFromOptions(check_ksp);CHKERRQ(ierr);
170*92cccca0SStefano 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 */
171c7017625SStefano Zampini     ierr = KSPSetUp(check_ksp);CHKERRQ(ierr);
172*92cccca0SStefano Zampini     ierr = KSPGetPC(local_ksp,&local_pc);CHKERRQ(ierr);
173*92cccca0SStefano Zampini     ierr = KSPSetPC(check_ksp,local_pc);CHKERRQ(ierr);
174*92cccca0SStefano Zampini     ierr = KSPGetTolerances(check_ksp,NULL,NULL,NULL,&maxit);CHKERRQ(ierr);
175*92cccca0SStefano Zampini     ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PetscMin(10,maxit));CHKERRQ(ierr);
176*92cccca0SStefano Zampini     ierr = VecSetRandom(work2,NULL);CHKERRQ(ierr);
177*92cccca0SStefano Zampini     ierr = MatMult(local_mat,work2,work1);CHKERRQ(ierr);
178*92cccca0SStefano Zampini     ierr = KSPSolve(check_ksp,work1,work1);CHKERRQ(ierr);
179*92cccca0SStefano Zampini     ierr = KSPCheckSolve(check_ksp,pc,work1);CHKERRQ(ierr);
180*92cccca0SStefano Zampini     ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr);
181*92cccca0SStefano Zampini     ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr);
182c7017625SStefano Zampini     ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr);
183c7017625SStefano Zampini     ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr);
184c7017625SStefano Zampini     if (pcbddc->dbg_flag) {
185c7017625SStefano Zampini       if (isdir) {
186*92cccca0SStefano 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);
187c7017625SStefano Zampini       } else {
188*92cccca0SStefano 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);
189c7017625SStefano Zampini       }
190c7017625SStefano Zampini     }
191*92cccca0SStefano Zampini     if (needscaling) shell_ctx->scale = 1.0/lambda_max;
192*92cccca0SStefano Zampini 
193*92cccca0SStefano Zampini     if (needscaling && pcbddc->dbg_flag) { /* test for scaling factor */
194*92cccca0SStefano Zampini       PC new_pc;
195*92cccca0SStefano Zampini 
196*92cccca0SStefano Zampini       ierr = VecSetRandom(work2,NULL);CHKERRQ(ierr);
197*92cccca0SStefano Zampini       ierr = MatMult(local_mat,work2,work1);CHKERRQ(ierr);
198*92cccca0SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)check_ksp),&new_pc);CHKERRQ(ierr);
199*92cccca0SStefano Zampini       ierr = PCSetType(new_pc,PCKSP);CHKERRQ(ierr);
200*92cccca0SStefano Zampini       ierr = PCSetOperators(new_pc,local_mat,local_pmat);CHKERRQ(ierr);
201*92cccca0SStefano Zampini       ierr = PCKSPSetKSP(new_pc,local_ksp);CHKERRQ(ierr);
202*92cccca0SStefano Zampini       ierr = KSPSetPC(check_ksp,new_pc);CHKERRQ(ierr);
203*92cccca0SStefano Zampini       ierr = PCDestroy(&new_pc);CHKERRQ(ierr);
204*92cccca0SStefano Zampini       ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit);CHKERRQ(ierr);
205*92cccca0SStefano Zampini       ierr = KSPSetPreSolve(check_ksp,NULL,NULL);CHKERRQ(ierr);
206*92cccca0SStefano Zampini       ierr = KSPSetPostSolve(check_ksp,NULL,NULL);CHKERRQ(ierr);
207*92cccca0SStefano Zampini       ierr = KSPSolve(check_ksp,work1,work1);CHKERRQ(ierr);
208*92cccca0SStefano Zampini       ierr = KSPCheckSolve(check_ksp,pc,work1);CHKERRQ(ierr);
209*92cccca0SStefano Zampini       ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr);
210*92cccca0SStefano Zampini       ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr);
211*92cccca0SStefano Zampini       ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr);
212*92cccca0SStefano Zampini       ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr);
213*92cccca0SStefano Zampini       if (pcbddc->dbg_flag) {
214*92cccca0SStefano Zampini         if (isdir) {
215*92cccca0SStefano 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);
216*92cccca0SStefano Zampini         } else {
217*92cccca0SStefano 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);
218*92cccca0SStefano Zampini         }
219*92cccca0SStefano Zampini       }
220*92cccca0SStefano Zampini     }
221c7017625SStefano Zampini     ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr);
222*92cccca0SStefano Zampini     ierr = VecDestroy(&work1);CHKERRQ(ierr);
223*92cccca0SStefano Zampini     ierr = VecDestroy(&work2);CHKERRQ(ierr);
224c7017625SStefano Zampini   }
225c7017625SStefano Zampini 
226*92cccca0SStefano Zampini   if (pcbddc->dbg_flag) {
227c7017625SStefano Zampini     Vec       work1,work2,work3;
228*92cccca0SStefano Zampini     PetscReal test_err;
229674ae819SStefano Zampini 
230*92cccca0SStefano Zampini     /* check nullspace basis is solved exactly */
231*92cccca0SStefano Zampini     ierr = VecDuplicate(shell_ctx->fw[0],&work1);CHKERRQ(ierr);
232*92cccca0SStefano Zampini     ierr = VecDuplicate(shell_ctx->fw[0],&work2);CHKERRQ(ierr);
233*92cccca0SStefano Zampini     ierr = VecDuplicate(shell_ctx->fw[0],&work3);CHKERRQ(ierr);
234*92cccca0SStefano Zampini     ierr = VecSetRandom(shell_ctx->sw[0],NULL);CHKERRQ(ierr);
235*92cccca0SStefano Zampini     ierr = MatMult(shell_ctx->basis_mat,shell_ctx->sw[0],work1);CHKERRQ(ierr);
236674ae819SStefano Zampini     ierr = VecCopy(work1,work2);CHKERRQ(ierr);
237674ae819SStefano Zampini     ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr);
238*92cccca0SStefano Zampini     ierr = KSPSolve(local_ksp,work3,work1);CHKERRQ(ierr);
239c7017625SStefano Zampini     ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr);
240674ae819SStefano Zampini     ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr);
241185763b3SStefano Zampini     if (isdir) {
242c7017625SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
243674ae819SStefano Zampini     } else {
244c7017625SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
245674ae819SStefano Zampini     }
246674ae819SStefano Zampini     ierr = VecDestroy(&work1);CHKERRQ(ierr);
247674ae819SStefano Zampini     ierr = VecDestroy(&work2);CHKERRQ(ierr);
248674ae819SStefano Zampini     ierr = VecDestroy(&work3);CHKERRQ(ierr);
249*92cccca0SStefano Zampini   }
250674ae819SStefano Zampini   PetscFunctionReturn(0);
251674ae819SStefano Zampini }
252