xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 0369aaf70b351ebf2b584f83653ca0d29b2f88d4)
1 /* TODOLIST
2 
3    Solvers
4    - Add support for cholesky for coarse solver (similar to local solvers)
5    - Propagate ksp prefixes for solvers to mat objects?
6 
7    User interface
8    - ** DM attached to pc?
9 
10    Debugging output
11    - * Better management of verbosity levels of debugging output
12 
13    Extra
14    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
15    - BDDC with MG framework?
16 
17    FETIDP
18    - Move FETIDP code to its own classes
19 
20    MATIS related operations contained in BDDC code
21    - Provide general case for subassembling
22 
23 */
24 
25 #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
26 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
27 #include <petscblaslapack.h>
28 
29 /* temporarily declare it */
30 PetscErrorCode PCApply_BDDC(PC,Vec,Vec);
31 
32 /* -------------------------------------------------------------------------- */
33 #undef __FUNCT__
34 #define __FUNCT__ "PCSetFromOptions_BDDC"
35 PetscErrorCode PCSetFromOptions_BDDC(PetscOptions *PetscOptionsObject,PC pc)
36 {
37   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
38   PetscErrorCode ierr;
39 
40   PetscFunctionBegin;
41   ierr = PetscOptionsHead(PetscOptionsObject,"BDDC options");CHKERRQ(ierr);
42   /* Verbose debugging */
43   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
44   /* Primal space cumstomization */
45   ierr = PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);CHKERRQ(ierr);
46   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
47   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
48   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
49   ierr = PetscOptionsBool("-pc_bddc_use_true_nnsp","Use near null space attached to the matrix without modifications","none",pcbddc->use_nnsp_true,&pcbddc->use_nnsp_true,NULL);CHKERRQ(ierr);
50   ierr = PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is used when multiple constraints are present)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);CHKERRQ(ierr);
51   /* Change of basis */
52   ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not internal change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr);
53   ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not internal change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr);
54   if (!pcbddc->use_change_of_basis) {
55     pcbddc->use_change_on_faces = PETSC_FALSE;
56   }
57   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
58   ierr = PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);CHKERRQ(ierr);
59   ierr = PetscOptionsInt("-pc_bddc_coarse_redistribute","Number of procs where to redistribute coarse problem","none",pcbddc->redistribute_coarse,&pcbddc->redistribute_coarse,NULL);CHKERRQ(ierr);
60   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
61   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
62   ierr = PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);CHKERRQ(ierr);
63   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
64   ierr = PetscOptionsBool("-pc_bddc_schur_rebuild","Whether or not the interface graph for Schur principal minors has to be rebuilt (i.e. define the interface without any adjacency)","none",pcbddc->sub_schurs_rebuild,&pcbddc->sub_schurs_rebuild,NULL);CHKERRQ(ierr);
65   ierr = PetscOptionsInt("-pc_bddc_schur_layers","Number of dofs' layers for the computation of principal minors (i.e. -1 uses all dofs)","none",pcbddc->sub_schurs_layers,&pcbddc->sub_schurs_layers,NULL);CHKERRQ(ierr);
66   ierr = PetscOptionsBool("-pc_bddc_schur_use_useradj","Whether or not the CSR graph specified by the user should be used for computing successive layers (default is to use adj of local mat)","none",pcbddc->sub_schurs_use_useradj,&pcbddc->sub_schurs_use_useradj,NULL);CHKERRQ(ierr);
67   ierr = PetscOptionsBool("-pc_bddc_deluxe_faster","Faster application of deluxe scaling (requires extra work during setup)","none",pcbddc->faster_deluxe,&pcbddc->faster_deluxe,NULL);CHKERRQ(ierr);
68   ierr = PetscOptionsReal("-pc_bddc_adaptive_threshold","Threshold to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&pcbddc->adaptive_threshold,NULL);CHKERRQ(ierr);
69   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);CHKERRQ(ierr);
70   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);CHKERRQ(ierr);
71   ierr = PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);CHKERRQ(ierr);
72   ierr = PetscOptionsInt("-pc_bddc_coarse_adj","Number of processors where to map the coarse adjacency list","none",pcbddc->coarse_adj_red,&pcbddc->coarse_adj_red,NULL);CHKERRQ(ierr);
73   ierr = PetscOptionsBool("-pc_bddc_benign_trick","Apply the benign subspace trick to a class of saddle point problems","none",pcbddc->benign_saddle_point,&pcbddc->benign_saddle_point,NULL);CHKERRQ(ierr);
74   ierr = PetscOptionsTail();CHKERRQ(ierr);
75   PetscFunctionReturn(0);
76 }
77 /* -------------------------------------------------------------------------- */
78 #undef __FUNCT__
79 #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
80 static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
81 {
82   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
83   PetscErrorCode ierr;
84 
85   PetscFunctionBegin;
86   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
87   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
88   pcbddc->user_ChangeOfBasisMatrix = change;
89   PetscFunctionReturn(0);
90 }
91 #undef __FUNCT__
92 #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
93 /*@
94  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
95 
96    Collective on PC
97 
98    Input Parameters:
99 +  pc - the preconditioning context
100 -  change - the change of basis matrix
101 
102    Level: intermediate
103 
104    Notes:
105 
106 .seealso: PCBDDC
107 @*/
108 PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
109 {
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
114   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
115   PetscCheckSameComm(pc,1,change,2);
116   if (pc->mat) {
117     PetscInt rows_c,cols_c,rows,cols;
118     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
119     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
120     if (rows_c != rows) {
121       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
122     }
123     if (cols_c != cols) {
124       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
125     }
126     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
127     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
128     if (rows_c != rows) {
129       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
130     }
131     if (cols_c != cols) {
132       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
133     }
134   }
135   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
136   PetscFunctionReturn(0);
137 }
138 /* -------------------------------------------------------------------------- */
139 #undef __FUNCT__
140 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
141 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
142 {
143   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
144   PetscErrorCode ierr;
145 
146   PetscFunctionBegin;
147   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
148   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
149   pcbddc->user_primal_vertices = PrimalVertices;
150   PetscFunctionReturn(0);
151 }
152 #undef __FUNCT__
153 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
154 /*@
155  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
156 
157    Collective
158 
159    Input Parameters:
160 +  pc - the preconditioning context
161 -  PrimalVertices - index set of primal vertices in local numbering (can be empty)
162 
163    Level: intermediate
164 
165    Notes:
166 
167 .seealso: PCBDDC
168 @*/
169 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
170 {
171   PetscErrorCode ierr;
172 
173   PetscFunctionBegin;
174   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
175   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
176   PetscCheckSameComm(pc,1,PrimalVertices,2);
177   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
178   PetscFunctionReturn(0);
179 }
180 /* -------------------------------------------------------------------------- */
181 #undef __FUNCT__
182 #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
183 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
184 {
185   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
186 
187   PetscFunctionBegin;
188   pcbddc->coarsening_ratio = k;
189   PetscFunctionReturn(0);
190 }
191 
192 #undef __FUNCT__
193 #define __FUNCT__ "PCBDDCSetCoarseningRatio"
194 /*@
195  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
196 
197    Logically collective on PC
198 
199    Input Parameters:
200 +  pc - the preconditioning context
201 -  k - coarsening ratio (H/h at the coarser level)
202 
203    Options Database Keys:
204 .    -pc_bddc_coarsening_ratio
205 
206    Level: intermediate
207 
208    Notes:
209      Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
210 
211 .seealso: PCBDDC, PCBDDCSetLevels()
212 @*/
213 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
214 {
215   PetscErrorCode ierr;
216 
217   PetscFunctionBegin;
218   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
219   PetscValidLogicalCollectiveInt(pc,k,2);
220   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
221   PetscFunctionReturn(0);
222 }
223 
224 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
225 #undef __FUNCT__
226 #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
227 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
228 {
229   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
230 
231   PetscFunctionBegin;
232   pcbddc->use_exact_dirichlet_trick = flg;
233   PetscFunctionReturn(0);
234 }
235 
236 #undef __FUNCT__
237 #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
238 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
239 {
240   PetscErrorCode ierr;
241 
242   PetscFunctionBegin;
243   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
244   PetscValidLogicalCollectiveBool(pc,flg,2);
245   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
246   PetscFunctionReturn(0);
247 }
248 
249 #undef __FUNCT__
250 #define __FUNCT__ "PCBDDCSetLevel_BDDC"
251 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
252 {
253   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
254 
255   PetscFunctionBegin;
256   pcbddc->current_level = level;
257   PetscFunctionReturn(0);
258 }
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "PCBDDCSetLevel"
262 PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
263 {
264   PetscErrorCode ierr;
265 
266   PetscFunctionBegin;
267   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
268   PetscValidLogicalCollectiveInt(pc,level,2);
269   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
270   PetscFunctionReturn(0);
271 }
272 
273 #undef __FUNCT__
274 #define __FUNCT__ "PCBDDCSetLevels_BDDC"
275 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
276 {
277   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
278 
279   PetscFunctionBegin;
280   pcbddc->max_levels = levels;
281   PetscFunctionReturn(0);
282 }
283 
284 #undef __FUNCT__
285 #define __FUNCT__ "PCBDDCSetLevels"
286 /*@
287  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
288 
289    Logically collective on PC
290 
291    Input Parameters:
292 +  pc - the preconditioning context
293 -  levels - the maximum number of levels (max 9)
294 
295    Options Database Keys:
296 .    -pc_bddc_levels
297 
298    Level: intermediate
299 
300    Notes:
301      Default value is 0, i.e. traditional one-level BDDC
302 
303 .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
304 @*/
305 PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
306 {
307   PetscErrorCode ierr;
308 
309   PetscFunctionBegin;
310   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
311   PetscValidLogicalCollectiveInt(pc,levels,2);
312   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
313   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
314   PetscFunctionReturn(0);
315 }
316 /* -------------------------------------------------------------------------- */
317 
318 #undef __FUNCT__
319 #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
320 static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
321 {
322   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
323   PetscErrorCode ierr;
324 
325   PetscFunctionBegin;
326   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
327   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
328   pcbddc->NullSpace = NullSpace;
329   PetscFunctionReturn(0);
330 }
331 
332 #undef __FUNCT__
333 #define __FUNCT__ "PCBDDCSetNullSpace"
334 /*@
335  PCBDDCSetNullSpace - Set nullspace for BDDC operator
336 
337    Logically collective on PC and MatNullSpace
338 
339    Input Parameters:
340 +  pc - the preconditioning context
341 -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
342 
343    Level: intermediate
344 
345    Notes:
346 
347 .seealso: PCBDDC
348 @*/
349 PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
350 {
351   PetscErrorCode ierr;
352 
353   PetscFunctionBegin;
354   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
355   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
356   PetscCheckSameComm(pc,1,NullSpace,2);
357   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
358   PetscFunctionReturn(0);
359 }
360 /* -------------------------------------------------------------------------- */
361 
362 #undef __FUNCT__
363 #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
364 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
365 {
366   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
367   PetscErrorCode ierr;
368 
369   PetscFunctionBegin;
370   /* last user setting takes precendence -> destroy any other customization */
371   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
372   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
373   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
374   pcbddc->DirichletBoundaries = DirichletBoundaries;
375   pcbddc->recompute_topography = PETSC_TRUE;
376   PetscFunctionReturn(0);
377 }
378 
379 #undef __FUNCT__
380 #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
381 /*@
382  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
383 
384    Collective
385 
386    Input Parameters:
387 +  pc - the preconditioning context
388 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
389 
390    Level: intermediate
391 
392    Notes:
393      Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node
394 
395 .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
396 @*/
397 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
398 {
399   PetscErrorCode ierr;
400 
401   PetscFunctionBegin;
402   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
403   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
404   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
405   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
406   PetscFunctionReturn(0);
407 }
408 /* -------------------------------------------------------------------------- */
409 
410 #undef __FUNCT__
411 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
412 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
413 {
414   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
415   PetscErrorCode ierr;
416 
417   PetscFunctionBegin;
418   /* last user setting takes precendence -> destroy any other customization */
419   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
420   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
421   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
422   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
423   pcbddc->recompute_topography = PETSC_TRUE;
424   PetscFunctionReturn(0);
425 }
426 
427 #undef __FUNCT__
428 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
429 /*@
430  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
431 
432    Collective
433 
434    Input Parameters:
435 +  pc - the preconditioning context
436 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
437 
438    Level: intermediate
439 
440    Notes:
441 
442 .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
443 @*/
444 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
445 {
446   PetscErrorCode ierr;
447 
448   PetscFunctionBegin;
449   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
450   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
451   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
452   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
453   PetscFunctionReturn(0);
454 }
455 /* -------------------------------------------------------------------------- */
456 
457 #undef __FUNCT__
458 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
459 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
460 {
461   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
462   PetscErrorCode ierr;
463 
464   PetscFunctionBegin;
465   /* last user setting takes precendence -> destroy any other customization */
466   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
467   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
468   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
469   pcbddc->NeumannBoundaries = NeumannBoundaries;
470   pcbddc->recompute_topography = PETSC_TRUE;
471   PetscFunctionReturn(0);
472 }
473 
474 #undef __FUNCT__
475 #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
476 /*@
477  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
478 
479    Collective
480 
481    Input Parameters:
482 +  pc - the preconditioning context
483 -  NeumannBoundaries - parallel IS defining the Neumann boundaries
484 
485    Level: intermediate
486 
487    Notes:
488      Any process can list any global node
489 
490 .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
491 @*/
492 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
493 {
494   PetscErrorCode ierr;
495 
496   PetscFunctionBegin;
497   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
498   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
499   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
500   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
501   PetscFunctionReturn(0);
502 }
503 /* -------------------------------------------------------------------------- */
504 
505 #undef __FUNCT__
506 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
507 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
508 {
509   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
510   PetscErrorCode ierr;
511 
512   PetscFunctionBegin;
513   /* last user setting takes precendence -> destroy any other customization */
514   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
515   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
516   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
517   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
518   pcbddc->recompute_topography = PETSC_TRUE;
519   PetscFunctionReturn(0);
520 }
521 
522 #undef __FUNCT__
523 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
524 /*@
525  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
526 
527    Collective
528 
529    Input Parameters:
530 +  pc - the preconditioning context
531 -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
532 
533    Level: intermediate
534 
535    Notes:
536 
537 .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
538 @*/
539 PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
540 {
541   PetscErrorCode ierr;
542 
543   PetscFunctionBegin;
544   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
545   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
546   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
547   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
548   PetscFunctionReturn(0);
549 }
550 /* -------------------------------------------------------------------------- */
551 
552 #undef __FUNCT__
553 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
554 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
555 {
556   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
557 
558   PetscFunctionBegin;
559   *DirichletBoundaries = pcbddc->DirichletBoundaries;
560   PetscFunctionReturn(0);
561 }
562 
563 #undef __FUNCT__
564 #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
565 /*@
566  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
567 
568    Collective
569 
570    Input Parameters:
571 .  pc - the preconditioning context
572 
573    Output Parameters:
574 .  DirichletBoundaries - index set defining the Dirichlet boundaries
575 
576    Level: intermediate
577 
578    Notes:
579      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
580 
581 .seealso: PCBDDC
582 @*/
583 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
584 {
585   PetscErrorCode ierr;
586 
587   PetscFunctionBegin;
588   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
589   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
590   PetscFunctionReturn(0);
591 }
592 /* -------------------------------------------------------------------------- */
593 
594 #undef __FUNCT__
595 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
596 static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
597 {
598   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
599 
600   PetscFunctionBegin;
601   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
602   PetscFunctionReturn(0);
603 }
604 
605 #undef __FUNCT__
606 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
607 /*@
608  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
609 
610    Collective
611 
612    Input Parameters:
613 .  pc - the preconditioning context
614 
615    Output Parameters:
616 .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
617 
618    Level: intermediate
619 
620    Notes:
621      The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetDirichletBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetDirichletBoundaries).
622           In the latter case, the IS will be available after PCSetUp.
623 
624 .seealso: PCBDDC
625 @*/
626 PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
627 {
628   PetscErrorCode ierr;
629 
630   PetscFunctionBegin;
631   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
632   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
633   PetscFunctionReturn(0);
634 }
635 /* -------------------------------------------------------------------------- */
636 
637 #undef __FUNCT__
638 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
639 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
640 {
641   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
642 
643   PetscFunctionBegin;
644   *NeumannBoundaries = pcbddc->NeumannBoundaries;
645   PetscFunctionReturn(0);
646 }
647 
648 #undef __FUNCT__
649 #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
650 /*@
651  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
652 
653    Collective
654 
655    Input Parameters:
656 .  pc - the preconditioning context
657 
658    Output Parameters:
659 .  NeumannBoundaries - index set defining the Neumann boundaries
660 
661    Level: intermediate
662 
663    Notes:
664      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
665 
666 .seealso: PCBDDC
667 @*/
668 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
669 {
670   PetscErrorCode ierr;
671 
672   PetscFunctionBegin;
673   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
674   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
675   PetscFunctionReturn(0);
676 }
677 /* -------------------------------------------------------------------------- */
678 
679 #undef __FUNCT__
680 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
681 static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
682 {
683   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
684 
685   PetscFunctionBegin;
686   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
687   PetscFunctionReturn(0);
688 }
689 
690 #undef __FUNCT__
691 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
692 /*@
693  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
694 
695    Collective
696 
697    Input Parameters:
698 .  pc - the preconditioning context
699 
700    Output Parameters:
701 .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
702 
703    Level: intermediate
704 
705    Notes:
706      The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetNeumannBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetNeumannBoundaries).
707           In the latter case, the IS will be available after PCSetUp.
708 
709 .seealso: PCBDDC
710 @*/
711 PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
712 {
713   PetscErrorCode ierr;
714 
715   PetscFunctionBegin;
716   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
717   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
718   PetscFunctionReturn(0);
719 }
720 /* -------------------------------------------------------------------------- */
721 
722 #undef __FUNCT__
723 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
724 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
725 {
726   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
727   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
728   PetscErrorCode ierr;
729 
730   PetscFunctionBegin;
731   /* free old CSR */
732   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
733   /* TODO: PCBDDCGraphSetAdjacency */
734   /* get CSR into graph structure */
735   if (copymode == PETSC_COPY_VALUES) {
736     ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr);
737     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
738     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
739     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
740   } else if (copymode == PETSC_OWN_POINTER) {
741     mat_graph->xadj = (PetscInt*)xadj;
742     mat_graph->adjncy = (PetscInt*)adjncy;
743   }
744   mat_graph->nvtxs_csr = nvtxs;
745   PetscFunctionReturn(0);
746 }
747 
748 #undef __FUNCT__
749 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
750 /*@
751  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local matrix
752 
753    Not collective
754 
755    Input Parameters:
756 +  pc - the preconditioning context
757 .  nvtxs - number of local vertices of the graph (i.e., the size of the local problem)
758 .  xadj, adjncy - the CSR graph
759 -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
760 
761    Level: intermediate
762 
763    Notes:
764 
765 .seealso: PCBDDC,PetscCopyMode
766 @*/
767 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
768 {
769   void (*f)(void) = 0;
770   PetscErrorCode ierr;
771 
772   PetscFunctionBegin;
773   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
774   PetscValidIntPointer(xadj,3);
775   PetscValidIntPointer(adjncy,4);
776   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
777     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
778   }
779   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
780   /* free arrays if PCBDDC is not the PC type */
781   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
782   if (!f && copymode == PETSC_OWN_POINTER) {
783     ierr = PetscFree(xadj);CHKERRQ(ierr);
784     ierr = PetscFree(adjncy);CHKERRQ(ierr);
785   }
786   PetscFunctionReturn(0);
787 }
788 /* -------------------------------------------------------------------------- */
789 
790 #undef __FUNCT__
791 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
792 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
793 {
794   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
795   PetscInt i;
796   PetscErrorCode ierr;
797 
798   PetscFunctionBegin;
799   /* Destroy ISes if they were already set */
800   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
801     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
802   }
803   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
804   /* last user setting takes precendence -> destroy any other customization */
805   for (i=0;i<pcbddc->n_ISForDofs;i++) {
806     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
807   }
808   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
809   pcbddc->n_ISForDofs = 0;
810   /* allocate space then set */
811   if (n_is) {
812     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
813   }
814   for (i=0;i<n_is;i++) {
815     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
816     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
817   }
818   pcbddc->n_ISForDofsLocal=n_is;
819   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
820   pcbddc->recompute_topography = PETSC_TRUE;
821   PetscFunctionReturn(0);
822 }
823 
824 #undef __FUNCT__
825 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
826 /*@
827  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
828 
829    Collective
830 
831    Input Parameters:
832 +  pc - the preconditioning context
833 .  n_is - number of index sets defining the fields
834 -  ISForDofs - array of IS describing the fields in local ordering
835 
836    Level: intermediate
837 
838    Notes:
839      n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
840 
841 .seealso: PCBDDC
842 @*/
843 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
844 {
845   PetscInt       i;
846   PetscErrorCode ierr;
847 
848   PetscFunctionBegin;
849   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
850   PetscValidLogicalCollectiveInt(pc,n_is,2);
851   for (i=0;i<n_is;i++) {
852     PetscCheckSameComm(pc,1,ISForDofs[i],3);
853     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
854   }
855   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
856   PetscFunctionReturn(0);
857 }
858 /* -------------------------------------------------------------------------- */
859 
860 #undef __FUNCT__
861 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
862 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
863 {
864   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
865   PetscInt i;
866   PetscErrorCode ierr;
867 
868   PetscFunctionBegin;
869   /* Destroy ISes if they were already set */
870   for (i=0;i<pcbddc->n_ISForDofs;i++) {
871     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
872   }
873   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
874   /* last user setting takes precendence -> destroy any other customization */
875   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
876     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
877   }
878   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
879   pcbddc->n_ISForDofsLocal = 0;
880   /* allocate space then set */
881   if (n_is) {
882     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
883   }
884   for (i=0;i<n_is;i++) {
885     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
886     pcbddc->ISForDofs[i]=ISForDofs[i];
887   }
888   pcbddc->n_ISForDofs=n_is;
889   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
890   pcbddc->recompute_topography = PETSC_TRUE;
891   PetscFunctionReturn(0);
892 }
893 
894 #undef __FUNCT__
895 #define __FUNCT__ "PCBDDCSetDofsSplitting"
896 /*@
897  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
898 
899    Collective
900 
901    Input Parameters:
902 +  pc - the preconditioning context
903 .  n_is - number of index sets defining the fields
904 -  ISForDofs - array of IS describing the fields in global ordering
905 
906    Level: intermediate
907 
908    Notes:
909      Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
910 
911 .seealso: PCBDDC
912 @*/
913 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
914 {
915   PetscInt       i;
916   PetscErrorCode ierr;
917 
918   PetscFunctionBegin;
919   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
920   PetscValidLogicalCollectiveInt(pc,n_is,2);
921   for (i=0;i<n_is;i++) {
922     PetscCheckSameComm(pc,1,ISForDofs[i],3);
923     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
924   }
925   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
926   PetscFunctionReturn(0);
927 }
928 
929 /* -------------------------------------------------------------------------- */
930 #undef __FUNCT__
931 #define __FUNCT__ "PCPreSolve_BDDC"
932 /* -------------------------------------------------------------------------- */
933 /*
934    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
935                      guess if a transformation of basis approach has been selected.
936 
937    Input Parameter:
938 +  pc - the preconditioner contex
939 
940    Application Interface Routine: PCPreSolve()
941 
942    Notes:
943      The interface routine PCPreSolve() is not usually called directly by
944    the user, but instead is called by KSPSolve().
945 */
946 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
947 {
948   PetscErrorCode ierr;
949   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
950   PC_IS          *pcis = (PC_IS*)(pc->data);
951   Vec            used_vec;
952   PetscBool      copy_rhs = PETSC_TRUE;
953 
954   PetscFunctionBegin;
955   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
956   if (ksp) {
957     PetscBool iscg;
958     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
959     if (!iscg) {
960       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
961     }
962   }
963   /* Creates parallel work vectors used in presolve */
964   if (!pcbddc->original_rhs) {
965     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
966   }
967   if (!pcbddc->temp_solution) {
968     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
969   }
970 
971   if (x) {
972     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
973     used_vec = x;
974   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
975     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
976     used_vec = pcbddc->temp_solution;
977     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
978   }
979 
980   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
981   if (ksp) {
982     /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
983     ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
984     if (!pcbddc->ksp_guess_nonzero) {
985       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
986     }
987   }
988 
989   pcbddc->rhs_change = PETSC_FALSE;
990   /* Take into account zeroed rows -> change rhs and store solution removed */
991   if (rhs) {
992     IS dirIS = NULL;
993 
994     /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
995     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
996     if (dirIS) {
997       Mat_IS            *matis = (Mat_IS*)pc->pmat->data;
998       PetscInt          dirsize,i,*is_indices;
999       PetscScalar       *array_x;
1000       const PetscScalar *array_diagonal;
1001 
1002       ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
1003       ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
1004       ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1005       ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1006       ierr = VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1007       ierr = VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1008       ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
1009       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
1010       ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
1011       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
1012       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
1013       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
1014       ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
1015       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
1016       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1017       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1018       pcbddc->rhs_change = PETSC_TRUE;
1019       ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
1020     }
1021   }
1022 
1023   /* remove the computed solution or the initial guess from the rhs */
1024   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
1025     /* store the original rhs */
1026     if (copy_rhs) {
1027       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1028       copy_rhs = PETSC_FALSE;
1029     }
1030     pcbddc->rhs_change = PETSC_TRUE;
1031     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
1032     ierr = MatMultAdd(pc->mat,used_vec,rhs,rhs);CHKERRQ(ierr);
1033     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
1034     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
1035     if (ksp) {
1036       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr);
1037     }
1038   }
1039   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1040 
1041   /* When using the benign trick: (TODO: what about FETI-DP?)
1042      - change rhs on pressures (iteration matrix is surely of type MATIS)
1043      - compute initial vector in benign space
1044   */
1045   if (rhs && pcbddc->benign_saddle_point) {
1046     Mat_IS *matis = (Mat_IS*)(pc->mat->data);
1047 
1048     /* store the original rhs */
1049     if (copy_rhs) {
1050       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1051       copy_rhs = PETSC_FALSE;
1052     }
1053 
1054     /* now change (locally) the basis */
1055     ierr = VecScatterBegin(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1056     ierr = VecScatterEnd(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1057     if (pcbddc->benign_change) {
1058       ierr = MatMultTranspose(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
1059       ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1060       ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1061       /* change local iteration matrix */
1062       ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_REUSE_MATRIX,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
1063       ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr);
1064       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1065       pcbddc->benign_original_mat = matis->A;
1066       ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
1067       ierr = PetscObjectReference((PetscObject)pcbddc->local_mat);CHKERRQ(ierr);
1068       matis->A = pcbddc->local_mat;
1069       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1070     } else {
1071       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1072       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1073     }
1074     pcbddc->rhs_change = PETSC_TRUE;
1075 
1076     /* compute u^*_h as in Xuemin Tu's thesis (see Section 4.8.1) */
1077     /* TODO: what about Stokes? */
1078     if (!pcbddc->benign_vec) {
1079       ierr = VecDuplicate(rhs,&pcbddc->benign_vec);CHKERRQ(ierr);
1080     }
1081     ierr = VecSet(pcis->vec1_global,0.);CHKERRQ(ierr);
1082     pcbddc->benign_p0 = 0.;
1083     if (pcbddc->benign_p0_lidx >= 0) {
1084       const PetscScalar *array;
1085 
1086       ierr = VecGetArrayRead(pcis->vec2_N,&array);CHKERRQ(ierr);
1087       pcbddc->benign_p0 = array[pcbddc->benign_p0_lidx];
1088       ierr = VecRestoreArrayRead(pcis->vec2_N,&array);CHKERRQ(ierr);
1089     }
1090     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] benign u* should be %1.4e\n",PetscGlobalRank,pcbddc->benign_p0);CHKERRQ(ierr);
1091     ierr = PCBDDCBenignPopOrPushP0(pc,pcis->vec1_global,PETSC_FALSE);CHKERRQ(ierr);
1092     ierr = PCApply_BDDC(pc,pcis->vec1_global,pcbddc->benign_vec);CHKERRQ(ierr);
1093     pcbddc->benign_p0 = 0.;
1094     ierr = PCBDDCBenignPopOrPushP0(pc,pcbddc->benign_vec,PETSC_FALSE);CHKERRQ(ierr);
1095     if (pcbddc->benign_saddle_point) {
1096       Mat_IS* matis = (Mat_IS*)(pc->mat->data);
1097       ierr = VecScatterBegin(matis->rctx,pcbddc->benign_vec,matis->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1098       ierr = VecScatterEnd(matis->rctx,pcbddc->benign_vec,matis->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1099 
1100       if (pcbddc->benign_p0_lidx >=0) {
1101         Mat B0;
1102         Vec dummy_vec;
1103         PetscScalar *array;
1104         PetscInt ii[2];
1105 
1106         ii[0] = 0;
1107         ii[1] = pcbddc->B0_ncol;
1108         ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,1,pcis->n,ii,pcbddc->B0_cols,pcbddc->B0_vals,&B0);CHKERRQ(ierr);
1109         ierr = MatCreateVecs(B0,NULL,&dummy_vec);CHKERRQ(ierr);
1110         ierr = MatMult(B0,matis->x,dummy_vec);CHKERRQ(ierr);
1111         ierr = VecGetArray(dummy_vec,&array);CHKERRQ(ierr);
1112         ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] benign u* is %1.4e\n",PetscGlobalRank,array[0]);CHKERRQ(ierr);
1113         ierr = VecRestoreArray(dummy_vec,&array);CHKERRQ(ierr);
1114         ierr = MatDestroy(&B0);CHKERRQ(ierr);
1115         ierr = VecDestroy(&dummy_vec);CHKERRQ(ierr);
1116       }
1117     }
1118   }
1119 
1120   /* change rhs and iteration matrix if using the change of basis */
1121   if (pcbddc->ChangeOfBasisMatrix) {
1122     PCBDDCChange_ctx change_ctx;
1123 
1124     /* get change ctx */
1125     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1126 
1127     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1128     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1129     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1130     change_ctx->original_mat = pc->mat;
1131 
1132     /* change iteration matrix */
1133     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1134     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1135     pc->mat = pcbddc->new_global_mat;
1136 
1137     /* store the original rhs */
1138     if (copy_rhs) {
1139       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1140       copy_rhs = PETSC_FALSE;
1141     }
1142 
1143     /* change rhs */
1144     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1145     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
1146     pcbddc->rhs_change = PETSC_TRUE;
1147   }
1148 
1149   /* remove non-benign solution from the rhs */
1150   if (pcbddc->benign_saddle_point) {
1151     /* store the original rhs */
1152     if (copy_rhs) {
1153       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1154       copy_rhs = PETSC_FALSE;
1155     }
1156     ierr = VecScale(pcbddc->benign_vec,-1.0);CHKERRQ(ierr);
1157     ierr = MatMultAdd(pc->mat,pcbddc->benign_vec,rhs,rhs);CHKERRQ(ierr);
1158     ierr = VecScale(pcbddc->benign_vec,-1.0);CHKERRQ(ierr);
1159     pcbddc->rhs_change = PETSC_TRUE;
1160   }
1161 
1162   /* set initial guess if using PCG */
1163   if (x && pcbddc->use_exact_dirichlet_trick) {
1164     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1165     ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1166     ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1167     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1168     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1169     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1170     if (ksp) {
1171       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1172     }
1173   }
1174 
1175   /* remove nullspace if present */
1176   if (ksp && x && pcbddc->NullSpace) {
1177     ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr);
1178     /* store the original rhs */
1179     if (copy_rhs) {
1180       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1181       copy_rhs = PETSC_FALSE;
1182     }
1183     pcbddc->rhs_change = PETSC_TRUE;
1184     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1185   }
1186   PetscFunctionReturn(0);
1187 }
1188 
1189 /* -------------------------------------------------------------------------- */
1190 #undef __FUNCT__
1191 #define __FUNCT__ "PCPostSolve_BDDC"
1192 /* -------------------------------------------------------------------------- */
1193 /*
1194    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1195                      approach has been selected. Also, restores rhs to its original state.
1196 
1197    Input Parameter:
1198 +  pc - the preconditioner contex
1199 
1200    Application Interface Routine: PCPostSolve()
1201 
1202    Notes:
1203      The interface routine PCPostSolve() is not usually called directly by
1204      the user, but instead is called by KSPSolve().
1205 */
1206 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1207 {
1208   PetscErrorCode ierr;
1209   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1210 
1211   PetscFunctionBegin;
1212   if (pcbddc->ChangeOfBasisMatrix) {
1213     PCBDDCChange_ctx change_ctx;
1214 
1215     /* get change ctx */
1216     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1217 
1218     /* restore iteration matrix */
1219     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1220     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1221     pc->mat = change_ctx->original_mat;
1222 
1223     /* get solution in original basis */
1224     if (x) {
1225       PC_IS *pcis = (PC_IS*)(pc->data);
1226 
1227       /* restore solution on pressures */
1228       if (pcbddc->benign_saddle_point) {
1229         Mat_IS *matis = (Mat_IS*)pc->mat->data;
1230 
1231         /* add non-benign solution */
1232         ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
1233 
1234         /* change basis on pressures for x */
1235         ierr = VecScatterBegin(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1236         ierr = VecScatterEnd(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1237         if (pcbddc->benign_change) {
1238 
1239           ierr = MatMult(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
1240           ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1241           ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1242         } else {
1243           ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1244           ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1245         }
1246       }
1247       /* change basis on x */
1248       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1249       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
1250     }
1251   }
1252 
1253   /* add solution removed in presolve */
1254   if (x && pcbddc->rhs_change) {
1255     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
1256   }
1257 
1258   /* restore rhs to its original state */
1259   if (rhs && pcbddc->rhs_change) {
1260     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1261   }
1262   pcbddc->rhs_change = PETSC_FALSE;
1263 
1264   /* restore ksp guess state */
1265   if (ksp) {
1266     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
1267   }
1268   PetscFunctionReturn(0);
1269 }
1270 /* -------------------------------------------------------------------------- */
1271 #undef __FUNCT__
1272 #define __FUNCT__ "PCSetUp_BDDC"
1273 /* -------------------------------------------------------------------------- */
1274 /*
1275    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1276                   by setting data structures and options.
1277 
1278    Input Parameter:
1279 +  pc - the preconditioner context
1280 
1281    Application Interface Routine: PCSetUp()
1282 
1283    Notes:
1284      The interface routine PCSetUp() is not usually called directly by
1285      the user, but instead is called by PCApply() if necessary.
1286 */
1287 PetscErrorCode PCSetUp_BDDC(PC pc)
1288 {
1289   PetscErrorCode ierr;
1290   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
1291   Mat_IS*        matis;
1292   MatNullSpace   nearnullspace;
1293   IS             zerodiag = NULL;
1294   PetscInt       nrows,ncols;
1295   PetscBool      computetopography,computesolvers,computesubschurs;
1296   PetscBool      computeconstraintsmatrix;
1297   PetscBool      new_nearnullspace_provided,ismatis;
1298 
1299   PetscFunctionBegin;
1300   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr);
1301   if (!ismatis) {
1302     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
1303   }
1304   ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr);
1305   if (nrows != ncols) {
1306     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
1307   }
1308   matis = (Mat_IS*)pc->pmat->data;
1309   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1310   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1311      Also, BDDC directly build the Dirichlet problem */
1312   /* split work */
1313   if (pc->setupcalled) {
1314     if (pc->flag == SAME_NONZERO_PATTERN) {
1315       computetopography = PETSC_FALSE;
1316       computesolvers = PETSC_TRUE;
1317     } else { /* DIFFERENT_NONZERO_PATTERN */
1318       computetopography = PETSC_TRUE;
1319       computesolvers = PETSC_TRUE;
1320     }
1321   } else {
1322     computetopography = PETSC_TRUE;
1323     computesolvers = PETSC_TRUE;
1324   }
1325   if (pcbddc->recompute_topography) {
1326     computetopography = PETSC_TRUE;
1327   }
1328   computeconstraintsmatrix = PETSC_FALSE;
1329   if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) {
1330     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling");
1331   }
1332   pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling);
1333   if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;
1334 
1335   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
1336   if (pcbddc->faster_deluxe && pcbddc->adaptive_selection && pcbddc->use_change_of_basis) {
1337     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute faster deluxe if adaptivity and change of basis are both requested. Rerun with -pc_bddc_deluxe_faster false");
1338   }
1339 
1340   /* check if the iteration matrix is of type MATIS in case the benign trick has been requested */
1341   ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&ismatis);CHKERRQ(ierr);
1342   if (pcbddc->benign_saddle_point && !ismatis) {
1343     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner with benign subspace trick requires the iteration matrix to be of type MATIS");
1344   }
1345 
1346   /* Get stdout for dbg */
1347   if (pcbddc->dbg_flag) {
1348     if (!pcbddc->dbg_viewer) {
1349       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1350       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1351     }
1352     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1353   }
1354 
1355   if (pcbddc->user_ChangeOfBasisMatrix) {
1356     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1357     pcbddc->use_change_of_basis = PETSC_FALSE;
1358     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
1359   } else {
1360     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1361     if (!pcbddc->benign_saddle_point) {
1362       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1363       pcbddc->local_mat = matis->A;
1364     } else { /* TODO: handle user change of basis */
1365       PetscInt  nz;
1366       PetscBool sorted;
1367 
1368       ierr = MatFindZeroDiagonals(matis->A,&zerodiag);CHKERRQ(ierr);
1369       ierr = ISSorted(zerodiag,&sorted);CHKERRQ(ierr);
1370       if (!sorted) {
1371         ierr = ISSort(zerodiag);CHKERRQ(ierr);
1372       }
1373       ierr = ISGetLocalSize(zerodiag,&nz);CHKERRQ(ierr);
1374       if (nz) {
1375         IS                zerodiagc;
1376         PetscScalar       *array;
1377         const PetscInt    *idxs,*idxsc;
1378         PetscInt          i,n,*nnz;
1379 
1380         /* TODO: add check for shared dofs */
1381         pcbddc->use_change_of_basis = PETSC_TRUE;
1382         pcbddc->use_change_on_faces = PETSC_TRUE;
1383         ierr = MatGetLocalSize(matis->A,&n,NULL);CHKERRQ(ierr);
1384         ierr = ISComplement(zerodiag,0,n,&zerodiagc);CHKERRQ(ierr);
1385         ierr = ISGetIndices(zerodiag,&idxs);CHKERRQ(ierr);
1386         ierr = ISGetIndices(zerodiagc,&idxsc);CHKERRQ(ierr);
1387         /* local change of basis for pressures */
1388         ierr = MatDestroy(&pcbddc->benign_change);CHKERRQ(ierr);
1389         ierr = MatCreate(PetscObjectComm((PetscObject)matis->A),&pcbddc->benign_change);CHKERRQ(ierr);
1390         ierr = MatSetType(pcbddc->benign_change,MATAIJ);CHKERRQ(ierr);
1391         ierr = MatSetSizes(pcbddc->benign_change,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
1392         ierr = PetscMalloc1(n,&nnz);CHKERRQ(ierr);
1393         for (i=0;i<n-nz;i++) nnz[idxsc[i]] = 1; /* identity on velocities */
1394         for (i=0;i<nz-1;i++) nnz[idxs[i]] = 2; /* change on pressures */
1395         nnz[idxs[nz-1]] = nz; /* last local pressure dof: _0 set */
1396         ierr = MatSeqAIJSetPreallocation(pcbddc->benign_change,0,nnz);CHKERRQ(ierr);
1397         ierr = PetscFree(nnz);CHKERRQ(ierr);
1398         /* set identity on velocities */
1399         for (i=0;i<n-nz;i++) {
1400           ierr = MatSetValue(pcbddc->benign_change,idxsc[i],idxsc[i],1.,INSERT_VALUES);CHKERRQ(ierr);
1401         }
1402         /* set change on pressures */
1403         for (i=0;i<nz-1;i++) {
1404           PetscScalar vals[2];
1405           PetscInt    cols[2];
1406 
1407           /* TODO: add quadrature */
1408           cols[0] = idxs[i];
1409           cols[1] = idxs[nz-1];
1410           vals[0] = 1.;
1411           vals[1] = 1./nz;
1412           ierr = MatSetValues(pcbddc->benign_change,1,idxs+i,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
1413         }
1414         ierr = PetscMalloc1(nz,&array);CHKERRQ(ierr);
1415         for (i=0;i<nz-1;i++) array[i] = -1.;
1416         array[nz-1] = 1./nz;
1417         ierr = MatSetValues(pcbddc->benign_change,1,idxs+nz-1,nz,idxs,array,INSERT_VALUES);CHKERRQ(ierr);
1418         ierr = PetscFree(array);CHKERRQ(ierr);
1419         ierr = MatAssemblyBegin(pcbddc->benign_change,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1420         ierr = MatAssemblyEnd(pcbddc->benign_change,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1421         /* TODO: need optimization? */
1422         ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_INITIAL_MATRIX,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
1423         /* store local and global idxs for p0 */
1424         pcbddc->benign_p0_lidx = idxs[nz-1];
1425         ierr = ISLocalToGlobalMappingApply(pc->pmat->rmap->mapping,1,&idxs[nz-1],&pcbddc->benign_p0_gidx);CHKERRQ(ierr);
1426         ierr = ISRestoreIndices(zerodiag,&idxs);CHKERRQ(ierr);
1427         ierr = ISRestoreIndices(zerodiagc,&idxsc);CHKERRQ(ierr);
1428         ierr = ISDestroy(&zerodiagc);CHKERRQ(ierr);
1429         /* pop B0 mat from pcbddc->local_mat */
1430         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
1431       } else { /* this is unlikely to happen but, just in case, destroy the empty IS */
1432         pcbddc->benign_p0_lidx = -1;
1433         pcbddc->benign_p0_gidx = -1;
1434         ierr = ISDestroy(&zerodiag);CHKERRQ(ierr);
1435         ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1436         pcbddc->local_mat = matis->A;
1437       }
1438     }
1439   }
1440 
1441   /* workaround for reals */
1442 #if !defined(PETSC_USE_COMPLEX)
1443   if (matis->A->symmetric_set) {
1444     ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr);
1445   }
1446 #endif
1447   if (matis->A->symmetric_set) {
1448     ierr = MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);CHKERRQ(ierr);
1449   }
1450   if (matis->A->spd_set) {
1451     ierr = MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);CHKERRQ(ierr);
1452   }
1453 
1454   /* Set up all the "iterative substructuring" common block without computing solvers */
1455   {
1456     Mat temp_mat;
1457 
1458     temp_mat = matis->A;
1459     matis->A = pcbddc->local_mat;
1460     ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr);
1461     pcbddc->local_mat = matis->A;
1462     matis->A = temp_mat;
1463   }
1464 
1465   /* Analyze interface */
1466   if (computetopography) {
1467     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1468     computeconstraintsmatrix = PETSC_TRUE;
1469   }
1470 
1471   /* check b(v_I,p_0) = 0 for all v_I */
1472   if (pcbddc->dbg_flag && zerodiag) {
1473     PC_IS          *pcis = (PC_IS*)(pc->data);
1474     IS             dirIS = NULL;
1475     PetscScalar    *vals;
1476     const PetscInt *idxs;
1477     PetscInt       i,nz;
1478 
1479     /* p0 */
1480     ierr = VecSet(pcis->vec1_N,0.);CHKERRQ(ierr);
1481     ierr = PetscMalloc1(pcis->n,&vals);CHKERRQ(ierr);
1482     ierr = ISGetLocalSize(zerodiag,&nz);CHKERRQ(ierr);
1483     ierr = ISGetIndices(zerodiag,&idxs);CHKERRQ(ierr);
1484     for (i=0;i<nz;i++) vals[i] = 1.;
1485     ierr = VecSetValues(pcis->vec1_N,nz,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1486     ierr = VecAssemblyBegin(pcis->vec1_N);CHKERRQ(ierr);
1487     ierr = VecAssemblyEnd(pcis->vec1_N);CHKERRQ(ierr);
1488     /* v_I */
1489     ierr = VecSetRandom(pcis->vec2_N,NULL);CHKERRQ(ierr);
1490     for (i=0;i<nz;i++) vals[i] = 0.;
1491     ierr = VecSetValues(pcis->vec2_N,nz,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1492     ierr = ISRestoreIndices(zerodiag,&idxs);CHKERRQ(ierr);
1493     ierr = ISGetIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
1494     for (i=0;i<pcis->n_B;i++) vals[i] = 0.;
1495     ierr = VecSetValues(pcis->vec2_N,pcis->n_B,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1496     ierr = ISRestoreIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
1497     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
1498     if (dirIS) {
1499       PetscInt n;
1500 
1501       ierr = ISGetLocalSize(dirIS,&n);CHKERRQ(ierr);
1502       ierr = ISGetIndices(dirIS,&idxs);CHKERRQ(ierr);
1503       for (i=0;i<n;i++) vals[i] = 0.;
1504       ierr = VecSetValues(pcis->vec2_N,n,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1505       ierr = ISRestoreIndices(dirIS,&idxs);CHKERRQ(ierr);
1506     }
1507     ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
1508     ierr = VecAssemblyBegin(pcis->vec2_N);CHKERRQ(ierr);
1509     ierr = VecAssemblyEnd(pcis->vec2_N);CHKERRQ(ierr);
1510     ierr = VecSet(matis->x,0.);CHKERRQ(ierr);
1511     ierr = MatMult(matis->A,pcis->vec1_N,matis->x);CHKERRQ(ierr);
1512     ierr = VecDot(matis->x,pcis->vec2_N,&vals[0]);CHKERRQ(ierr);
1513     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] check original mat: %1.4e\n",PetscGlobalRank,vals[0]);CHKERRQ(ierr);
1514     ierr = PetscFree(vals);CHKERRQ(ierr);
1515   }
1516 
1517   /* check PCBDDCBenignPopOrPush */
1518   if (pcbddc->dbg_flag && pcbddc->benign_saddle_point) {
1519     PC_IS  *pcis = (PC_IS*)(pc->data);
1520 
1521     ierr = VecSetRandom(pcis->vec1_global,NULL);CHKERRQ(ierr);
1522     pcbddc->benign_p0 = -PetscGlobalRank;
1523     ierr = PCBDDCBenignPopOrPushP0(pc,pcis->vec1_global,PETSC_FALSE);CHKERRQ(ierr);
1524     pcbddc->benign_p0 = 1;
1525     ierr = PCBDDCBenignPopOrPushP0(pc,pcis->vec1_global,PETSC_TRUE);CHKERRQ(ierr);
1526     if (pcbddc->benign_p0_gidx >=0 && pcbddc->benign_p0 != -PetscGlobalRank) {
1527       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error testing PCBDDCBenignPopOrPushP0! Found %1.4e instead of %1.4e\n",pcbddc->benign_p0,-PetscGlobalRank);CHKERRQ(ierr);
1528     }
1529   }
1530 
1531   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1532   if (computesolvers) {
1533     PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs;
1534 
1535     if (computesubschurs && computetopography) {
1536       ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr);
1537     }
1538     if (sub_schurs->use_mumps) {
1539       if (computesubschurs) {
1540         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1541       }
1542       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1543     } else {
1544       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1545       if (computesubschurs) {
1546         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1547       }
1548     }
1549     if (pcbddc->adaptive_selection) {
1550       ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr);
1551       computeconstraintsmatrix = PETSC_TRUE;
1552     }
1553   }
1554 
1555   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1556   new_nearnullspace_provided = PETSC_FALSE;
1557   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1558   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1559     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1560       new_nearnullspace_provided = PETSC_TRUE;
1561     } else {
1562       /* determine if the two nullspaces are different (should be lightweight) */
1563       if (nearnullspace != pcbddc->onearnullspace) {
1564         new_nearnullspace_provided = PETSC_TRUE;
1565       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1566         PetscInt         i;
1567         const Vec        *nearnullvecs;
1568         PetscObjectState state;
1569         PetscInt         nnsp_size;
1570         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1571         for (i=0;i<nnsp_size;i++) {
1572           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1573           if (pcbddc->onearnullvecs_state[i] != state) {
1574             new_nearnullspace_provided = PETSC_TRUE;
1575             break;
1576           }
1577         }
1578       }
1579     }
1580   } else {
1581     if (!nearnullspace) { /* both nearnullspaces are null */
1582       new_nearnullspace_provided = PETSC_FALSE;
1583     } else { /* nearnullspace attached later */
1584       new_nearnullspace_provided = PETSC_TRUE;
1585     }
1586   }
1587 
1588   /* Setup constraints and related work vectors */
1589   /* reset primal space flags */
1590   pcbddc->new_primal_space = PETSC_FALSE;
1591   pcbddc->new_primal_space_local = PETSC_FALSE;
1592   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1593     /* It also sets the primal space flags */
1594     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1595     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1596     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
1597   }
1598 
1599   if (computesolvers || pcbddc->new_primal_space) {
1600     if (pcbddc->use_change_of_basis) {
1601       PC_IS *pcis = (PC_IS*)(pc->data);
1602       Mat   temp_mat = NULL;
1603 
1604       if (zerodiag) {
1605         /* insert B0 in pcbddc->local_mat */
1606         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_FALSE);CHKERRQ(ierr);
1607         if (pcbddc->dbg_flag) {
1608           PC_IS          *pcis = (PC_IS*)(pc->data);
1609           IS             dirIS = NULL;
1610           PetscScalar    *vals;
1611           const PetscInt *idxs;
1612           PetscInt       i,nz;
1613 
1614           /* p0 */
1615           ierr = VecSet(pcis->vec1_N,0.);CHKERRQ(ierr);
1616           ierr = PetscMalloc1(pcis->n,&vals);CHKERRQ(ierr);
1617           ierr = ISGetLocalSize(zerodiag,&nz);CHKERRQ(ierr);
1618           ierr = ISGetIndices(zerodiag,&idxs);CHKERRQ(ierr);
1619           vals[0] = 1.;
1620           ierr = VecSetValues(pcis->vec1_N,1,&idxs[nz-1],vals,INSERT_VALUES);CHKERRQ(ierr);
1621           ierr = VecAssemblyBegin(pcis->vec1_N);CHKERRQ(ierr);
1622           ierr = VecAssemblyEnd(pcis->vec1_N);CHKERRQ(ierr);
1623           /* v_I */
1624           ierr = VecSetRandom(pcis->vec2_N,NULL);CHKERRQ(ierr);
1625           for (i=0;i<nz;i++) vals[i] = 0.;
1626           ierr = VecSetValues(pcis->vec2_N,nz,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1627           ierr = ISRestoreIndices(zerodiag,&idxs);CHKERRQ(ierr);
1628           ierr = ISGetIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
1629           for (i=0;i<pcis->n_B;i++) vals[i] = 0.;
1630           ierr = VecSetValues(pcis->vec2_N,pcis->n_B,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1631           ierr = ISRestoreIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
1632           ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
1633           if (dirIS) {
1634             PetscInt n;
1635 
1636             ierr = ISGetLocalSize(dirIS,&n);CHKERRQ(ierr);
1637             ierr = ISGetIndices(dirIS,&idxs);CHKERRQ(ierr);
1638             for (i=0;i<n;i++) vals[i] = 0.;
1639             ierr = VecSetValues(pcis->vec2_N,n,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1640             ierr = ISRestoreIndices(dirIS,&idxs);CHKERRQ(ierr);
1641           }
1642           ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
1643           ierr = VecAssemblyBegin(pcis->vec2_N);CHKERRQ(ierr);
1644           ierr = VecAssemblyEnd(pcis->vec2_N);CHKERRQ(ierr);
1645           ierr = VecSet(matis->x,0.);CHKERRQ(ierr);
1646           ierr = MatMult(pcbddc->local_mat,pcis->vec1_N,matis->x);CHKERRQ(ierr);
1647           ierr = VecDot(matis->x,pcis->vec2_N,&vals[0]);CHKERRQ(ierr);
1648           ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] check new mat: %1.4e\n",PetscGlobalRank,vals[0]);CHKERRQ(ierr);
1649           ierr = PetscFree(vals);CHKERRQ(ierr);
1650         }
1651         /* hack: swap pointers */
1652         temp_mat = matis->A;
1653         matis->A = pcbddc->local_mat;
1654         pcbddc->local_mat = NULL;
1655       }
1656       ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
1657       if (zerodiag) {
1658         /* restore original matrix */
1659         ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
1660         matis->A = temp_mat;
1661         /* pop B0 from pcbddc->local_mat */
1662         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
1663       }
1664       /* get submatrices */
1665       ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
1666       ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
1667       ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
1668       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
1669       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
1670       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
1671       /* set flag in pcis to not reuse submatrices during PCISCreate */
1672       pcis->reusesubmatrices = PETSC_FALSE;
1673     } else if (!pcbddc->user_ChangeOfBasisMatrix) {
1674       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1675       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1676       pcbddc->local_mat = matis->A;
1677     }
1678     /* SetUp coarse and local Neumann solvers */
1679     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1680     /* SetUp Scaling operator */
1681     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
1682   }
1683 
1684   /* check BDDC operator */
1685   if (pcbddc->dbg_flag && (
1686       (pcbddc->n_vertices == pcbddc->local_primal_size) || pcbddc->benign_saddle_point) ) {
1687     PC_IS          *pcis = (PC_IS*)(pc->data);
1688     Mat            S_j,B0=NULL,B0_B=NULL;
1689     Vec            dummy_vec=NULL,vec_check_B,vec_scale_P;
1690     PetscScalar    norm,p0_check,*array,*array2;
1691     PetscInt       i;
1692 
1693     /* B0 and B0_B */
1694     if (zerodiag) {
1695       IS       dummy;
1696       PetscInt ii[2];
1697 
1698       ii[0] = 0;
1699       ii[1] = pcbddc->B0_ncol;
1700       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,1,pcis->n,ii,pcbddc->B0_cols,pcbddc->B0_vals,&B0);CHKERRQ(ierr);
1701       ierr = ISCreateStride(PETSC_COMM_SELF,1,0,1,&dummy);CHKERRQ(ierr);
1702       ierr = MatGetSubMatrix(B0,dummy,pcis->is_B_local,MAT_INITIAL_MATRIX,&B0_B);CHKERRQ(ierr);
1703       ierr = MatCreateVecs(B0_B,NULL,&dummy_vec);CHKERRQ(ierr);
1704       ierr = ISDestroy(&dummy);CHKERRQ(ierr);
1705     }
1706     /* I need a primal vector to scale primal nodes since BDDC sums contibutions */
1707     ierr = VecDuplicate(pcbddc->vec1_P,&vec_scale_P);CHKERRQ(ierr);
1708     ierr = VecSet(pcbddc->vec1_P,1.0);CHKERRQ(ierr);
1709     ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,pcbddc->vec1_P,pcbddc->coarse_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1710     ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,pcbddc->vec1_P,pcbddc->coarse_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1711     ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,pcbddc->coarse_vec,vec_scale_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1712     ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,pcbddc->coarse_vec,vec_scale_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1713     ierr = VecReciprocal(vec_scale_P);CHKERRQ(ierr);
1714     /* S_j */
1715     ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&S_j);CHKERRQ(ierr);
1716     ierr = MatSchurComplementSetKSP(S_j,pcbddc->ksp_D);CHKERRQ(ierr);
1717 
1718     /* mimic vector in \widetilde{W}_\Gamma */
1719     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
1720     /* continuous in primal space */
1721     ierr = VecSetRandom(pcbddc->coarse_vec,NULL);CHKERRQ(ierr);
1722     ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1723     ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1724     ierr = VecGetArray(pcbddc->vec1_P,&array);CHKERRQ(ierr);
1725     if (zerodiag) {
1726       p0_check = array[pcbddc->local_primal_size-1];
1727     } else {
1728       p0_check = 0;
1729     }
1730     ierr = VecSetValues(pcis->vec1_N,pcbddc->local_primal_size,pcbddc->local_primal_ref_node,array,INSERT_VALUES);CHKERRQ(ierr);
1731     ierr = VecRestoreArray(pcbddc->vec1_P,&array);CHKERRQ(ierr);
1732     ierr = VecAssemblyBegin(pcis->vec1_N);CHKERRQ(ierr);
1733     ierr = VecAssemblyEnd(pcis->vec1_N);CHKERRQ(ierr);
1734     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1735     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1736     ierr = VecDuplicate(pcis->vec2_B,&vec_check_B);CHKERRQ(ierr);
1737     ierr = VecCopy(pcis->vec2_B,vec_check_B);CHKERRQ(ierr);
1738 
1739     /* assemble rhs for coarse problem */
1740     /* widetilde{S}_\Gamma w_\Gamma + \widetilde{B0}^T_B p0 */
1741     /* local with Schur */
1742     ierr = MatMult(S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
1743     if (zerodiag) {
1744       ierr = VecGetArray(dummy_vec,&array);CHKERRQ(ierr);
1745       array[0] = p0_check;
1746       ierr = VecRestoreArray(dummy_vec,&array);CHKERRQ(ierr);
1747       ierr = MatMultTransposeAdd(B0_B,dummy_vec,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
1748     }
1749     /* sum on primal nodes the local contributions */
1750     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1751     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1752     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
1753     ierr = VecGetArray(pcbddc->vec1_P,&array2);CHKERRQ(ierr);
1754     for (i=0;i<pcbddc->local_primal_size;i++) array2[i] = array[pcbddc->local_primal_ref_node[i]];
1755     ierr = VecRestoreArray(pcbddc->vec1_P,&array2);CHKERRQ(ierr);
1756     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
1757     ierr = VecSet(pcbddc->coarse_vec,0.);CHKERRQ(ierr);
1758     ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,pcbddc->vec1_P,pcbddc->coarse_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1759     ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,pcbddc->vec1_P,pcbddc->coarse_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1760     ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1761     ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1762     ierr = VecGetArray(pcbddc->vec1_P,&array);CHKERRQ(ierr);
1763     /* scale primal nodes (BDDC sums contibutions) */
1764     ierr = VecPointwiseMult(pcbddc->vec1_P,vec_scale_P,pcbddc->vec1_P);CHKERRQ(ierr);
1765     ierr = VecSetValues(pcis->vec1_N,pcbddc->local_primal_size,pcbddc->local_primal_ref_node,array,INSERT_VALUES);CHKERRQ(ierr);
1766     ierr = VecRestoreArray(pcbddc->vec1_P,&array);CHKERRQ(ierr);
1767     ierr = VecAssemblyBegin(pcis->vec1_N);CHKERRQ(ierr);
1768     ierr = VecAssemblyEnd(pcis->vec1_N);CHKERRQ(ierr);
1769     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1770     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1771     /* global: \widetilde{B0}_B w_\Gamma */
1772     if (zerodiag) {
1773       ierr = MatMult(B0_B,pcis->vec2_B,dummy_vec);CHKERRQ(ierr);
1774       ierr = VecGetArray(dummy_vec,&array);CHKERRQ(ierr);
1775       pcbddc->benign_p0 = array[0];
1776       ierr = VecRestoreArray(dummy_vec,&array);CHKERRQ(ierr);
1777     } else {
1778       pcbddc->benign_p0 = 0.;
1779     }
1780     /* BDDC */
1781     ierr = VecSet(pcis->vec1_D,0.);CHKERRQ(ierr);
1782     ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
1783 
1784     ierr = VecCopy(pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
1785     ierr = VecAXPY(pcis->vec1_B,-1.0,vec_check_B);CHKERRQ(ierr);
1786     ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&norm);CHKERRQ(ierr);
1787     PetscPrintf(PETSC_COMM_SELF,"[%d] BDDC local error is %1.4e\n",PetscGlobalRank,norm);
1788     if (pcbddc->benign_p0_lidx >= 0) {
1789       PetscPrintf(PETSC_COMM_SELF,"[%d] BDDC p0 error is %1.4e\n",PetscGlobalRank,PetscAbsScalar(pcbddc->benign_p0-p0_check));
1790     }
1791 
1792     ierr = VecDestroy(&vec_scale_P);CHKERRQ(ierr);
1793     ierr = VecDestroy(&vec_check_B);CHKERRQ(ierr);
1794     ierr = VecDestroy(&dummy_vec);CHKERRQ(ierr);
1795     ierr = MatDestroy(&S_j);CHKERRQ(ierr);
1796     ierr = MatDestroy(&B0);CHKERRQ(ierr);
1797     ierr = MatDestroy(&B0_B);CHKERRQ(ierr);
1798   }
1799   ierr = ISDestroy(&zerodiag);CHKERRQ(ierr);
1800 
1801   if (pcbddc->dbg_flag) {
1802     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1803   }
1804   PetscFunctionReturn(0);
1805 }
1806 
1807 /* -------------------------------------------------------------------------- */
1808 /*
1809    PCApply_BDDC - Applies the BDDC operator to a vector.
1810 
1811    Input Parameters:
1812 +  pc - the preconditioner context
1813 -  r - input vector (global)
1814 
1815    Output Parameter:
1816 .  z - output vector (global)
1817 
1818    Application Interface Routine: PCApply()
1819  */
1820 #undef __FUNCT__
1821 #define __FUNCT__ "PCApply_BDDC"
1822 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1823 {
1824   PC_IS             *pcis = (PC_IS*)(pc->data);
1825   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1826   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
1827   PetscErrorCode    ierr;
1828   const PetscScalar one = 1.0;
1829   const PetscScalar m_one = -1.0;
1830   const PetscScalar zero = 0.0;
1831 
1832 /* This code is similar to that provided in nn.c for PCNN
1833    NN interface preconditioner changed to BDDC
1834    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
1835 
1836   PetscFunctionBegin;
1837   if (pcbddc->benign_saddle_point) { /* extract p0 from r */
1838     ierr = PCBDDCBenignPopOrPushP0(pc,r,PETSC_TRUE);CHKERRQ(ierr);
1839   }
1840   if (!pcbddc->use_exact_dirichlet_trick) {
1841     ierr = VecCopy(r,z);CHKERRQ(ierr);
1842     /* First Dirichlet solve */
1843     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1844     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1845     /*
1846       Assembling right hand side for BDDC operator
1847       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1848       - pcis->vec1_B the interface part of the global vector z
1849     */
1850     if (n_D) {
1851       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1852       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1853       if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1854       ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1855     } else {
1856       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1857     }
1858     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1859     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1860     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1861   } else {
1862     if (pcbddc->switch_static) {
1863       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1864     }
1865     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1866   }
1867 
1868   /* Apply interface preconditioner
1869      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1870   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
1871 
1872   /* Apply transpose of partition of unity operator */
1873   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1874 
1875   /* Second Dirichlet solve and assembling of output */
1876   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1877   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1878   if (n_B) {
1879     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1880     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1881   } else if (pcbddc->switch_static) {
1882     ierr = MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1883   }
1884   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1885 
1886   if (!pcbddc->use_exact_dirichlet_trick) {
1887     if (pcbddc->switch_static) {
1888       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1889     } else {
1890       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1891     }
1892     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1893     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1894   } else {
1895     if (pcbddc->switch_static) {
1896       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1897     } else {
1898       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1899     }
1900     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1901     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1902   }
1903 
1904   if (pcbddc->benign_saddle_point) { /* push p0 (computed in PCBDDCApplyInterface) */
1905     ierr = PCBDDCBenignPopOrPushP0(pc,z,PETSC_FALSE);CHKERRQ(ierr);
1906   }
1907   PetscFunctionReturn(0);
1908 }
1909 
1910 /* -------------------------------------------------------------------------- */
1911 /*
1912    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
1913 
1914    Input Parameters:
1915 +  pc - the preconditioner context
1916 -  r - input vector (global)
1917 
1918    Output Parameter:
1919 .  z - output vector (global)
1920 
1921    Application Interface Routine: PCApplyTranspose()
1922  */
1923 #undef __FUNCT__
1924 #define __FUNCT__ "PCApplyTranspose_BDDC"
1925 PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
1926 {
1927   PC_IS             *pcis = (PC_IS*)(pc->data);
1928   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1929   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
1930   PetscErrorCode    ierr;
1931   const PetscScalar one = 1.0;
1932   const PetscScalar m_one = -1.0;
1933   const PetscScalar zero = 0.0;
1934 
1935   PetscFunctionBegin;
1936   if (!pcbddc->use_exact_dirichlet_trick) {
1937     ierr = VecCopy(r,z);CHKERRQ(ierr);
1938     /* First Dirichlet solve */
1939     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1940     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1941     /*
1942       Assembling right hand side for BDDC operator
1943       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1944       - pcis->vec1_B the interface part of the global vector z
1945     */
1946     if (n_D) {
1947       ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1948       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1949       if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1950       ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1951     } else {
1952       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1953     }
1954     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1955     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1956     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1957   } else {
1958     if (pcbddc->switch_static) {
1959       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1960     }
1961     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1962   }
1963 
1964   /* Apply interface preconditioner
1965      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1966   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
1967 
1968   /* Apply transpose of partition of unity operator */
1969   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1970 
1971   /* Second Dirichlet solve and assembling of output */
1972   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1973   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1974   if (n_B) {
1975     ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1976     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1977   } else if (pcbddc->switch_static) {
1978     ierr = MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1979   }
1980   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1981   if (!pcbddc->use_exact_dirichlet_trick) {
1982     if (pcbddc->switch_static) {
1983       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1984     } else {
1985       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1986     }
1987     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1988     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1989   } else {
1990     if (pcbddc->switch_static) {
1991       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1992     } else {
1993       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1994     }
1995     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1996     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1997   }
1998   PetscFunctionReturn(0);
1999 }
2000 /* -------------------------------------------------------------------------- */
2001 
2002 #undef __FUNCT__
2003 #define __FUNCT__ "PCDestroy_BDDC"
2004 PetscErrorCode PCDestroy_BDDC(PC pc)
2005 {
2006   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
2007   PetscErrorCode ierr;
2008 
2009   PetscFunctionBegin;
2010   /* free data created by PCIS */
2011   ierr = PCISDestroy(pc);CHKERRQ(ierr);
2012   /* free BDDC custom data  */
2013   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
2014   /* destroy objects related to topography */
2015   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
2016   /* free allocated graph structure */
2017   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
2018   /* free allocated sub schurs structure */
2019   ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr);
2020   /* destroy objects for scaling operator */
2021   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
2022   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
2023   /* free solvers stuff */
2024   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
2025   /* free global vectors needed in presolve */
2026   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
2027   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
2028   /* free stuff for change of basis hooks */
2029   if (pcbddc->new_global_mat) {
2030     PCBDDCChange_ctx change_ctx;
2031     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
2032     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
2033     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
2034     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
2035     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
2036   }
2037   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
2038   /* remove functions */
2039   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
2040   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
2041   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
2042   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
2043   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
2044   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
2045   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
2046   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
2047   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
2048   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
2049   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
2050   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
2051   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
2052   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
2053   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
2054   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
2055   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
2056   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
2057   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
2058   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
2059   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
2060   /* Free the private data structure */
2061   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2062   PetscFunctionReturn(0);
2063 }
2064 /* -------------------------------------------------------------------------- */
2065 
2066 #undef __FUNCT__
2067 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
2068 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
2069 {
2070   FETIDPMat_ctx  mat_ctx;
2071   Vec            copy_standard_rhs;
2072   PC_IS*         pcis;
2073   PC_BDDC*       pcbddc;
2074   PetscErrorCode ierr;
2075 
2076   PetscFunctionBegin;
2077   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
2078   pcis = (PC_IS*)mat_ctx->pc->data;
2079   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
2080 
2081   /*
2082      change of basis for physical rhs if needed
2083      It also changes the rhs in case of dirichlet boundaries
2084      TODO: better management when FETIDP will have its own class
2085   */
2086   ierr = VecDuplicate(standard_rhs,&copy_standard_rhs);CHKERRQ(ierr);
2087   ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr);
2088   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr);
2089   /* store vectors for computation of fetidp final solution */
2090   ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2091   ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2092   /* scale rhs since it should be unassembled */
2093   /* TODO use counter scaling? (also below) */
2094   ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2095   ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2096   /* Apply partition of unity */
2097   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
2098   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
2099   if (!pcbddc->switch_static) {
2100     /* compute partially subassembled Schur complement right-hand side */
2101     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
2102     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
2103     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
2104     ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr);
2105     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2106     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2107     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
2108     ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2109     ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2110     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
2111   }
2112   ierr = VecDestroy(&copy_standard_rhs);CHKERRQ(ierr);
2113   /* BDDC rhs */
2114   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
2115   if (pcbddc->switch_static) {
2116     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
2117   }
2118   /* apply BDDC */
2119   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
2120   /* Application of B_delta and assembling of rhs for fetidp fluxes */
2121   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
2122   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
2123   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2124   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2125   PetscFunctionReturn(0);
2126 }
2127 
2128 #undef __FUNCT__
2129 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
2130 /*@
2131  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side
2132 
2133    Collective
2134 
2135    Input Parameters:
2136 +  fetidp_mat      - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
2137 -  standard_rhs    - the right-hand side of the original linear system
2138 
2139    Output Parameters:
2140 .  fetidp_flux_rhs - the right-hand side for the FETI-DP linear system
2141 
2142    Level: developer
2143 
2144    Notes:
2145 
2146 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
2147 @*/
2148 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
2149 {
2150   FETIDPMat_ctx  mat_ctx;
2151   PetscErrorCode ierr;
2152 
2153   PetscFunctionBegin;
2154   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
2155   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
2156   PetscFunctionReturn(0);
2157 }
2158 /* -------------------------------------------------------------------------- */
2159 
2160 #undef __FUNCT__
2161 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
2162 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
2163 {
2164   FETIDPMat_ctx  mat_ctx;
2165   PC_IS*         pcis;
2166   PC_BDDC*       pcbddc;
2167   PetscErrorCode ierr;
2168 
2169   PetscFunctionBegin;
2170   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
2171   pcis = (PC_IS*)mat_ctx->pc->data;
2172   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
2173 
2174   /* apply B_delta^T */
2175   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2176   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2177   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
2178   /* compute rhs for BDDC application */
2179   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
2180   if (pcbddc->switch_static) {
2181     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
2182   }
2183   /* apply BDDC */
2184   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
2185   /* put values into standard global vector */
2186   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2187   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2188   if (!pcbddc->switch_static) {
2189     /* compute values into the interior if solved for the partially subassembled Schur complement */
2190     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
2191     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
2192     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
2193   }
2194   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2195   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2196   /* final change of basis if needed
2197      Is also sums the dirichlet part removed during RHS assembling */
2198   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
2199   PetscFunctionReturn(0);
2200 }
2201 
2202 #undef __FUNCT__
2203 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
2204 /*@
2205  PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system
2206 
2207    Collective
2208 
2209    Input Parameters:
2210 +  fetidp_mat      - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
2211 -  fetidp_flux_sol - the solution of the FETI-DP linear system
2212 
2213    Output Parameters:
2214 .  standard_sol    - the solution defined on the physical domain
2215 
2216    Level: developer
2217 
2218    Notes:
2219 
2220 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
2221 @*/
2222 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
2223 {
2224   FETIDPMat_ctx  mat_ctx;
2225   PetscErrorCode ierr;
2226 
2227   PetscFunctionBegin;
2228   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
2229   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
2230   PetscFunctionReturn(0);
2231 }
2232 /* -------------------------------------------------------------------------- */
2233 
2234 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
2235 extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
2236 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
2237 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
2238 extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
2239 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
2240 
2241 #undef __FUNCT__
2242 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
2243 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
2244 {
2245 
2246   FETIDPMat_ctx  fetidpmat_ctx;
2247   Mat            newmat;
2248   FETIDPPC_ctx   fetidppc_ctx;
2249   PC             newpc;
2250   MPI_Comm       comm;
2251   PetscErrorCode ierr;
2252 
2253   PetscFunctionBegin;
2254   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
2255   /* FETIDP linear matrix */
2256   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
2257   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
2258   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
2259   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
2260   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
2261   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
2262   ierr = MatSetUp(newmat);CHKERRQ(ierr);
2263   /* FETIDP preconditioner */
2264   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
2265   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
2266   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
2267   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
2268   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
2269   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
2270   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
2271   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
2272   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
2273   ierr = PCSetUp(newpc);CHKERRQ(ierr);
2274   /* return pointers for objects created */
2275   *fetidp_mat=newmat;
2276   *fetidp_pc=newpc;
2277   PetscFunctionReturn(0);
2278 }
2279 
2280 #undef __FUNCT__
2281 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
2282 /*@
2283  PCBDDCCreateFETIDPOperators - Create FETI-DP operators
2284 
2285    Collective
2286 
2287    Input Parameters:
2288 .  pc - the BDDC preconditioning context (setup should have been called before)
2289 
2290    Output Parameters:
2291 +  fetidp_mat - shell FETI-DP matrix object
2292 -  fetidp_pc  - shell Dirichlet preconditioner for FETI-DP matrix
2293 
2294    Options Database Keys:
2295 .    -fetidp_fullyredundant <false> - use or not a fully redundant set of Lagrange multipliers
2296 
2297    Level: developer
2298 
2299    Notes:
2300      Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose
2301 
2302 .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
2303 @*/
2304 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
2305 {
2306   PetscErrorCode ierr;
2307 
2308   PetscFunctionBegin;
2309   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2310   if (pc->setupcalled) {
2311     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
2312   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
2313   PetscFunctionReturn(0);
2314 }
2315 /* -------------------------------------------------------------------------- */
2316 /*MC
2317    PCBDDC - Balancing Domain Decomposition by Constraints.
2318 
2319    An implementation of the BDDC preconditioner based on
2320 
2321 .vb
2322    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
2323    [2] A. Klawonn and O. B. Widlund. "Dual-Primal FETI Methods for Linear Elasticity", http://cs.nyu.edu/csweb/Research/TechReports/TR2004-855/TR2004-855.pdf
2324    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
2325    [4] C. Pechstein and C. R. Dohrmann. "Modern domain decomposition methods BDDC, deluxe scaling, and an algebraic approach", Seminar talk, Linz, December 2013, http://people.ricam.oeaw.ac.at/c.pechstein/pechstein-bddc2013.pdf
2326 .ve
2327 
2328    The matrix to be preconditioned (Pmat) must be of type MATIS.
2329 
2330    Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
2331 
2332    It also works with unsymmetric and indefinite problems.
2333 
2334    Unlike 'conventional' interface preconditioners, PCBDDC iterates over all degrees of freedom, not just those on the interface. This allows the use of approximate solvers on the subdomains.
2335 
2336    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace()
2337 
2338    Boundary nodes are split in vertices, edges and faces classes using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph()
2339    Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesLocalIS()
2340 
2341    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.
2342 
2343    Change of basis is performed similarly to [2] when requested. When more than one constraint is present on a single connected component (i.e. an edge or a face), a robust method based on local QR factorizations is used.
2344    User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()
2345 
2346    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.
2347 
2348    Adaptive selection of primal constraints [4] is supported for SPD systems with high-contrast in the coefficients if MUMPS is present. Future versions of the code will also consider using MKL_PARDISO or PASTIX.
2349 
2350    An experimental interface to the FETI-DP method is available. FETI-DP operators could be created using PCBDDCCreateFETIDPOperators(). A stand-alone class for the FETI-DP method will be provided in the next releases.
2351    Deluxe scaling is not supported yet for FETI-DP.
2352 
2353    Options Database Keys (some of them, run with -h for a complete list):
2354 
2355 .    -pc_bddc_use_vertices <true> - use or not vertices in primal space
2356 .    -pc_bddc_use_edges <true> - use or not edges in primal space
2357 .    -pc_bddc_use_faces <false> - use or not faces in primal space
2358 .    -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
2359 .    -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
2360 .    -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
2361 .    -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
2362 .    -pc_bddc_levels <0> - maximum number of levels for multilevel
2363 .    -pc_bddc_coarsening_ratio <8> - number of subdomains which will be aggregated together at the coarser level (e.g. H/h ratio at the coarser level, significative only in the multilevel case)
2364 .    -pc_bddc_redistribute <0> - size of a subset of processors where the coarse problem will be remapped (the value is ignored if not at the coarsest level)
2365 .    -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
2366 .    -pc_bddc_schur_layers <-1> - select the economic version of deluxe scaling by specifying the number of layers (-1 corresponds to the original deluxe scaling)
2367 .    -pc_bddc_adaptive_threshold <0.0> - when a value greater than one is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS installed)
2368 -    -pc_bddc_check_level <0> - set verbosity level of debugging output
2369 
2370    Options for Dirichlet, Neumann or coarse solver can be set with
2371 .vb
2372       -pc_bddc_dirichlet_
2373       -pc_bddc_neumann_
2374       -pc_bddc_coarse_
2375 .ve
2376    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.
2377 
2378    When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
2379 .vb
2380       -pc_bddc_dirichlet_lN_
2381       -pc_bddc_neumann_lN_
2382       -pc_bddc_coarse_lN_
2383 .ve
2384    Note that level number ranges from the finest (0) to the coarsest (N).
2385    In order to specify options for the BDDC operators at the coarser levels (and not for the solvers), prepend -pc_bddc_coarse_ or -pc_bddc_coarse_l to the option, e.g.
2386 .vb
2387      -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
2388 .ve
2389    will use a threshold of 5 for constraints' selection at the first coarse level and will redistribute the coarse problem of the first coarse level on 3 processors
2390 
2391    Level: intermediate
2392 
2393    Developer notes:
2394 
2395    Contributed by Stefano Zampini
2396 
2397 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
2398 M*/
2399 
2400 #undef __FUNCT__
2401 #define __FUNCT__ "PCCreate_BDDC"
2402 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
2403 {
2404   PetscErrorCode      ierr;
2405   PC_BDDC             *pcbddc;
2406 
2407   PetscFunctionBegin;
2408   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
2409   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
2410   pc->data  = (void*)pcbddc;
2411 
2412   /* create PCIS data structure */
2413   ierr = PCISCreate(pc);CHKERRQ(ierr);
2414 
2415   /* BDDC customization */
2416   pcbddc->use_local_adj       = PETSC_TRUE;
2417   pcbddc->use_vertices        = PETSC_TRUE;
2418   pcbddc->use_edges           = PETSC_TRUE;
2419   pcbddc->use_faces           = PETSC_FALSE;
2420   pcbddc->use_change_of_basis = PETSC_FALSE;
2421   pcbddc->use_change_on_faces = PETSC_FALSE;
2422   pcbddc->switch_static       = PETSC_FALSE;
2423   pcbddc->use_nnsp_true       = PETSC_FALSE;
2424   pcbddc->use_qr_single       = PETSC_FALSE;
2425   pcbddc->symmetric_primal    = PETSC_TRUE;
2426   pcbddc->benign_saddle_point = PETSC_FALSE;
2427   pcbddc->dbg_flag            = 0;
2428   /* private */
2429   pcbddc->local_primal_size          = 0;
2430   pcbddc->local_primal_size_cc       = 0;
2431   pcbddc->local_primal_ref_node      = 0;
2432   pcbddc->local_primal_ref_mult      = 0;
2433   pcbddc->n_vertices                 = 0;
2434   pcbddc->primal_indices_local_idxs  = 0;
2435   pcbddc->recompute_topography       = PETSC_FALSE;
2436   pcbddc->coarse_size                = -1;
2437   pcbddc->new_primal_space           = PETSC_FALSE;
2438   pcbddc->new_primal_space_local     = PETSC_FALSE;
2439   pcbddc->global_primal_indices      = 0;
2440   pcbddc->onearnullspace             = 0;
2441   pcbddc->onearnullvecs_state        = 0;
2442   pcbddc->user_primal_vertices       = 0;
2443   pcbddc->NullSpace                  = 0;
2444   pcbddc->temp_solution              = 0;
2445   pcbddc->original_rhs               = 0;
2446   pcbddc->local_mat                  = 0;
2447   pcbddc->ChangeOfBasisMatrix        = 0;
2448   pcbddc->user_ChangeOfBasisMatrix   = 0;
2449   pcbddc->new_global_mat             = 0;
2450   pcbddc->coarse_vec                 = 0;
2451   pcbddc->coarse_ksp                 = 0;
2452   pcbddc->coarse_phi_B               = 0;
2453   pcbddc->coarse_phi_D               = 0;
2454   pcbddc->coarse_psi_B               = 0;
2455   pcbddc->coarse_psi_D               = 0;
2456   pcbddc->vec1_P                     = 0;
2457   pcbddc->vec1_R                     = 0;
2458   pcbddc->vec2_R                     = 0;
2459   pcbddc->local_auxmat1              = 0;
2460   pcbddc->local_auxmat2              = 0;
2461   pcbddc->R_to_B                     = 0;
2462   pcbddc->R_to_D                     = 0;
2463   pcbddc->ksp_D                      = 0;
2464   pcbddc->ksp_R                      = 0;
2465   pcbddc->NeumannBoundaries          = 0;
2466   pcbddc->NeumannBoundariesLocal     = 0;
2467   pcbddc->DirichletBoundaries        = 0;
2468   pcbddc->DirichletBoundariesLocal   = 0;
2469   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
2470   pcbddc->n_ISForDofs                = 0;
2471   pcbddc->n_ISForDofsLocal           = 0;
2472   pcbddc->ISForDofs                  = 0;
2473   pcbddc->ISForDofsLocal             = 0;
2474   pcbddc->ConstraintMatrix           = 0;
2475   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
2476   pcbddc->coarse_loc_to_glob         = 0;
2477   pcbddc->coarsening_ratio           = 8;
2478   pcbddc->coarse_adj_red             = 0;
2479   pcbddc->current_level              = 0;
2480   pcbddc->max_levels                 = 0;
2481   pcbddc->use_coarse_estimates       = PETSC_FALSE;
2482   pcbddc->redistribute_coarse        = 0;
2483   pcbddc->coarse_subassembling       = 0;
2484   pcbddc->coarse_subassembling_init  = 0;
2485 
2486   /* benign subspace trick */
2487   pcbddc->B0_ncol                    = 0;
2488   pcbddc->B0_cols                    = NULL;
2489   pcbddc->B0_vals                    = NULL;
2490   pcbddc->benign_change              = 0;
2491   pcbddc->benign_vec                 = 0;
2492   pcbddc->benign_original_mat        = 0;
2493   pcbddc->benign_sf                  = 0;
2494 
2495   /* create local graph structure */
2496   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
2497 
2498   /* scaling */
2499   pcbddc->work_scaling          = 0;
2500   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
2501   pcbddc->faster_deluxe         = PETSC_FALSE;
2502 
2503   /* create sub schurs structure */
2504   ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr);
2505   pcbddc->sub_schurs_rebuild     = PETSC_FALSE;
2506   pcbddc->sub_schurs_layers      = -1;
2507   pcbddc->sub_schurs_use_useradj = PETSC_FALSE;
2508 
2509   pcbddc->computed_rowadj = PETSC_FALSE;
2510 
2511   /* adaptivity */
2512   pcbddc->adaptive_threshold      = 0.0;
2513   pcbddc->adaptive_nmax           = 0;
2514   pcbddc->adaptive_nmin           = 0;
2515 
2516   /* function pointers */
2517   pc->ops->apply               = PCApply_BDDC;
2518   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2519   pc->ops->setup               = PCSetUp_BDDC;
2520   pc->ops->destroy             = PCDestroy_BDDC;
2521   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
2522   pc->ops->view                = 0;
2523   pc->ops->applyrichardson     = 0;
2524   pc->ops->applysymmetricleft  = 0;
2525   pc->ops->applysymmetricright = 0;
2526   pc->ops->presolve            = PCPreSolve_BDDC;
2527   pc->ops->postsolve           = PCPostSolve_BDDC;
2528 
2529   /* composing function */
2530   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
2531   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
2532   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
2533   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
2534   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
2535   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
2536   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
2537   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
2538   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2539   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
2540   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2541   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
2542   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2543   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
2544   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2545   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
2546   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
2547   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
2548   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
2549   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
2550   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
2551   PetscFunctionReturn(0);
2552 }
2553 
2554