xref: /petsc/src/ksp/pc/impls/bddc/bddcnullspace.c (revision 3cbbc5ff037ae03fd9010ee4a8e4d7c21ef7cb71) !
1 #include <../src/ksp/pc/impls/bddc/bddc.h>
2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3 
4 static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y)
5 {
6   NullSpaceCorrection_ctx pc_ctx;
7   PetscErrorCode          ierr;
8 
9   PetscFunctionBegin;
10   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
11   /* E */
12   ierr = MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);CHKERRQ(ierr);
13   ierr = MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);CHKERRQ(ierr);
14   /* P^-1 */
15   ierr = PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);CHKERRQ(ierr);
16   /* E^T */
17   ierr = MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);CHKERRQ(ierr);
18   ierr = VecScale(pc_ctx->work_small_1,-1.0);CHKERRQ(ierr);
19   ierr = MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);CHKERRQ(ierr);
20   /* Sum contributions */
21   ierr = MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);CHKERRQ(ierr);
22   if (pc_ctx->apply_scaling) {
23     ierr = VecScale(y,pc_ctx->scale);CHKERRQ(ierr);
24   }
25   PetscFunctionReturn(0);
26 }
27 
28 static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc)
29 {
30   NullSpaceCorrection_ctx pc_ctx;
31   PetscErrorCode          ierr;
32 
33   PetscFunctionBegin;
34   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
35   ierr = VecDestroy(&pc_ctx->work_small_1);CHKERRQ(ierr);
36   ierr = VecDestroy(&pc_ctx->work_small_2);CHKERRQ(ierr);
37   ierr = VecDestroy(&pc_ctx->work_full_1);CHKERRQ(ierr);
38   ierr = VecDestroy(&pc_ctx->work_full_2);CHKERRQ(ierr);
39   ierr = MatDestroy(&pc_ctx->basis_mat);CHKERRQ(ierr);
40   ierr = MatDestroy(&pc_ctx->Lbasis_mat);CHKERRQ(ierr);
41   ierr = MatDestroy(&pc_ctx->Kbasis_mat);CHKERRQ(ierr);
42   ierr = PCDestroy(&pc_ctx->local_pc);CHKERRQ(ierr);
43   ierr = PetscFree(pc_ctx);CHKERRQ(ierr);
44   PetscFunctionReturn(0);
45 }
46 
47 PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling)
48 {
49   PC_BDDC                  *pcbddc = (PC_BDDC*)pc->data;
50   PC_IS                    *pcis = (PC_IS*)pc->data;
51   Mat_IS                   *matis = (Mat_IS*)pc->pmat->data;
52   MatNullSpace             NullSpace = NULL;
53   KSP                      local_ksp;
54   PC                       newpc;
55   NullSpaceCorrection_ctx  shell_ctx;
56   Mat                      local_mat,local_pmat,small_mat,inv_small_mat;
57   const Vec                *nullvecs;
58   VecScatter               scatter_ctx;
59   IS                       is_aux,local_dofs;
60   MatFactorInfo            matinfo;
61   PetscScalar              *basis_mat,*Kbasis_mat,*array,*array_mat;
62   PetscScalar              one = 1.0,zero = 0.0, m_one = -1.0;
63   PetscInt                 basis_dofs,basis_size,nnsp_size,i,k;
64   PetscBool                nnsp_has_cnst;
65   PetscReal                test_err,lambda_min,lambda_max;
66   PetscBool                setsym,issym=PETSC_FALSE;
67   PetscErrorCode           ierr;
68 
69   PetscFunctionBegin;
70   ierr = MatGetNullSpace(matis->A,&NullSpace);CHKERRQ(ierr);
71   if (!NullSpace) {
72     if (pcbddc->dbg_flag) {
73       ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d doesn't have local nullspace: no need for correction in %s solver \n",PetscGlobalRank,isdir ? "Dirichlet" : "Neumann");CHKERRQ(ierr);
74     }
75     PetscFunctionReturn(0);
76   }
77   /* Infer the local solver */
78   if (isdir) {
79     /* Dirichlet solver */
80     local_ksp = pcbddc->ksp_D;
81     local_dofs = pcis->is_I_local;
82   } else {
83     /* Neumann solver */
84     local_ksp = pcbddc->ksp_R;
85     local_dofs = pcbddc->is_R_local;
86   }
87   ierr = ISGetSize(local_dofs,&basis_dofs);CHKERRQ(ierr);
88   ierr = KSPGetOperators(local_ksp,&local_mat,&local_pmat);CHKERRQ(ierr);
89 
90   /* Get null space vecs */
91   ierr = MatNullSpaceGetVecs(NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);CHKERRQ(ierr);
92   basis_size = nnsp_size;
93   if (nnsp_has_cnst) basis_size++;
94 
95    /* Create shell ctx */
96   ierr = PetscNew(&shell_ctx);CHKERRQ(ierr);
97   shell_ctx->apply_scaling = needscaling;
98 
99   /* Create work vectors in shell context */
100   ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);CHKERRQ(ierr);
101   ierr = VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);CHKERRQ(ierr);
102   ierr = VecSetType(shell_ctx->work_small_1,VECSEQ);CHKERRQ(ierr);
103   ierr = VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);CHKERRQ(ierr);
104   ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);CHKERRQ(ierr);
105   ierr = VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);CHKERRQ(ierr);
106   ierr = VecSetType(shell_ctx->work_full_1,VECSEQ);CHKERRQ(ierr);
107   ierr = VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);CHKERRQ(ierr);
108 
109   /* Allocate workspace */
110   ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat);CHKERRQ(ierr);
111   ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);CHKERRQ(ierr);
112   ierr = MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr);
113   ierr = MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr);
114 
115   /* Restrict local null space on selected dofs
116      and compute matrices N and K*N */
117   ierr = VecScatterCreate(pcis->vec1_N,local_dofs,shell_ctx->work_full_1,(IS)0,&scatter_ctx);CHKERRQ(ierr);
118   for (k=0;k<nnsp_size;k++) {
119     ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr);
120     ierr = VecScatterBegin(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
121     ierr = VecScatterEnd(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
122     ierr = VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr);
123     ierr = MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);CHKERRQ(ierr);
124     ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr);
125     ierr = VecResetArray(shell_ctx->work_full_2);CHKERRQ(ierr);
126   }
127   if (nnsp_has_cnst) {
128     ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr);
129     ierr = VecSet(shell_ctx->work_full_1,one);CHKERRQ(ierr);
130     ierr = VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr);
131     ierr = MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);CHKERRQ(ierr);
132     ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr);
133     ierr = VecResetArray(shell_ctx->work_full_2);CHKERRQ(ierr);
134   }
135   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
136   ierr = MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr);
137   ierr = MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr);
138 
139   /* Assemble another Mat object in shell context */
140   ierr = MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);CHKERRQ(ierr);
141   ierr = MatFactorInfoInitialize(&matinfo);CHKERRQ(ierr);
142   ierr = ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);CHKERRQ(ierr);
143   ierr = MatLUFactor(small_mat,is_aux,is_aux,&matinfo);CHKERRQ(ierr);
144   ierr = ISDestroy(&is_aux);CHKERRQ(ierr);
145   ierr = PetscMalloc1(basis_size*basis_size,&array_mat);CHKERRQ(ierr);
146   for (k=0;k<basis_size;k++) {
147     ierr = VecSet(shell_ctx->work_small_1,zero);CHKERRQ(ierr);
148     ierr = VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);CHKERRQ(ierr);
149     ierr = VecAssemblyBegin(shell_ctx->work_small_1);CHKERRQ(ierr);
150     ierr = VecAssemblyEnd(shell_ctx->work_small_1);CHKERRQ(ierr);
151     ierr = MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);CHKERRQ(ierr);
152     ierr = VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr);
153     for (i=0;i<basis_size;i++) {
154       array_mat[i*basis_size+k]=array[i];
155     }
156     ierr = VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr);
157   }
158   ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);CHKERRQ(ierr);
159   ierr = MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);CHKERRQ(ierr);
160   ierr = PetscFree(array_mat);CHKERRQ(ierr);
161   ierr = MatDestroy(&inv_small_mat);CHKERRQ(ierr);
162   ierr = MatDestroy(&small_mat);CHKERRQ(ierr);
163   ierr = MatScale(shell_ctx->Kbasis_mat,m_one);CHKERRQ(ierr);
164 
165   /* Rebuild local PC */
166   ierr = KSPGetPC(local_ksp,&shell_ctx->local_pc);CHKERRQ(ierr);
167   ierr = PetscObjectReference((PetscObject)shell_ctx->local_pc);CHKERRQ(ierr);
168   ierr = PCCreate(PETSC_COMM_SELF,&newpc);CHKERRQ(ierr);
169   ierr = PCSetOperators(newpc,local_mat,local_mat);CHKERRQ(ierr);
170   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
171   ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr);
172   ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr);
173   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr);
174   ierr = PCSetUp(newpc);CHKERRQ(ierr);
175   ierr = KSPSetPC(local_ksp,newpc);CHKERRQ(ierr);
176   ierr = KSPSetUp(local_ksp);CHKERRQ(ierr);
177 
178   /* Create ksp object suitable for extreme eigenvalues' estimation */
179   if (needscaling) {
180     KSP check_ksp;
181     Vec *workv;
182 
183     ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr);
184     ierr = KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);CHKERRQ(ierr);
185     ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr);
186     ierr = KSPSetTolerances(check_ksp,1.e-4,1.e-12,PETSC_DEFAULT,basis_dofs);CHKERRQ(ierr);
187     ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr);
188     ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr);
189     if (issym) {
190       ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr);
191     }
192     ierr = VecDuplicateVecs(shell_ctx->work_full_1,2,&workv);CHKERRQ(ierr);
193     ierr = KSPSetPC(check_ksp,newpc);CHKERRQ(ierr);
194     ierr = KSPSetUp(check_ksp);CHKERRQ(ierr);
195     ierr = VecSetRandom(workv[1],NULL);CHKERRQ(ierr);
196     ierr = MatMult(local_mat,workv[1],workv[0]);CHKERRQ(ierr);
197     ierr = KSPSolve(check_ksp,workv[0],workv[0]);CHKERRQ(ierr);
198     ierr = VecAXPY(workv[0],-1.,workv[1]);CHKERRQ(ierr);
199     ierr = VecNorm(workv[0],NORM_INFINITY,&test_err);CHKERRQ(ierr);
200     ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr);
201     ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr);
202     if (pcbddc->dbg_flag) {
203       if (isdir) {
204         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);
205       } else {
206         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);
207       }
208     }
209     shell_ctx->scale = 1./lambda_max;
210     ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr);
211     ierr = VecDestroyVecs(2,&workv);CHKERRQ(ierr);
212   }
213   ierr = PCDestroy(&newpc);CHKERRQ(ierr);
214   PetscFunctionReturn(0);
215 }
216 
217 PetscErrorCode PCBDDCNullSpaceCheckCorrection(PC pc, PetscBool isdir)
218 {
219   PC_BDDC                  *pcbddc = (PC_BDDC*)pc->data;
220   Mat_IS                   *matis = (Mat_IS*)pc->pmat->data;
221   KSP                      check_ksp,local_ksp;
222   MatNullSpace             NullSpace = NULL;
223   NullSpaceCorrection_ctx  shell_ctx;
224   PC                       check_pc;
225   Mat                      test_mat,local_mat;
226   PetscReal                test_err,lambda_min,lambda_max;
227   PetscBool                setsym,issym=PETSC_FALSE;
228   Vec                      work1,work2,work3;
229   PetscInt                 k;
230   PetscErrorCode           ierr;
231 
232   PetscFunctionBegin;
233   ierr = MatGetNullSpace(matis->A,&NullSpace);CHKERRQ(ierr);
234   if (!NullSpace) PetscFunctionReturn(0);
235   if (!pcbddc->dbg_flag) PetscFunctionReturn(0);
236   if (isdir) local_ksp = pcbddc->ksp_D;
237   else local_ksp = pcbddc->ksp_R;
238   ierr = KSPGetOperators(local_ksp,&local_mat,NULL);CHKERRQ(ierr);
239   ierr = KSPGetPC(local_ksp,&check_pc);CHKERRQ(ierr);
240   ierr = PCShellGetContext(check_pc,(void**)&shell_ctx);CHKERRQ(ierr);
241   ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr);
242   ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr);
243   ierr = VecDuplicate(shell_ctx->work_full_1,&work3);CHKERRQ(ierr);
244 
245   ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr);
246   ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr);
247   ierr = VecCopy(work1,work2);CHKERRQ(ierr);
248   ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr);
249   ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr);
250   ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr);
251   ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr);
252   if (isdir) {
253     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
254   } else {
255     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
256   }
257 
258   ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr);
259   ierr = MatShift(test_mat,1.);CHKERRQ(ierr);
260   ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr);
261   ierr = MatDestroy(&test_mat);CHKERRQ(ierr);
262   if (isdir) {
263     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
264   } else {
265     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
266   }
267 
268   /* Create ksp object suitable for extreme eigenvalues' estimation */
269   ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr);
270   ierr = KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);CHKERRQ(ierr);
271   ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr);
272   ierr = KSPSetTolerances(check_ksp,1.e-10,1.e-8,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
273   ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr);
274   ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr);
275   if (issym) {
276     ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr);
277   }
278   ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr);
279   ierr = KSPSetUp(check_ksp);CHKERRQ(ierr);
280   ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr);
281   ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr);
282   ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr);
283   ierr = VecAXPY(work2,-1.,work1);CHKERRQ(ierr);
284   ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr);
285   ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr);
286   ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr);
287   if (isdir) {
288     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet adapted KSP (scale %d) %1.14e (it %d, eigs %1.6e %1.6e)\n",PetscGlobalRank,shell_ctx->apply_scaling,test_err,k,lambda_min,lambda_max);CHKERRQ(ierr);
289   } else {
290     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann adapted KSP (scale %d) %1.14e (it %d, eigs %1.6e %1.6e)\n",PetscGlobalRank,shell_ctx->apply_scaling,test_err,k,lambda_min,lambda_max);CHKERRQ(ierr);
291   }
292   ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr);
293   ierr = VecDestroy(&work1);CHKERRQ(ierr);
294   ierr = VecDestroy(&work2);CHKERRQ(ierr);
295   ierr = VecDestroy(&work3);CHKERRQ(ierr);
296   PetscFunctionReturn(0);
297 }
298