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