xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 7d24d15547f91c03459755ce0a63d25236fb26a1)
1 /* TODOLIST
2 
3    ConstraintsSetup
4    - tolerances for constraints as an option (take care of single precision!)
5    - Can MAT_IGNORE_ZERO_ENTRIES be used for Constraints Matrix?
6 
7    Solvers
8    - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers)
9    - Propagate ksp prefixes for solvers to mat objects?
10    - Propagate nearnullspace info among levels
11 
12    User interface
13    - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access)
14    - Provide PCApplyTranpose_BDDC
15    - DofSplitting and DM attached to pc?
16 
17    Debugging output
18    - Better management of verbosity levels of debugging output
19 
20    Build
21    - make runexe59
22 
23    Extra
24    - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
25    - Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver?
26    - add support for computing h,H and related using coordinates?
27    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
28    - Better management in PCIS code
29    - BDDC with MG framework?
30 
31    FETIDP
32    - Move FETIDP code to its own classes
33 
34    MATIS related operations contained in BDDC code
35    - Provide general case for subassembling
36    - Preallocation routines in MatISGetMPIAXAIJ
37 
38 */
39 
40 /* ----------------------------------------------------------------------------------------------------------------------------------------------
41    Implementation of BDDC preconditioner based on:
42    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
43    ---------------------------------------------------------------------------------------------------------------------------------------------- */
44 
45 #include "bddc.h" /*I "petscpc.h" I*/  /* includes for fortran wrappers */
46 #include "bddcprivate.h"
47 #include <petscblaslapack.h>
48 
49 /* -------------------------------------------------------------------------- */
50 #undef __FUNCT__
51 #define __FUNCT__ "PCSetFromOptions_BDDC"
52 PetscErrorCode PCSetFromOptions_BDDC(PC pc)
53 {
54   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
55   PetscErrorCode ierr;
56 
57   PetscFunctionBegin;
58   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
59   /* Verbose debugging */
60   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
61   /* Primal space cumstomization */
62   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
63   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
64   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
65   /* Change of basis */
66   ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr);
67   ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr);
68   if (!pcbddc->use_change_of_basis) {
69     pcbddc->use_change_on_faces = PETSC_FALSE;
70   }
71   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
72   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);
73   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
74   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
75   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
76   ierr = PetscOptionsTail();CHKERRQ(ierr);
77   PetscFunctionReturn(0);
78 }
79 /* -------------------------------------------------------------------------- */
80 #undef __FUNCT__
81 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
82 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
83 {
84   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
85   PetscErrorCode ierr;
86 
87   PetscFunctionBegin;
88   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
89   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
90   pcbddc->user_primal_vertices = PrimalVertices;
91   PetscFunctionReturn(0);
92 }
93 #undef __FUNCT__
94 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
95 /*@
96  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
97 
98    Not collective
99 
100    Input Parameters:
101 +  pc - the preconditioning context
102 -  PrimalVertices - index set of primal vertices in local numbering
103 
104    Level: intermediate
105 
106    Notes:
107 
108 .seealso: PCBDDC
109 @*/
110 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
111 {
112   PetscErrorCode ierr;
113 
114   PetscFunctionBegin;
115   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
116   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
117   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
118   PetscFunctionReturn(0);
119 }
120 /* -------------------------------------------------------------------------- */
121 #undef __FUNCT__
122 #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
123 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
124 {
125   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
126 
127   PetscFunctionBegin;
128   pcbddc->coarsening_ratio = k;
129   PetscFunctionReturn(0);
130 }
131 
132 #undef __FUNCT__
133 #define __FUNCT__ "PCBDDCSetCoarseningRatio"
134 /*@
135  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
136 
137    Logically collective on PC
138 
139    Input Parameters:
140 +  pc - the preconditioning context
141 -  k - coarsening ratio (H/h at the coarser level)
142 
143    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
144 
145    Level: intermediate
146 
147    Notes:
148 
149 .seealso: PCBDDC
150 @*/
151 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
152 {
153   PetscErrorCode ierr;
154 
155   PetscFunctionBegin;
156   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
157   PetscValidLogicalCollectiveInt(pc,k,2);
158   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
163 #undef __FUNCT__
164 #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
165 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
166 {
167   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
168 
169   PetscFunctionBegin;
170   pcbddc->use_exact_dirichlet_trick = flg;
171   PetscFunctionReturn(0);
172 }
173 
174 #undef __FUNCT__
175 #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
176 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
177 {
178   PetscErrorCode ierr;
179 
180   PetscFunctionBegin;
181   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
182   PetscValidLogicalCollectiveBool(pc,flg,2);
183   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
184   PetscFunctionReturn(0);
185 }
186 
187 #undef __FUNCT__
188 #define __FUNCT__ "PCBDDCSetLevel_BDDC"
189 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
190 {
191   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
192 
193   PetscFunctionBegin;
194   pcbddc->current_level = level;
195   PetscFunctionReturn(0);
196 }
197 
198 #undef __FUNCT__
199 #define __FUNCT__ "PCBDDCSetLevel"
200 PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
201 {
202   PetscErrorCode ierr;
203 
204   PetscFunctionBegin;
205   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
206   PetscValidLogicalCollectiveInt(pc,level,2);
207   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "PCBDDCSetLevels_BDDC"
213 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
214 {
215   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
216 
217   PetscFunctionBegin;
218   pcbddc->max_levels = levels;
219   PetscFunctionReturn(0);
220 }
221 
222 #undef __FUNCT__
223 #define __FUNCT__ "PCBDDCSetLevels"
224 /*@
225  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
226 
227    Logically collective on PC
228 
229    Input Parameters:
230 +  pc - the preconditioning context
231 -  levels - the maximum number of levels (max 9)
232 
233    Default value is 0, i.e. traditional one-level BDDC
234 
235    Level: intermediate
236 
237    Notes:
238 
239 .seealso: PCBDDC
240 @*/
241 PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
242 {
243   PetscErrorCode ierr;
244 
245   PetscFunctionBegin;
246   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
247   PetscValidLogicalCollectiveInt(pc,levels,2);
248   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
249   PetscFunctionReturn(0);
250 }
251 /* -------------------------------------------------------------------------- */
252 
253 #undef __FUNCT__
254 #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
255 static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
256 {
257   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
258   PetscErrorCode ierr;
259 
260   PetscFunctionBegin;
261   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
262   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
263   pcbddc->NullSpace = NullSpace;
264   PetscFunctionReturn(0);
265 }
266 
267 #undef __FUNCT__
268 #define __FUNCT__ "PCBDDCSetNullSpace"
269 /*@
270  PCBDDCSetNullSpace - Set nullspace for BDDC operator
271 
272    Logically collective on PC and MatNullSpace
273 
274    Input Parameters:
275 +  pc - the preconditioning context
276 -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
277 
278    Level: intermediate
279 
280    Notes:
281 
282 .seealso: PCBDDC
283 @*/
284 PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
285 {
286   PetscErrorCode ierr;
287 
288   PetscFunctionBegin;
289   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
290   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
291   PetscCheckSameComm(pc,1,NullSpace,2);
292   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
293   PetscFunctionReturn(0);
294 }
295 /* -------------------------------------------------------------------------- */
296 
297 #undef __FUNCT__
298 #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
299 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
300 {
301   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
302   PetscErrorCode ierr;
303 
304   PetscFunctionBegin;
305   /* last user setting takes precendence -> destroy any other customization */
306   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
307   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
308   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
309   pcbddc->DirichletBoundaries = DirichletBoundaries;
310   pcbddc->recompute_topography = PETSC_TRUE;
311   PetscFunctionReturn(0);
312 }
313 
314 #undef __FUNCT__
315 #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
316 /*@
317  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
318 
319    Collective
320 
321    Input Parameters:
322 +  pc - the preconditioning context
323 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
324 
325    Level: intermediate
326 
327    Notes: Any process can list any global node
328 
329 .seealso: PCBDDC
330 @*/
331 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
332 {
333   PetscErrorCode ierr;
334 
335   PetscFunctionBegin;
336   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
337   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
338   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
339   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
340   PetscFunctionReturn(0);
341 }
342 /* -------------------------------------------------------------------------- */
343 
344 #undef __FUNCT__
345 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
346 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
347 {
348   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
349   PetscErrorCode ierr;
350 
351   PetscFunctionBegin;
352   /* last user setting takes precendence -> destroy any other customization */
353   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
354   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
355   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
356   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
357   pcbddc->recompute_topography = PETSC_TRUE;
358   PetscFunctionReturn(0);
359 }
360 
361 #undef __FUNCT__
362 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
363 /*@
364  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
365 
366    Collective
367 
368    Input Parameters:
369 +  pc - the preconditioning context
370 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
371 
372    Level: intermediate
373 
374    Notes:
375 
376 .seealso: PCBDDC
377 @*/
378 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
379 {
380   PetscErrorCode ierr;
381 
382   PetscFunctionBegin;
383   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
384   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
385   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
386   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
387   PetscFunctionReturn(0);
388 }
389 /* -------------------------------------------------------------------------- */
390 
391 #undef __FUNCT__
392 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
393 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
394 {
395   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
396   PetscErrorCode ierr;
397 
398   PetscFunctionBegin;
399   /* last user setting takes precendence -> destroy any other customization */
400   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
401   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
402   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
403   pcbddc->NeumannBoundaries = NeumannBoundaries;
404   pcbddc->recompute_topography = PETSC_TRUE;
405   PetscFunctionReturn(0);
406 }
407 
408 #undef __FUNCT__
409 #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
410 /*@
411  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
412 
413    Collective
414 
415    Input Parameters:
416 +  pc - the preconditioning context
417 -  NeumannBoundaries - parallel IS defining the Neumann boundaries
418 
419    Level: intermediate
420 
421    Notes: Any process can list any global node
422 
423 .seealso: PCBDDC
424 @*/
425 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
426 {
427   PetscErrorCode ierr;
428 
429   PetscFunctionBegin;
430   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
431   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
432   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
433   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
434   PetscFunctionReturn(0);
435 }
436 /* -------------------------------------------------------------------------- */
437 
438 #undef __FUNCT__
439 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
440 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
441 {
442   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
443   PetscErrorCode ierr;
444 
445   PetscFunctionBegin;
446   /* last user setting takes precendence -> destroy any other customization */
447   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
448   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
449   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
450   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
451   pcbddc->recompute_topography = PETSC_TRUE;
452   PetscFunctionReturn(0);
453 }
454 
455 #undef __FUNCT__
456 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
457 /*@
458  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
459 
460    Collective
461 
462    Input Parameters:
463 +  pc - the preconditioning context
464 -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
465 
466    Level: intermediate
467 
468    Notes:
469 
470 .seealso: PCBDDC
471 @*/
472 PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
473 {
474   PetscErrorCode ierr;
475 
476   PetscFunctionBegin;
477   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
478   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
479   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
480   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
481   PetscFunctionReturn(0);
482 }
483 /* -------------------------------------------------------------------------- */
484 
485 #undef __FUNCT__
486 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
487 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
488 {
489   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
490 
491   PetscFunctionBegin;
492   *DirichletBoundaries = pcbddc->DirichletBoundaries;
493   PetscFunctionReturn(0);
494 }
495 
496 #undef __FUNCT__
497 #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
498 /*@
499  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
500 
501    Collective
502 
503    Input Parameters:
504 .  pc - the preconditioning context
505 
506    Output Parameters:
507 .  DirichletBoundaries - index set defining the Dirichlet boundaries
508 
509    Level: intermediate
510 
511    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
512 
513 .seealso: PCBDDC
514 @*/
515 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
516 {
517   PetscErrorCode ierr;
518 
519   PetscFunctionBegin;
520   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
521   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
522   PetscFunctionReturn(0);
523 }
524 /* -------------------------------------------------------------------------- */
525 
526 #undef __FUNCT__
527 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
528 static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
529 {
530   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
531 
532   PetscFunctionBegin;
533   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
534   PetscFunctionReturn(0);
535 }
536 
537 #undef __FUNCT__
538 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
539 /*@
540  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
541 
542    Collective
543 
544    Input Parameters:
545 .  pc - the preconditioning context
546 
547    Output Parameters:
548 .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
549 
550    Level: intermediate
551 
552    Notes:
553 
554 .seealso: PCBDDC
555 @*/
556 PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
557 {
558   PetscErrorCode ierr;
559 
560   PetscFunctionBegin;
561   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
562   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
563   PetscFunctionReturn(0);
564 }
565 /* -------------------------------------------------------------------------- */
566 
567 #undef __FUNCT__
568 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
569 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
570 {
571   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
572 
573   PetscFunctionBegin;
574   *NeumannBoundaries = pcbddc->NeumannBoundaries;
575   PetscFunctionReturn(0);
576 }
577 
578 #undef __FUNCT__
579 #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
580 /*@
581  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
582 
583    Collective
584 
585    Input Parameters:
586 .  pc - the preconditioning context
587 
588    Output Parameters:
589 .  NeumannBoundaries - index set defining the Neumann boundaries
590 
591    Level: intermediate
592 
593    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
594 
595 .seealso: PCBDDC
596 @*/
597 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
598 {
599   PetscErrorCode ierr;
600 
601   PetscFunctionBegin;
602   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
603   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
604   PetscFunctionReturn(0);
605 }
606 /* -------------------------------------------------------------------------- */
607 
608 #undef __FUNCT__
609 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
610 static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
611 {
612   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
613 
614   PetscFunctionBegin;
615   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
616   PetscFunctionReturn(0);
617 }
618 
619 #undef __FUNCT__
620 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
621 /*@
622  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
623 
624    Collective
625 
626    Input Parameters:
627 .  pc - the preconditioning context
628 
629    Output Parameters:
630 .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
631 
632    Level: intermediate
633 
634    Notes:
635 
636 .seealso: PCBDDC
637 @*/
638 PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
639 {
640   PetscErrorCode ierr;
641 
642   PetscFunctionBegin;
643   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
644   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
645   PetscFunctionReturn(0);
646 }
647 /* -------------------------------------------------------------------------- */
648 
649 #undef __FUNCT__
650 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
651 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
652 {
653   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
654   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
655   PetscErrorCode ierr;
656 
657   PetscFunctionBegin;
658   /* free old CSR */
659   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
660   /* TODO: PCBDDCGraphSetAdjacency */
661   /* get CSR into graph structure */
662   if (copymode == PETSC_COPY_VALUES) {
663     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
664     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
665     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
666     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
667   } else if (copymode == PETSC_OWN_POINTER) {
668     mat_graph->xadj = (PetscInt*)xadj;
669     mat_graph->adjncy = (PetscInt*)adjncy;
670   }
671   mat_graph->nvtxs_csr = nvtxs;
672   PetscFunctionReturn(0);
673 }
674 
675 #undef __FUNCT__
676 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
677 /*@
678  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
679 
680    Not collective
681 
682    Input Parameters:
683 +  pc - the preconditioning context
684 .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
685 .  xadj, adjncy - the CSR graph
686 -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
687 
688    Level: intermediate
689 
690    Notes:
691 
692 .seealso: PCBDDC,PetscCopyMode
693 @*/
694 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
695 {
696   void (*f)(void) = 0;
697   PetscErrorCode ierr;
698 
699   PetscFunctionBegin;
700   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
701   PetscValidIntPointer(xadj,3);
702   PetscValidIntPointer(xadj,4);
703   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
704     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
705   }
706   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
707   /* free arrays if PCBDDC is not the PC type */
708   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
709   if (!f && copymode == PETSC_OWN_POINTER) {
710     ierr = PetscFree(xadj);CHKERRQ(ierr);
711     ierr = PetscFree(adjncy);CHKERRQ(ierr);
712   }
713   PetscFunctionReturn(0);
714 }
715 /* -------------------------------------------------------------------------- */
716 
717 #undef __FUNCT__
718 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
719 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
720 {
721   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
722   PetscInt i;
723   PetscErrorCode ierr;
724 
725   PetscFunctionBegin;
726   /* Destroy ISes if they were already set */
727   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
728     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
729   }
730   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
731   /* last user setting takes precendence -> destroy any other customization */
732   for (i=0;i<pcbddc->n_ISForDofs;i++) {
733     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
734   }
735   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
736   pcbddc->n_ISForDofs = 0;
737   /* allocate space then set */
738   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
739   for (i=0;i<n_is;i++) {
740     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
741     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
742   }
743   pcbddc->n_ISForDofsLocal=n_is;
744   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
745   PetscFunctionReturn(0);
746 }
747 
748 #undef __FUNCT__
749 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
750 /*@
751  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
752 
753    Collective
754 
755    Input Parameters:
756 +  pc - the preconditioning context
757 -  n_is - number of index sets defining the fields
758 .  ISForDofs - array of IS describing the fields in local ordering
759 
760    Level: intermediate
761 
762    Notes: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to a different field.
763 
764 .seealso: PCBDDC
765 @*/
766 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
767 {
768   PetscInt       i;
769   PetscErrorCode ierr;
770 
771   PetscFunctionBegin;
772   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
773   PetscValidLogicalCollectiveInt(pc,n_is,2);
774   for (i=0;i<n_is;i++) {
775     PetscCheckSameComm(pc,1,ISForDofs[i],3);
776     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
777   }
778   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
779   PetscFunctionReturn(0);
780 }
781 /* -------------------------------------------------------------------------- */
782 
783 #undef __FUNCT__
784 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
785 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
786 {
787   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
788   PetscInt i;
789   PetscErrorCode ierr;
790 
791   PetscFunctionBegin;
792   /* Destroy ISes if they were already set */
793   for (i=0;i<pcbddc->n_ISForDofs;i++) {
794     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
795   }
796   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
797   /* last user setting takes precendence -> destroy any other customization */
798   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
799     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
800   }
801   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
802   pcbddc->n_ISForDofsLocal = 0;
803   /* allocate space then set */
804   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
805   for (i=0;i<n_is;i++) {
806     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
807     pcbddc->ISForDofs[i]=ISForDofs[i];
808   }
809   pcbddc->n_ISForDofs=n_is;
810   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
811   PetscFunctionReturn(0);
812 }
813 
814 #undef __FUNCT__
815 #define __FUNCT__ "PCBDDCSetDofsSplitting"
816 /*@
817  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
818 
819    Collective
820 
821    Input Parameters:
822 +  pc - the preconditioning context
823 -  n_is - number of index sets defining the fields
824 .  ISForDofs - array of IS describing the fields in global ordering
825 
826    Level: intermediate
827 
828    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.
829 
830 .seealso: PCBDDC
831 @*/
832 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
833 {
834   PetscInt       i;
835   PetscErrorCode ierr;
836 
837   PetscFunctionBegin;
838   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
839   PetscValidLogicalCollectiveInt(pc,n_is,2);
840   for (i=0;i<n_is;i++) {
841     PetscCheckSameComm(pc,1,ISForDofs[i],3);
842     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
843   }
844   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
845   PetscFunctionReturn(0);
846 }
847 /* -------------------------------------------------------------------------- */
848 #undef __FUNCT__
849 #define __FUNCT__ "PCPreSolve_BDDC"
850 /* -------------------------------------------------------------------------- */
851 /*
852    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
853                      guess if a transformation of basis approach has been selected.
854 
855    Input Parameter:
856 +  pc - the preconditioner contex
857 
858    Application Interface Routine: PCPreSolve()
859 
860    Notes:
861    The interface routine PCPreSolve() is not usually called directly by
862    the user, but instead is called by KSPSolve().
863 */
864 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
865 {
866   PetscErrorCode ierr;
867   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
868   PC_IS          *pcis = (PC_IS*)(pc->data);
869   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
870   Mat            temp_mat;
871   IS             dirIS;
872   Vec            used_vec;
873   PetscBool      guess_nonzero;
874 
875   PetscFunctionBegin;
876   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
877   if (ksp) {
878     PetscBool iscg;
879     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
880     if (!iscg) {
881       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
882     }
883   }
884   /* Creates parallel work vectors used in presolve */
885   if (!pcbddc->original_rhs) {
886     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
887   }
888   if (!pcbddc->temp_solution) {
889     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
890   }
891   if (x) {
892     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
893     used_vec = x;
894   } else {
895     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
896     used_vec = pcbddc->temp_solution;
897     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
898   }
899   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
900   if (ksp) {
901     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
902     if (!guess_nonzero) {
903       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
904     }
905   }
906 
907   /* store the original rhs */
908   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
909 
910   /* Take into account zeroed rows -> change rhs and store solution removed */
911   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
912   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
913   if (rhs && dirIS) {
914     PetscInt    dirsize,i,*is_indices;
915     PetscScalar *array_x,*array_diagonal;
916 
917     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
918     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
919     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
920     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
921     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
922     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
923     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
924     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
925     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
926     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
927     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
928     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
929     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
930     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
931     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
932     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
933 
934     /* remove the computed solution from the rhs */
935     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
936     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
937     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
938   }
939 
940   /* store partially computed solution and set initial guess */
941   if (x) {
942     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
943     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
944     if (pcbddc->use_exact_dirichlet_trick) {
945       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
946       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
947       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
948       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
949       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
950       if (ksp) {
951         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
952       }
953     }
954   }
955 
956   /* prepare MatMult and rhs for solver */
957   if (pcbddc->use_change_of_basis) {
958     /* swap pointers for local matrices */
959     temp_mat = matis->A;
960     matis->A = pcbddc->local_mat;
961     pcbddc->local_mat = temp_mat;
962     if (rhs) {
963       /* Get local rhs and apply transformation of basis */
964       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
965       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
966       /* from original basis to modified basis */
967       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
968       /* put back modified values into the global vec using INSERT_VALUES copy mode */
969       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
970       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
971     }
972   }
973 
974   /* remove nullspace if present */
975   if (ksp && pcbddc->NullSpace) {
976     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
977     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
978   }
979   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
980   PetscFunctionReturn(0);
981 }
982 /* -------------------------------------------------------------------------- */
983 #undef __FUNCT__
984 #define __FUNCT__ "PCPostSolve_BDDC"
985 /* -------------------------------------------------------------------------- */
986 /*
987    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
988                      approach has been selected. Also, restores rhs to its original state.
989 
990    Input Parameter:
991 +  pc - the preconditioner contex
992 
993    Application Interface Routine: PCPostSolve()
994 
995    Notes:
996    The interface routine PCPostSolve() is not usually called directly by
997    the user, but instead is called by KSPSolve().
998 */
999 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1000 {
1001   PetscErrorCode ierr;
1002   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1003   PC_IS          *pcis   = (PC_IS*)(pc->data);
1004   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
1005   Mat            temp_mat;
1006 
1007   PetscFunctionBegin;
1008   if (pcbddc->use_change_of_basis) {
1009     /* swap pointers for local matrices */
1010     temp_mat = matis->A;
1011     matis->A = pcbddc->local_mat;
1012     pcbddc->local_mat = temp_mat;
1013   }
1014   if (pcbddc->use_change_of_basis && x) {
1015     /* Get Local boundary and apply transformation of basis to solution vector */
1016     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1017     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1018     /* from modified basis to original basis */
1019     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
1020     /* put back modified values into the global vec using INSERT_VALUES copy mode */
1021     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1022     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1023   }
1024   /* add solution removed in presolve */
1025   if (x) {
1026     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
1027   }
1028   /* restore rhs to its original state */
1029   if (rhs) {
1030     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1031   }
1032   PetscFunctionReturn(0);
1033 }
1034 /* -------------------------------------------------------------------------- */
1035 #undef __FUNCT__
1036 #define __FUNCT__ "PCSetUp_BDDC"
1037 /* -------------------------------------------------------------------------- */
1038 /*
1039    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1040                   by setting data structures and options.
1041 
1042    Input Parameter:
1043 +  pc - the preconditioner context
1044 
1045    Application Interface Routine: PCSetUp()
1046 
1047    Notes:
1048    The interface routine PCSetUp() is not usually called directly by
1049    the user, but instead is called by PCApply() if necessary.
1050 */
1051 PetscErrorCode PCSetUp_BDDC(PC pc)
1052 {
1053   PetscErrorCode   ierr;
1054   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
1055   MatNullSpace     nearnullspace;
1056   MatStructure     flag;
1057   PetscBool        computeis,computetopography,computesolvers;
1058   PetscBool        new_nearnullspace_provided;
1059 
1060   PetscFunctionBegin;
1061   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1062   /* PCIS does not support MatStructure flags different from SAME_PRECONDITIONER */
1063   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1064      Also, BDDC directly build the Dirichlet problem */
1065 
1066   /* split work */
1067   if (pc->setupcalled) {
1068     computeis = PETSC_FALSE;
1069     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
1070     if (flag == SAME_PRECONDITIONER) {
1071       PetscFunctionReturn(0);
1072     } else if (flag == SAME_NONZERO_PATTERN) {
1073       computetopography = PETSC_FALSE;
1074       computesolvers = PETSC_TRUE;
1075     } else { /* DIFFERENT_NONZERO_PATTERN */
1076       computetopography = PETSC_TRUE;
1077       computesolvers = PETSC_TRUE;
1078     }
1079   } else {
1080     computeis = PETSC_TRUE;
1081     computetopography = PETSC_TRUE;
1082     computesolvers = PETSC_TRUE;
1083   }
1084   if (pcbddc->recompute_topography) {
1085     computetopography = PETSC_TRUE;
1086   }
1087 
1088   /* Get stdout for dbg */
1089   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
1090     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
1091     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1092     if (pcbddc->current_level) {
1093       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
1094     }
1095   }
1096 
1097   /* Set up all the "iterative substructuring" common block without computing solvers */
1098   if (computeis) {
1099     /* HACK INTO PCIS */
1100     PC_IS* pcis = (PC_IS*)pc->data;
1101     pcis->computesolvers = PETSC_FALSE;
1102     ierr = PCISSetUp(pc);CHKERRQ(ierr);
1103     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr);
1104   }
1105 
1106   /* Analyze interface */
1107   if (computetopography) {
1108     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1109   }
1110 
1111   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1112   new_nearnullspace_provided = PETSC_FALSE;
1113   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1114   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1115     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1116       new_nearnullspace_provided = PETSC_TRUE;
1117     } else {
1118       /* determine if the two nullspaces are different (should be lightweight) */
1119       if (nearnullspace != pcbddc->onearnullspace) {
1120         new_nearnullspace_provided = PETSC_TRUE;
1121       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1122         PetscInt         i;
1123         const Vec        *nearnullvecs;
1124         PetscObjectState state;
1125         PetscInt         nnsp_size;
1126         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1127         for (i=0;i<nnsp_size;i++) {
1128           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1129           if (pcbddc->onearnullvecs_state[i] != state) {
1130             new_nearnullspace_provided = PETSC_TRUE;
1131             break;
1132           }
1133         }
1134       }
1135     }
1136   } else {
1137     if (!nearnullspace) { /* both nearnullspaces are null */
1138       new_nearnullspace_provided = PETSC_FALSE;
1139     } else { /* nearnullspace attached later */
1140       new_nearnullspace_provided = PETSC_TRUE;
1141     }
1142   }
1143 
1144   /* Setup constraints and related work vectors */
1145   /* reset primal space flags */
1146   pcbddc->new_primal_space = PETSC_FALSE;
1147   pcbddc->new_primal_space_local = PETSC_FALSE;
1148   if (computetopography || new_nearnullspace_provided) {
1149     /* It also sets the primal space flags */
1150     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1151     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1152     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
1153   }
1154 
1155   if (computesolvers || pcbddc->new_primal_space) {
1156     /* reset data */
1157     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1158     /* Create coarse and local stuffs */
1159     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1160     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
1161   }
1162   if (pcbddc->dbg_flag && pcbddc->current_level) {
1163     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
1164   }
1165   PetscFunctionReturn(0);
1166 }
1167 
1168 /* -------------------------------------------------------------------------- */
1169 /*
1170    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
1171 
1172    Input Parameters:
1173 .  pc - the preconditioner context
1174 .  r - input vector (global)
1175 
1176    Output Parameter:
1177 .  z - output vector (global)
1178 
1179    Application Interface Routine: PCApply()
1180  */
1181 #undef __FUNCT__
1182 #define __FUNCT__ "PCApply_BDDC"
1183 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1184 {
1185   PC_IS             *pcis = (PC_IS*)(pc->data);
1186   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1187   PetscErrorCode    ierr;
1188   const PetscScalar one = 1.0;
1189   const PetscScalar m_one = -1.0;
1190   const PetscScalar zero = 0.0;
1191 
1192 /* This code is similar to that provided in nn.c for PCNN
1193    NN interface preconditioner changed to BDDC
1194    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
1195 
1196   PetscFunctionBegin;
1197   if (!pcbddc->use_exact_dirichlet_trick) {
1198     /* First Dirichlet solve */
1199     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1200     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1201     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1202     /*
1203       Assembling right hand side for BDDC operator
1204       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1205       - pcis->vec1_B the interface part of the global vector z
1206     */
1207     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1208     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1209     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1210     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1211     ierr = VecCopy(r,z);CHKERRQ(ierr);
1212     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1213     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1214     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1215   } else {
1216     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1217     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1218     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1219   }
1220 
1221   /* Apply interface preconditioner
1222      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1223   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
1224 
1225   /* Apply transpose of partition of unity operator */
1226   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1227 
1228   /* Second Dirichlet solve and assembling of output */
1229   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1230   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1231   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1232   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1233   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1234   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1235   if (pcbddc->switch_static) { ierr = VecAXPY (pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1236   ierr = VecAXPY (pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
1237   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1238   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1239   PetscFunctionReturn(0);
1240 }
1241 /* -------------------------------------------------------------------------- */
1242 
1243 #undef __FUNCT__
1244 #define __FUNCT__ "PCDestroy_BDDC"
1245 PetscErrorCode PCDestroy_BDDC(PC pc)
1246 {
1247   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1248   PetscErrorCode ierr;
1249 
1250   PetscFunctionBegin;
1251   /* free data created by PCIS */
1252   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1253   /* free BDDC custom data  */
1254   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1255   /* destroy objects related to topography */
1256   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1257   /* free allocated graph structure */
1258   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1259   /* free data for scaling operator */
1260   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1261   /* free solvers stuff */
1262   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1263   /* free global vectors needed in presolve */
1264   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
1265   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1266   ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr);
1267   /* remove functions */
1268   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1269   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
1270   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1271   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
1272   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1273   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1274   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1275   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1276   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1277   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1278   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1279   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1280   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1281   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1282   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1283   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1284   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1285   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1286   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1287   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1288   /* Free the private data structure */
1289   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1290   PetscFunctionReturn(0);
1291 }
1292 /* -------------------------------------------------------------------------- */
1293 
1294 #undef __FUNCT__
1295 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
1296 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1297 {
1298   FETIDPMat_ctx  mat_ctx;
1299   PC_IS*         pcis;
1300   PC_BDDC*       pcbddc;
1301   PetscErrorCode ierr;
1302 
1303   PetscFunctionBegin;
1304   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1305   pcis = (PC_IS*)mat_ctx->pc->data;
1306   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1307 
1308   /* change of basis for physical rhs if needed
1309      It also changes the rhs in case of dirichlet boundaries */
1310   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
1311   /* store vectors for computation of fetidp final solution */
1312   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1313   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1314   /* scale rhs since it should be unassembled */
1315   /* TODO use counter scaling? (also below) */
1316   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1317   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1318   /* Apply partition of unity */
1319   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1320   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1321   if (!pcbddc->switch_static) {
1322     /* compute partially subassembled Schur complement right-hand side */
1323     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1324     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
1325     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1326     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
1327     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1328     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1329     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1330     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1331     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1332     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1333   }
1334   /* BDDC rhs */
1335   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
1336   if (pcbddc->switch_static) {
1337     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1338   }
1339   /* apply BDDC */
1340   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
1341   /* Application of B_delta and assembling of rhs for fetidp fluxes */
1342   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
1343   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
1344   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1345   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1346   /* restore original rhs */
1347   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
1348   PetscFunctionReturn(0);
1349 }
1350 
1351 #undef __FUNCT__
1352 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
1353 /*@
1354  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
1355 
1356    Collective
1357 
1358    Input Parameters:
1359 +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
1360 .  standard_rhs - the right-hand side for your linear system
1361 
1362    Output Parameters:
1363 -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
1364 
1365    Level: developer
1366 
1367    Notes:
1368 
1369 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1370 @*/
1371 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1372 {
1373   FETIDPMat_ctx  mat_ctx;
1374   PetscErrorCode ierr;
1375 
1376   PetscFunctionBegin;
1377   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1378   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
1379   PetscFunctionReturn(0);
1380 }
1381 /* -------------------------------------------------------------------------- */
1382 
1383 #undef __FUNCT__
1384 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
1385 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1386 {
1387   FETIDPMat_ctx  mat_ctx;
1388   PC_IS*         pcis;
1389   PC_BDDC*       pcbddc;
1390   PetscErrorCode ierr;
1391 
1392   PetscFunctionBegin;
1393   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1394   pcis = (PC_IS*)mat_ctx->pc->data;
1395   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1396 
1397   /* apply B_delta^T */
1398   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1399   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1400   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
1401   /* compute rhs for BDDC application */
1402   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1403   if (pcbddc->switch_static) {
1404     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1405   }
1406   /* apply BDDC */
1407   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
1408   /* put values into standard global vector */
1409   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1410   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1411   if (!pcbddc->switch_static) {
1412     /* compute values into the interior if solved for the partially subassembled Schur complement */
1413     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
1414     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
1415     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1416   }
1417   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1418   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1419   /* final change of basis if needed
1420      Is also sums the dirichlet part removed during RHS assembling */
1421   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
1422   PetscFunctionReturn(0);
1423 }
1424 
1425 #undef __FUNCT__
1426 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
1427 /*@
1428  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
1429 
1430    Collective
1431 
1432    Input Parameters:
1433 +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
1434 .  fetidp_flux_sol - the solution of the FETIDP linear system
1435 
1436    Output Parameters:
1437 -  standard_sol      - the solution defined on the physical domain
1438 
1439    Level: developer
1440 
1441    Notes:
1442 
1443 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1444 @*/
1445 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1446 {
1447   FETIDPMat_ctx  mat_ctx;
1448   PetscErrorCode ierr;
1449 
1450   PetscFunctionBegin;
1451   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1452   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
1453   PetscFunctionReturn(0);
1454 }
1455 /* -------------------------------------------------------------------------- */
1456 
1457 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1458 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1459 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1460 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1461 
1462 #undef __FUNCT__
1463 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
1464 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1465 {
1466 
1467   FETIDPMat_ctx  fetidpmat_ctx;
1468   Mat            newmat;
1469   FETIDPPC_ctx   fetidppc_ctx;
1470   PC             newpc;
1471   MPI_Comm       comm;
1472   PetscErrorCode ierr;
1473 
1474   PetscFunctionBegin;
1475   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1476   /* FETIDP linear matrix */
1477   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
1478   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
1479   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
1480   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1481   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
1482   ierr = MatSetUp(newmat);CHKERRQ(ierr);
1483   /* FETIDP preconditioner */
1484   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
1485   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
1486   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
1487   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
1488   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
1489   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1490   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
1491   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
1492   ierr = PCSetUp(newpc);CHKERRQ(ierr);
1493   /* return pointers for objects created */
1494   *fetidp_mat=newmat;
1495   *fetidp_pc=newpc;
1496   PetscFunctionReturn(0);
1497 }
1498 
1499 #undef __FUNCT__
1500 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
1501 /*@
1502  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
1503 
1504    Collective
1505 
1506    Input Parameters:
1507 +  pc - the BDDC preconditioning context already setup
1508 
1509    Output Parameters:
1510 .  fetidp_mat - shell FETIDP matrix object
1511 .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
1512 
1513    Options Database Keys:
1514 -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
1515 
1516    Level: developer
1517 
1518    Notes:
1519      Currently the only operation provided for FETIDP matrix is MatMult
1520 
1521 .seealso: PCBDDC
1522 @*/
1523 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1524 {
1525   PetscErrorCode ierr;
1526 
1527   PetscFunctionBegin;
1528   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1529   if (pc->setupcalled) {
1530     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1531   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1532   PetscFunctionReturn(0);
1533 }
1534 /* -------------------------------------------------------------------------- */
1535 /*MC
1536    PCBDDC - Balancing Domain Decomposition by Constraints.
1537 
1538    An implementation of the BDDC preconditioner based on
1539 
1540 .vb
1541    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
1542    [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
1543    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
1544 .ve
1545 
1546    The matrix to be preconditioned (Pmat) must be of type MATIS.
1547 
1548    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
1549 
1550    It also works with unsymmetric and indefinite problems.
1551 
1552    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.
1553 
1554    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
1555 
1556    Boundary nodes are split in vertices, edges and faces 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
1557 
1558    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace.
1559 
1560    Change of basis is performed similarly to [2] when requested. When more the 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.
1561 
1562    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
1563 
1564    Options Database Keys:
1565 
1566 .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
1567 .    -pc_bddc_use_edges <1> - use or not edges in primal space
1568 .    -pc_bddc_use_faces <0> - use or not faces in primal space
1569 .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
1570 .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
1571 .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
1572 .    -pc_bddc_levels <0> - maximum number of levels for multilevel
1573 .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
1574 -    -pc_bddc_check_level <0> - set verbosity level of debugging output
1575 
1576    Options for Dirichlet, Neumann or coarse solver can be set with
1577 .vb
1578       -pc_bddc_dirichlet_
1579       -pc_bddc_neumann_
1580       -pc_bddc_coarse_
1581 .ve
1582    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
1583 
1584    When using a multilevel approach, solvers' options at the N-th level can be specified as
1585 .vb
1586       -pc_bddc_dirichlet_N_
1587       -pc_bddc_neumann_N_
1588       -pc_bddc_coarse_N_
1589 .ve
1590    Note that level number ranges from the finest 0 to the coarsest N
1591 
1592    Level: intermediate
1593 
1594    Developer notes:
1595      Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose
1596 
1597      New deluxe scaling operator will be available soon.
1598 
1599    Contributed by Stefano Zampini
1600 
1601 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1602 M*/
1603 
1604 #undef __FUNCT__
1605 #define __FUNCT__ "PCCreate_BDDC"
1606 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1607 {
1608   PetscErrorCode      ierr;
1609   PC_BDDC             *pcbddc;
1610 
1611   PetscFunctionBegin;
1612   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1613   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1614   pc->data  = (void*)pcbddc;
1615 
1616   /* create PCIS data structure */
1617   ierr = PCISCreate(pc);CHKERRQ(ierr);
1618 
1619   /* BDDC customization */
1620   pcbddc->use_vertices        = PETSC_TRUE;
1621   pcbddc->use_edges           = PETSC_TRUE;
1622   pcbddc->use_faces           = PETSC_FALSE;
1623   pcbddc->use_change_of_basis = PETSC_FALSE;
1624   pcbddc->use_change_on_faces = PETSC_FALSE;
1625   pcbddc->switch_static       = PETSC_FALSE;
1626   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
1627   pcbddc->dbg_flag            = 0;
1628 
1629   pcbddc->BtoNmap                    = 0;
1630   pcbddc->local_primal_size          = 0;
1631   pcbddc->n_vertices                 = 0;
1632   pcbddc->n_actual_vertices          = 0;
1633   pcbddc->n_constraints              = 0;
1634   pcbddc->primal_indices_local_idxs  = 0;
1635   pcbddc->recompute_topography       = PETSC_FALSE;
1636   pcbddc->coarse_size                = 0;
1637   pcbddc->new_primal_space           = PETSC_FALSE;
1638   pcbddc->new_primal_space_local     = PETSC_FALSE;
1639   pcbddc->global_primal_indices      = 0;
1640   pcbddc->onearnullspace             = 0;
1641   pcbddc->onearnullvecs_state        = 0;
1642   pcbddc->user_primal_vertices       = 0;
1643   pcbddc->NullSpace                  = 0;
1644   pcbddc->temp_solution              = 0;
1645   pcbddc->original_rhs               = 0;
1646   pcbddc->local_mat                  = 0;
1647   pcbddc->ChangeOfBasisMatrix        = 0;
1648   pcbddc->coarse_vec                 = 0;
1649   pcbddc->coarse_rhs                 = 0;
1650   pcbddc->coarse_ksp                 = 0;
1651   pcbddc->coarse_phi_B               = 0;
1652   pcbddc->coarse_phi_D               = 0;
1653   pcbddc->coarse_psi_B               = 0;
1654   pcbddc->coarse_psi_D               = 0;
1655   pcbddc->vec1_P                     = 0;
1656   pcbddc->vec1_R                     = 0;
1657   pcbddc->vec2_R                     = 0;
1658   pcbddc->local_auxmat1              = 0;
1659   pcbddc->local_auxmat2              = 0;
1660   pcbddc->R_to_B                     = 0;
1661   pcbddc->R_to_D                     = 0;
1662   pcbddc->ksp_D                      = 0;
1663   pcbddc->ksp_R                      = 0;
1664   pcbddc->NeumannBoundaries          = 0;
1665   pcbddc->NeumannBoundariesLocal     = 0;
1666   pcbddc->DirichletBoundaries        = 0;
1667   pcbddc->DirichletBoundariesLocal   = 0;
1668   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
1669   pcbddc->n_ISForDofs                = 0;
1670   pcbddc->n_ISForDofsLocal           = 0;
1671   pcbddc->ISForDofs                  = 0;
1672   pcbddc->ISForDofsLocal             = 0;
1673   pcbddc->ConstraintMatrix           = 0;
1674   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
1675   pcbddc->coarse_loc_to_glob         = 0;
1676   pcbddc->coarsening_ratio           = 8;
1677   pcbddc->current_level              = 0;
1678   pcbddc->max_levels                 = 0;
1679 
1680   /* create local graph structure */
1681   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1682 
1683   /* scaling */
1684   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1685   pcbddc->work_scaling               = 0;
1686 
1687   /* function pointers */
1688   pc->ops->apply               = PCApply_BDDC;
1689   pc->ops->applytranspose      = 0;
1690   pc->ops->setup               = PCSetUp_BDDC;
1691   pc->ops->destroy             = PCDestroy_BDDC;
1692   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1693   pc->ops->view                = 0;
1694   pc->ops->applyrichardson     = 0;
1695   pc->ops->applysymmetricleft  = 0;
1696   pc->ops->applysymmetricright = 0;
1697   pc->ops->presolve            = PCPreSolve_BDDC;
1698   pc->ops->postsolve           = PCPostSolve_BDDC;
1699 
1700   /* composing function */
1701   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1702   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
1703   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1704   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
1705   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1706   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1707   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1708   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1709   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1710   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1711   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1712   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1713   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1714   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1715   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1716   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
1717   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1718   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1719   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1720   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1721   PetscFunctionReturn(0);
1722 }
1723 
1724