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