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