xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 9881197aa671a785a78c2af771fc6e38049d9b71)
1 /* TODOLIST
2    DofSplitting and DM attached to pc?
3    Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet)
4    change how to deal with the coarse problem (PCBDDCSetCoarseEnvironment):
5      - simplify coarse problem structure -> PCBDDC or PCREDUDANT, nothing else -> same comm for all levels?
6      - remove coarse enums and allow use of PCBDDCGetCoarseKSP
7      - remove metis dependency -> use MatPartitioning for multilevel -> Assemble serial adjacency in PCBDDCAnalyzeInterface?
8    code refactoring:
9      - pick up better names for static functions
10    change options structure:
11      - insert BDDC into MG framework?
12    provide other ops? Ask to developers
13    remove all unused printf
14    man pages
15 */
16 
17 /* ----------------------------------------------------------------------------------------------------------------------------------------------
18    Implementation of BDDC preconditioner based on:
19    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
20    ---------------------------------------------------------------------------------------------------------------------------------------------- */
21 
22 #include "bddc.h" /*I "petscpc.h" I*/  /* includes for fortran wrappers */
23 #include "bddcprivate.h"
24 #include <petscblaslapack.h>
25 
26 /* -------------------------------------------------------------------------- */
27 #undef __FUNCT__
28 #define __FUNCT__ "PCSetFromOptions_BDDC"
29 PetscErrorCode PCSetFromOptions_BDDC(PC pc)
30 {
31   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
32   PetscErrorCode ierr;
33 
34   PetscFunctionBegin;
35   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
36   /* Verbose debugging */
37   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
38   /* Primal space cumstomization */
39   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
40   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
41   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
42   /* Change of basis */
43   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);
44   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);
45   if (!pcbddc->use_change_of_basis) {
46     pcbddc->use_change_on_faces = PETSC_FALSE;
47   }
48   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
49   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);
50   /* Coarse solver context */
51   static const char * const avail_coarse_problems[] = {"sequential","replicated","parallel","multilevel","CoarseProblemType","PC_BDDC_",0}; /*order of choiches depends on ENUM defined in bddc.h */
52   ierr = PetscOptionsEnum("-pc_bddc_coarse_problem_type","Set coarse problem type","none",avail_coarse_problems,(PetscEnum)pcbddc->coarse_problem_type,(PetscEnum*)&pcbddc->coarse_problem_type,NULL);CHKERRQ(ierr);
53   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
54   ierr = PetscOptionsInt("-pc_bddc_max_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
55   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
56   ierr = PetscOptionsTail();CHKERRQ(ierr);
57   PetscFunctionReturn(0);
58 }
59 /* -------------------------------------------------------------------------- */
60 #undef __FUNCT__
61 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
62 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
63 {
64   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
65   PetscErrorCode ierr;
66 
67   PetscFunctionBegin;
68   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
69   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
70   pcbddc->user_primal_vertices = PrimalVertices;
71   PetscFunctionReturn(0);
72 }
73 #undef __FUNCT__
74 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
75 /*@
76  PCBDDCSetPrimalVerticesLocalIS - Set user defined primal vertices in PCBDDC.
77 
78    Not collective
79 
80    Input Parameters:
81 +  pc - the preconditioning context
82 -  PrimalVertices - index sets of primal vertices in local numbering
83 
84    Level: intermediate
85 
86    Notes:
87 
88 .seealso: PCBDDC
89 @*/
90 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
91 {
92   PetscErrorCode ierr;
93 
94   PetscFunctionBegin;
95   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
96   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
97   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
98   PetscFunctionReturn(0);
99 }
100 /* -------------------------------------------------------------------------- */
101 #undef __FUNCT__
102 #define __FUNCT__ "PCBDDCSetCoarseProblemType_BDDC"
103 static PetscErrorCode PCBDDCSetCoarseProblemType_BDDC(PC pc, CoarseProblemType CPT)
104 {
105   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
106 
107   PetscFunctionBegin;
108   pcbddc->coarse_problem_type = CPT;
109   PetscFunctionReturn(0);
110 }
111 
112 #undef __FUNCT__
113 #define __FUNCT__ "PCBDDCSetCoarseProblemType"
114 /*@
115  PCBDDCSetCoarseProblemType - Set coarse problem type in PCBDDC.
116 
117    Not collective
118 
119    Input Parameters:
120 +  pc - the preconditioning context
121 -  CoarseProblemType - pick a better name and explain what this is
122 
123    Level: intermediate
124 
125    Notes:
126    Not collective but all procs must call with same arguments.
127 
128 .seealso: PCBDDC
129 @*/
130 PetscErrorCode PCBDDCSetCoarseProblemType(PC pc, CoarseProblemType CPT)
131 {
132   PetscErrorCode ierr;
133 
134   PetscFunctionBegin;
135   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
136   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseProblemType_C",(PC,CoarseProblemType),(pc,CPT));CHKERRQ(ierr);
137   PetscFunctionReturn(0);
138 }
139 /* -------------------------------------------------------------------------- */
140 #undef __FUNCT__
141 #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
142 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
143 {
144   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
145 
146   PetscFunctionBegin;
147   pcbddc->coarsening_ratio=k;
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "PCBDDCSetCoarseningRatio"
153 /*@
154  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel coarsening
155 
156    Logically collective on PC
157 
158    Input Parameters:
159 +  pc - the preconditioning context
160 -  k - coarsening ratio
161 
162    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level.
163 
164    Level: intermediate
165 
166    Notes:
167 
168 .seealso: PCBDDC
169 @*/
170 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
171 {
172   PetscErrorCode ierr;
173 
174   PetscFunctionBegin;
175   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
176   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
177   PetscFunctionReturn(0);
178 }
179 /* -------------------------------------------------------------------------- */
180 
181 #undef __FUNCT__
182 #define __FUNCT__ "PCBDDCSetMaxLevels_BDDC"
183 static PetscErrorCode PCBDDCSetMaxLevels_BDDC(PC pc,PetscInt max_levels)
184 {
185   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
186 
187   PetscFunctionBegin;
188   pcbddc->max_levels=max_levels;
189   PetscFunctionReturn(0);
190 }
191 
192 #undef __FUNCT__
193 #define __FUNCT__ "PCBDDCSetMaxLevels"
194 /*@
195  PCBDDCSetMaxLevels - Sets the maximum number of levels within the multilevel approach.
196 
197    Logically collective on PC
198 
199    Input Parameters:
200 +  pc - the preconditioning context
201 -  max_levels - the maximum number of levels
202 
203    Default value is 1, i.e. coarse problem will be solved inexactly with one application
204    of PCBDDC preconditioner if the multilevel approach is requested.
205 
206    Level: intermediate
207 
208    Notes:
209 
210 .seealso: PCBDDC
211 @*/
212 PetscErrorCode PCBDDCSetMaxLevels(PC pc,PetscInt max_levels)
213 {
214   PetscErrorCode ierr;
215 
216   PetscFunctionBegin;
217   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
218   ierr = PetscTryMethod(pc,"PCBDDCSetMaxLevels_C",(PC,PetscInt),(pc,max_levels));CHKERRQ(ierr);
219   PetscFunctionReturn(0);
220 }
221 /* -------------------------------------------------------------------------- */
222 
223 #undef __FUNCT__
224 #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
225 static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
226 {
227   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
228   PetscErrorCode ierr;
229 
230   PetscFunctionBegin;
231   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
232   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
233   pcbddc->NullSpace=NullSpace;
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "PCBDDCSetNullSpace"
239 /*@
240  PCBDDCSetNullSpace - Set NullSpace of global operator of BDDC preconditioned mat.
241 
242    Logically collective on PC and MatNullSpace
243 
244    Input Parameters:
245 +  pc - the preconditioning context
246 -  NullSpace - Null space of the linear operator to be preconditioned.
247 
248    Level: intermediate
249 
250    Notes:
251 
252 .seealso: PCBDDC
253 @*/
254 PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
255 {
256   PetscErrorCode ierr;
257 
258   PetscFunctionBegin;
259   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
260   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
261   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
262   PetscFunctionReturn(0);
263 }
264 /* -------------------------------------------------------------------------- */
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
268 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
269 {
270   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
271   PetscErrorCode ierr;
272 
273   PetscFunctionBegin;
274   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
275   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
276   pcbddc->DirichletBoundaries=DirichletBoundaries;
277   PetscFunctionReturn(0);
278 }
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
282 /*@
283  PCBDDCSetDirichletBoundaries - Set index set defining subdomain part (in local ordering)
284                               of Dirichlet boundaries for the global problem.
285 
286    Not collective
287 
288    Input Parameters:
289 +  pc - the preconditioning context
290 -  DirichletBoundaries - sequential index set defining the subdomain part of Dirichlet boundaries (can be NULL)
291 
292    Level: intermediate
293 
294    Notes:
295 
296 .seealso: PCBDDC
297 @*/
298 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
299 {
300   PetscErrorCode ierr;
301 
302   PetscFunctionBegin;
303   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
304   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
305   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
306   PetscFunctionReturn(0);
307 }
308 /* -------------------------------------------------------------------------- */
309 
310 #undef __FUNCT__
311 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
312 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
313 {
314   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
315   PetscErrorCode ierr;
316 
317   PetscFunctionBegin;
318   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
319   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
320   pcbddc->NeumannBoundaries=NeumannBoundaries;
321   PetscFunctionReturn(0);
322 }
323 
324 #undef __FUNCT__
325 #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
326 /*@
327  PCBDDCSetNeumannBoundaries - Set index set defining subdomain part (in local ordering)
328                               of Neumann boundaries for the global problem.
329 
330    Not collective
331 
332    Input Parameters:
333 +  pc - the preconditioning context
334 -  NeumannBoundaries - sequential index set defining the subdomain part of Neumann boundaries (can be NULL)
335 
336    Level: intermediate
337 
338    Notes:
339 
340 .seealso: PCBDDC
341 @*/
342 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
343 {
344   PetscErrorCode ierr;
345 
346   PetscFunctionBegin;
347   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
348   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
349   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
350   PetscFunctionReturn(0);
351 }
352 /* -------------------------------------------------------------------------- */
353 
354 #undef __FUNCT__
355 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
356 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
357 {
358   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
359 
360   PetscFunctionBegin;
361   *DirichletBoundaries = pcbddc->DirichletBoundaries;
362   PetscFunctionReturn(0);
363 }
364 
365 #undef __FUNCT__
366 #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
367 /*@
368  PCBDDCGetDirichletBoundaries - Get index set defining subdomain part (in local ordering)
369                                 of Dirichlet boundaries for the global problem.
370 
371    Not collective
372 
373    Input Parameters:
374 +  pc - the preconditioning context
375 
376    Output Parameters:
377 +  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
378 
379    Level: intermediate
380 
381    Notes:
382 
383 .seealso: PCBDDC
384 @*/
385 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
386 {
387   PetscErrorCode ierr;
388 
389   PetscFunctionBegin;
390   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
391   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
392   PetscFunctionReturn(0);
393 }
394 /* -------------------------------------------------------------------------- */
395 
396 #undef __FUNCT__
397 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
398 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
399 {
400   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
401 
402   PetscFunctionBegin;
403   *NeumannBoundaries = pcbddc->NeumannBoundaries;
404   PetscFunctionReturn(0);
405 }
406 
407 #undef __FUNCT__
408 #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
409 /*@
410  PCBDDCGetNeumannBoundaries - Get index set defining subdomain part (in local ordering)
411                               of Neumann boundaries for the global problem.
412 
413    Not collective
414 
415    Input Parameters:
416 +  pc - the preconditioning context
417 
418    Output Parameters:
419 +  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
420 
421    Level: intermediate
422 
423    Notes:
424 
425 .seealso: PCBDDC
426 @*/
427 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
428 {
429   PetscErrorCode ierr;
430 
431   PetscFunctionBegin;
432   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
433   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
434   PetscFunctionReturn(0);
435 }
436 /* -------------------------------------------------------------------------- */
437 
438 #undef __FUNCT__
439 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
440 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
441 {
442   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
443   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
444   PetscErrorCode ierr;
445 
446   PetscFunctionBegin;
447   /* free old CSR */
448   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
449   /* TODO: PCBDDCGraphSetAdjacency */
450   /* get CSR into graph structure */
451   if (copymode == PETSC_COPY_VALUES) {
452     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
453     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
454     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
455     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
456   } else if (copymode == PETSC_OWN_POINTER) {
457     mat_graph->xadj = (PetscInt*)xadj;
458     mat_graph->adjncy = (PetscInt*)adjncy;
459   }
460   mat_graph->nvtxs_csr = nvtxs;
461   PetscFunctionReturn(0);
462 }
463 
464 #undef __FUNCT__
465 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
466 /*@
467  PCBDDCSetLocalAdjacencyGraph - Set CSR graph of local matrix for use of PCBDDC.
468 
469    Not collective
470 
471    Input Parameters:
472 +  pc - the preconditioning context
473 -  nvtxs - number of local vertices of the graph
474 -  xadj, adjncy - the CSR graph
475 -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. In the former case the user must free the array passed in;
476                                                              in the latter case, memory must be obtained with PetscMalloc.
477 
478    Level: intermediate
479 
480    Notes:
481 
482 .seealso: PCBDDC
483 @*/
484 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
485 {
486   void (*f)(void) = 0;
487   PetscErrorCode ierr;
488 
489   PetscFunctionBegin;
490   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
491   PetscValidIntPointer(xadj,3);
492   PetscValidIntPointer(xadj,4);
493   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
494     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
495   }
496   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
497   /* free arrays if PCBDDC is not the PC type */
498   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
499   if (!f && copymode == PETSC_OWN_POINTER) {
500     ierr = PetscFree(xadj);CHKERRQ(ierr);
501     ierr = PetscFree(adjncy);CHKERRQ(ierr);
502   }
503   PetscFunctionReturn(0);
504 }
505 /* -------------------------------------------------------------------------- */
506 
507 #undef __FUNCT__
508 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
509 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
510 {
511   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
512   PetscInt i;
513   PetscErrorCode ierr;
514 
515   PetscFunctionBegin;
516   /* Destroy ISes if they were already set */
517   for (i=0;i<pcbddc->n_ISForDofs;i++) {
518     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
519   }
520   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
521   /* allocate space then set */
522   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
523   for (i=0;i<n_is;i++) {
524     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
525     pcbddc->ISForDofs[i]=ISForDofs[i];
526   }
527   pcbddc->n_ISForDofs=n_is;
528   PetscFunctionReturn(0);
529 }
530 
531 #undef __FUNCT__
532 #define __FUNCT__ "PCBDDCSetDofsSplitting"
533 /*@
534  PCBDDCSetDofsSplitting - Set index sets defining fields of local mat.
535 
536    Not collective
537 
538    Input Parameters:
539 +  pc - the preconditioning context
540 -  n - number of index sets defining the fields
541 -  IS[] - array of IS describing the fields
542 
543    Level: intermediate
544 
545    Notes:
546 
547 .seealso: PCBDDC
548 @*/
549 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
550 {
551   PetscErrorCode ierr;
552 
553   PetscFunctionBegin;
554   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
555   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
556   PetscFunctionReturn(0);
557 }
558 /* -------------------------------------------------------------------------- */
559 #undef __FUNCT__
560 #define __FUNCT__ "PCPreSolve_BDDC"
561 /* -------------------------------------------------------------------------- */
562 /*
563    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
564                      guess if a transformation of basis approach has been selected.
565 
566    Input Parameter:
567 +  pc - the preconditioner contex
568 
569    Application Interface Routine: PCPreSolve()
570 
571    Notes:
572    The interface routine PCPreSolve() is not usually called directly by
573    the user, but instead is called by KSPSolve().
574 */
575 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
576 {
577   PetscErrorCode ierr;
578   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
579   PC_IS          *pcis = (PC_IS*)(pc->data);
580   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
581   Mat            temp_mat;
582   IS             dirIS;
583   PetscInt       dirsize,i,*is_indices;
584   PetscScalar    *array_x,*array_diagonal;
585   Vec            used_vec;
586   PetscBool      guess_nonzero;
587 
588   PetscFunctionBegin;
589   /* Creates parallel work vectors used in presolve. */
590   if (!pcbddc->original_rhs) {
591     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
592   }
593   if (!pcbddc->temp_solution) {
594     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
595   }
596   if (x) {
597     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
598     used_vec = x;
599   } else {
600     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
601     used_vec = pcbddc->temp_solution;
602     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
603   }
604   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
605   if (ksp) {
606     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
607     if ( !guess_nonzero ) {
608       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
609     }
610   }
611 
612   if (rhs) { /* TODO: wiser handling of rhs removal, which is only needed in case of zeroed rows */
613     /* store the original rhs */
614     ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
615 
616     /* Take into account zeroed rows -> change rhs and store solution removed */
617     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
618     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
619     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
620     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
621     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
622     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
623     ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr);
624     if (dirIS) {
625       ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr);
626       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
627       ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
628       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
629       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
630       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
631       ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
632       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
633     }
634     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
635     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
636 
637     /* remove the computed solution from the rhs */
638     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
639     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
640     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
641   }
642 
643   /* store partially computed solution and set initial guess */
644   if (x) {
645     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
646     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
647     if (pcbddc->use_exact_dirichlet && !pcbddc->coarse_psi_B) {
648       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
649       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
650       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
651       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
652       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
653       if (ksp) {
654         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
655       }
656     }
657   }
658 
659   if (pcbddc->use_change_of_basis) {
660     /* swap pointers for local matrices */
661     temp_mat = matis->A;
662     matis->A = pcbddc->local_mat;
663     pcbddc->local_mat = temp_mat;
664   }
665   if (pcbddc->use_change_of_basis && rhs) {
666     /* Get local rhs and apply transformation of basis */
667     ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
668     ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
669     /* from original basis to modified basis */
670     ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
671     /* put back modified values into the global vec using INSERT_VALUES copy mode */
672     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
673     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
674   }
675   if (ksp && pcbddc->NullSpace) {
676     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
677     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
678   }
679   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
680   PetscFunctionReturn(0);
681 }
682 /* -------------------------------------------------------------------------- */
683 #undef __FUNCT__
684 #define __FUNCT__ "PCPostSolve_BDDC"
685 /* -------------------------------------------------------------------------- */
686 /*
687    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
688                      approach has been selected. Also, restores rhs to its original state.
689 
690    Input Parameter:
691 +  pc - the preconditioner contex
692 
693    Application Interface Routine: PCPostSolve()
694 
695    Notes:
696    The interface routine PCPostSolve() is not usually called directly by
697    the user, but instead is called by KSPSolve().
698 */
699 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
700 {
701   PetscErrorCode ierr;
702   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
703   PC_IS          *pcis   = (PC_IS*)(pc->data);
704   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
705   Mat            temp_mat;
706 
707   PetscFunctionBegin;
708   if (pcbddc->use_change_of_basis) {
709     /* swap pointers for local matrices */
710     temp_mat = matis->A;
711     matis->A = pcbddc->local_mat;
712     pcbddc->local_mat = temp_mat;
713   }
714   if (pcbddc->use_change_of_basis && x) {
715     /* Get Local boundary and apply transformation of basis to solution vector */
716     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
717     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
718     /* from modified basis to original basis */
719     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
720     /* put back modified values into the global vec using INSERT_VALUES copy mode */
721     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
722     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
723   }
724   /* add solution removed in presolve */
725   if (x) {
726     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
727   }
728   /* restore rhs to its original state */
729   if (rhs) {
730     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
731   }
732   PetscFunctionReturn(0);
733 }
734 /* -------------------------------------------------------------------------- */
735 #undef __FUNCT__
736 #define __FUNCT__ "PCSetUp_BDDC"
737 /* -------------------------------------------------------------------------- */
738 /*
739    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
740                   by setting data structures and options.
741 
742    Input Parameter:
743 +  pc - the preconditioner context
744 
745    Application Interface Routine: PCSetUp()
746 
747    Notes:
748    The interface routine PCSetUp() is not usually called directly by
749    the user, but instead is called by PCApply() if necessary.
750 */
751 PetscErrorCode PCSetUp_BDDC(PC pc)
752 {
753   PetscErrorCode ierr;
754   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
755   MatStructure   flag;
756   PetscBool      computeis,computetopography,computesolvers;
757 
758   PetscFunctionBegin;
759   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */
760   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
761      So, we set to pcnone the Neumann problem of pcis in order to avoid unneeded computation
762      Also, we decide to directly build the (same) Dirichlet problem */
763   ierr = PetscOptionsSetValue("-is_localN_pc_type","none");CHKERRQ(ierr);
764   ierr = PetscOptionsSetValue("-is_localD_pc_type","none");CHKERRQ(ierr);
765   /* Get stdout for dbg */
766   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
767     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
768     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
769   }
770   /* first attempt to split work */
771   if (pc->setupcalled) {
772     computeis = PETSC_FALSE;
773     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
774     if (flag == SAME_PRECONDITIONER) {
775       computetopography = PETSC_FALSE;
776       computesolvers = PETSC_FALSE;
777     } else if (flag == SAME_NONZERO_PATTERN) {
778       computetopography = PETSC_FALSE;
779       computesolvers = PETSC_TRUE;
780     } else { /* DIFFERENT_NONZERO_PATTERN */
781       computetopography = PETSC_TRUE;
782       computesolvers = PETSC_TRUE;
783     }
784   } else {
785     computeis = PETSC_TRUE;
786     computetopography = PETSC_TRUE;
787     computesolvers = PETSC_TRUE;
788   }
789   /* Set up all the "iterative substructuring" common block */
790   if (computeis) {
791     ierr = PCISSetUp(pc);CHKERRQ(ierr);
792   }
793   /* Analyze interface and set up local constraint and change of basis matrices */
794   if (computetopography) {
795     /* reset data */
796     ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
797     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
798     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
799   }
800   if (computesolvers) {
801     /* reset data */
802     ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
803     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
804     /* Create coarse and local stuffs */
805     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
806     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
807   }
808   PetscFunctionReturn(0);
809 }
810 
811 /* -------------------------------------------------------------------------- */
812 /*
813    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
814 
815    Input Parameters:
816 .  pc - the preconditioner context
817 .  r - input vector (global)
818 
819    Output Parameter:
820 .  z - output vector (global)
821 
822    Application Interface Routine: PCApply()
823  */
824 #undef __FUNCT__
825 #define __FUNCT__ "PCApply_BDDC"
826 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
827 {
828   PC_IS             *pcis = (PC_IS*)(pc->data);
829   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
830   PetscErrorCode    ierr;
831   const PetscScalar one = 1.0;
832   const PetscScalar m_one = -1.0;
833   const PetscScalar zero = 0.0;
834 
835 /* This code is similar to that provided in nn.c for PCNN
836    NN interface preconditioner changed to BDDC
837    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
838 
839   PetscFunctionBegin;
840   if (!pcbddc->use_exact_dirichlet || pcbddc->coarse_psi_B) {
841     /* First Dirichlet solve */
842     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
843     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
844     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
845     /*
846       Assembling right hand side for BDDC operator
847       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
848       - pcis->vec1_B the interface part of the global vector z
849     */
850     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
851     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
852     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
853     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
854     ierr = VecCopy(r,z);CHKERRQ(ierr);
855     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
856     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
857     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
858   } else {
859     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
860     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
861     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
862   }
863 
864   /* Apply interface preconditioner
865      input/output vecs: pcis->vec1_B and pcis->vec1_D */
866   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
867 
868   /* Apply transpose of partition of unity operator */
869   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
870 
871   /* Second Dirichlet solve and assembling of output */
872   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
873   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
874   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
875   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
876   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr);
877   ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr);
878   if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
879   ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr);
880   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
881   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
882   PetscFunctionReturn(0);
883 }
884 /* -------------------------------------------------------------------------- */
885 
886 #undef __FUNCT__
887 #define __FUNCT__ "PCDestroy_BDDC"
888 PetscErrorCode PCDestroy_BDDC(PC pc)
889 {
890   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
891   PetscErrorCode ierr;
892 
893   PetscFunctionBegin;
894   /* free data created by PCIS */
895   ierr = PCISDestroy(pc);CHKERRQ(ierr);
896   /* free BDDC custom data  */
897   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
898   /* destroy objects related to topography */
899   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
900   /* free allocated graph structure */
901   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
902   /* free data for scaling operator */
903   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
904   /* free solvers stuff */
905   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
906   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
907   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
908   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
909   ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
910   /* free global vectors needed in presolve */
911   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
912   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
913   /* remove functions */
914   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
915   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
916   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetMaxLevels_C",NULL);CHKERRQ(ierr);
917   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
918   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
919   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
920   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
921   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
922   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseProblemType_C",NULL);CHKERRQ(ierr);
923   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
924   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
925   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
926   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
927   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
928   /* Free the private data structure */
929   ierr = PetscFree(pc->data);CHKERRQ(ierr);
930   PetscFunctionReturn(0);
931 }
932 /* -------------------------------------------------------------------------- */
933 
934 #undef __FUNCT__
935 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
936 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
937 {
938   FETIDPMat_ctx  mat_ctx;
939   PC_IS*         pcis;
940   PC_BDDC*       pcbddc;
941   PetscErrorCode ierr;
942 
943   PetscFunctionBegin;
944   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
945   pcis = (PC_IS*)mat_ctx->pc->data;
946   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
947 
948   /* change of basis for physical rhs if needed
949      It also changes the rhs in case of dirichlet boundaries */
950   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
951   /* store vectors for computation of fetidp final solution */
952   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
953   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
954   /* scale rhs since it should be unassembled */
955   /* TODO use counter scaling? (also below) */
956   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
957   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
958   /* Apply partition of unity */
959   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
960   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
961   if (!pcbddc->switch_static) {
962     /* compute partially subassembled Schur complement right-hand side */
963     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
964     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
965     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
966     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
967     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
968     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
969     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
970     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
971     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
972     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
973   }
974   /* BDDC rhs */
975   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
976   if (pcbddc->switch_static) {
977     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
978   }
979   /* apply BDDC */
980   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
981   /* Application of B_delta and assembling of rhs for fetidp fluxes */
982   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
983   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
984   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
985   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
986   /* restore original rhs */
987   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
988   PetscFunctionReturn(0);
989 }
990 
991 #undef __FUNCT__
992 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
993 /*@
994  PCBDDCMatFETIDPGetRHS - Get rhs for FETIDP linear system.
995 
996    Collective
997 
998    Input Parameters:
999 +  fetidp_mat   - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators
1000 +  standard_rhs - the rhs of your linear system
1001 
1002    Output Parameters:
1003 +  fetidp_flux_rhs   - the rhs of the FETIDP linear system
1004 
1005    Level: developer
1006 
1007    Notes:
1008 
1009 .seealso: PCBDDC
1010 @*/
1011 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1012 {
1013   FETIDPMat_ctx  mat_ctx;
1014   PetscErrorCode ierr;
1015 
1016   PetscFunctionBegin;
1017   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1018   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
1019   PetscFunctionReturn(0);
1020 }
1021 /* -------------------------------------------------------------------------- */
1022 
1023 #undef __FUNCT__
1024 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
1025 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1026 {
1027   FETIDPMat_ctx  mat_ctx;
1028   PC_IS*         pcis;
1029   PC_BDDC*       pcbddc;
1030   PetscErrorCode ierr;
1031 
1032   PetscFunctionBegin;
1033   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1034   pcis = (PC_IS*)mat_ctx->pc->data;
1035   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1036 
1037   /* apply B_delta^T */
1038   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1039   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1040   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
1041   /* compute rhs for BDDC application */
1042   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1043   if (pcbddc->switch_static) {
1044     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1045   }
1046   /* apply BDDC */
1047   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
1048   /* put values into standard global vector */
1049   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1050   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1051   if (!pcbddc->switch_static) {
1052     /* compute values into the interior if solved for the partially subassembled Schur complement */
1053     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
1054     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
1055     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1056   }
1057   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1058   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1059   /* final change of basis if needed
1060      Is also sums the dirichlet part removed during RHS assembling */
1061   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
1062   PetscFunctionReturn(0);
1063 }
1064 
1065 #undef __FUNCT__
1066 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
1067 /*@
1068  PCBDDCMatFETIDPGetSolution - Get Solution for FETIDP linear system.
1069 
1070    Collective
1071 
1072    Input Parameters:
1073 +  fetidp_mat        - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators
1074 +  fetidp_flux_sol - the solution of the FETIDP linear system
1075 
1076    Output Parameters:
1077 +  standard_sol      - the solution on the global domain
1078 
1079    Level: developer
1080 
1081    Notes:
1082 
1083 .seealso: PCBDDC
1084 @*/
1085 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1086 {
1087   FETIDPMat_ctx  mat_ctx;
1088   PetscErrorCode ierr;
1089 
1090   PetscFunctionBegin;
1091   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1092   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
1093   PetscFunctionReturn(0);
1094 }
1095 /* -------------------------------------------------------------------------- */
1096 
1097 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1098 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1099 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1100 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1101 
1102 #undef __FUNCT__
1103 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
1104 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1105 {
1106 
1107   FETIDPMat_ctx  fetidpmat_ctx;
1108   Mat            newmat;
1109   FETIDPPC_ctx   fetidppc_ctx;
1110   PC             newpc;
1111   MPI_Comm       comm;
1112   PetscErrorCode ierr;
1113 
1114   PetscFunctionBegin;
1115   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1116   /* FETIDP linear matrix */
1117   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
1118   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
1119   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
1120   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1121   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
1122   ierr = MatSetUp(newmat);CHKERRQ(ierr);
1123   /* FETIDP preconditioner */
1124   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
1125   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
1126   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
1127   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
1128   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
1129   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1130   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
1131   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
1132   ierr = PCSetUp(newpc);CHKERRQ(ierr);
1133   /* return pointers for objects created */
1134   *fetidp_mat=newmat;
1135   *fetidp_pc=newpc;
1136   PetscFunctionReturn(0);
1137 }
1138 
1139 #undef __FUNCT__
1140 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
1141 /*@
1142  PCBDDCCreateFETIDPOperators - Create operators for FETIDP.
1143 
1144    Collective
1145 
1146    Input Parameters:
1147 +  pc - the BDDC preconditioning context (setup must be already called)
1148 
1149    Level: developer
1150 
1151    Notes:
1152 
1153 .seealso: PCBDDC
1154 @*/
1155 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1156 {
1157   PetscErrorCode ierr;
1158 
1159   PetscFunctionBegin;
1160   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1161   if (pc->setupcalled) {
1162     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1163   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1164   PetscFunctionReturn(0);
1165 }
1166 /* -------------------------------------------------------------------------- */
1167 /*MC
1168    PCBDDC - Balancing Domain Decomposition by Constraints.
1169 
1170    Options Database Keys:
1171 .    -pcbddc ??? -
1172 
1173    Level: intermediate
1174 
1175    Notes: The matrix used with this preconditioner must be of type MATIS
1176 
1177           Unlike more 'conventional' interface preconditioners, this iterates over ALL the
1178           degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers
1179           on the subdomains).
1180 
1181           Options for the coarse grid preconditioner can be set with -
1182           Options for the Dirichlet subproblem can be set with -
1183           Options for the Neumann subproblem can be set with -
1184 
1185    Contributed by Stefano Zampini
1186 
1187 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1188 M*/
1189 
1190 #undef __FUNCT__
1191 #define __FUNCT__ "PCCreate_BDDC"
1192 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1193 {
1194   PetscErrorCode      ierr;
1195   PC_BDDC             *pcbddc;
1196 
1197   PetscFunctionBegin;
1198   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1199   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1200   pc->data  = (void*)pcbddc;
1201 
1202   /* create PCIS data structure */
1203   ierr = PCISCreate(pc);CHKERRQ(ierr);
1204 
1205   /* BDDC customization */
1206   pcbddc->use_vertices        = PETSC_TRUE;
1207   pcbddc->use_edges           = PETSC_TRUE;
1208   pcbddc->use_faces           = PETSC_FALSE;
1209   pcbddc->use_change_of_basis = PETSC_FALSE;
1210   pcbddc->use_change_on_faces = PETSC_FALSE;
1211   pcbddc->switch_static       = PETSC_FALSE;
1212   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
1213   pcbddc->dbg_flag            = 0;
1214 
1215   pcbddc->user_primal_vertices       = 0;
1216   pcbddc->NullSpace                  = 0;
1217   pcbddc->temp_solution              = 0;
1218   pcbddc->original_rhs               = 0;
1219   pcbddc->local_mat                  = 0;
1220   pcbddc->ChangeOfBasisMatrix        = 0;
1221   pcbddc->coarse_vec                 = 0;
1222   pcbddc->coarse_rhs                 = 0;
1223   pcbddc->coarse_ksp                 = 0;
1224   pcbddc->coarse_phi_B               = 0;
1225   pcbddc->coarse_phi_D               = 0;
1226   pcbddc->coarse_psi_B               = 0;
1227   pcbddc->coarse_psi_D               = 0;
1228   pcbddc->vec1_P                     = 0;
1229   pcbddc->vec1_R                     = 0;
1230   pcbddc->vec2_R                     = 0;
1231   pcbddc->local_auxmat1              = 0;
1232   pcbddc->local_auxmat2              = 0;
1233   pcbddc->R_to_B                     = 0;
1234   pcbddc->R_to_D                     = 0;
1235   pcbddc->ksp_D                      = 0;
1236   pcbddc->ksp_R                      = 0;
1237   pcbddc->NeumannBoundaries          = 0;
1238   pcbddc->ISForDofs                  = 0;
1239   pcbddc->ConstraintMatrix           = 0;
1240   pcbddc->use_exact_dirichlet        = PETSC_TRUE;
1241   pcbddc->coarse_loc_to_glob         = 0;
1242   pcbddc->coarsening_ratio           = 8;
1243   pcbddc->current_level              = 0;
1244   pcbddc->max_levels                 = 1;
1245 
1246   /* create local graph structure */
1247   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1248 
1249   /* scaling */
1250   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1251   pcbddc->work_scaling               = 0;
1252 
1253   /* function pointers */
1254   pc->ops->apply               = PCApply_BDDC;
1255   pc->ops->applytranspose      = 0;
1256   pc->ops->setup               = PCSetUp_BDDC;
1257   pc->ops->destroy             = PCDestroy_BDDC;
1258   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1259   pc->ops->view                = 0;
1260   pc->ops->applyrichardson     = 0;
1261   pc->ops->applysymmetricleft  = 0;
1262   pc->ops->applysymmetricright = 0;
1263   pc->ops->presolve            = PCPreSolve_BDDC;
1264   pc->ops->postsolve           = PCPostSolve_BDDC;
1265 
1266   /* composing function */
1267   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1268   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
1269   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetMaxLevels_C",PCBDDCSetMaxLevels_BDDC);CHKERRQ(ierr);
1270   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1271   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1272   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1273   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1274   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1275   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseProblemType_C",PCBDDCSetCoarseProblemType_BDDC);CHKERRQ(ierr);
1276   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1277   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1278   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1279   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1280   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1281   PetscFunctionReturn(0);
1282 }
1283 
1284