xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision fcd409f55db887dc8ef901cc2123e1395f10ef8d)
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;
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   if (rhs) { /* TODO: wiser handling of rhs removal, which is only needed in case of zeroed rows */
618     /* store the original rhs */
619     ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
620 
621     /* Take into account zeroed rows -> change rhs and store solution removed */
622     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
623     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
624     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
625     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
626     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
627     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
628     ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr);
629     if (dirIS) {
630       ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr);
631       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
632       ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
633       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
634       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
635       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
636       ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
637       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
638     }
639     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
640     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
641 
642     /* remove the computed solution from the rhs */
643     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
644     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
645     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
646   }
647 
648   /* store partially computed solution and set initial guess */
649   if (x) {
650     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
651     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
652     if (pcbddc->use_exact_dirichlet && !pcbddc->coarse_psi_B) {
653       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
654       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
655       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
656       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
657       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
658       if (ksp) {
659         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
660       }
661     }
662   }
663 
664   if (pcbddc->use_change_of_basis) {
665     /* swap pointers for local matrices */
666     temp_mat = matis->A;
667     matis->A = pcbddc->local_mat;
668     pcbddc->local_mat = temp_mat;
669   }
670   if (pcbddc->use_change_of_basis && rhs) {
671     /* Get local rhs and apply transformation of basis */
672     ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
673     ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
674     /* from original basis to modified basis */
675     ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
676     /* put back modified values into the global vec using INSERT_VALUES copy mode */
677     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
678     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
679   }
680   if (ksp && pcbddc->NullSpace) {
681     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
682     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
683   }
684   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
685   PetscFunctionReturn(0);
686 }
687 /* -------------------------------------------------------------------------- */
688 #undef __FUNCT__
689 #define __FUNCT__ "PCPostSolve_BDDC"
690 /* -------------------------------------------------------------------------- */
691 /*
692    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
693                      approach has been selected. Also, restores rhs to its original state.
694 
695    Input Parameter:
696 +  pc - the preconditioner contex
697 
698    Application Interface Routine: PCPostSolve()
699 
700    Notes:
701    The interface routine PCPostSolve() is not usually called directly by
702    the user, but instead is called by KSPSolve().
703 */
704 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
705 {
706   PetscErrorCode ierr;
707   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
708   PC_IS          *pcis   = (PC_IS*)(pc->data);
709   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
710   Mat            temp_mat;
711 
712   PetscFunctionBegin;
713   if (pcbddc->use_change_of_basis) {
714     /* swap pointers for local matrices */
715     temp_mat = matis->A;
716     matis->A = pcbddc->local_mat;
717     pcbddc->local_mat = temp_mat;
718   }
719   if (pcbddc->use_change_of_basis && x) {
720     /* Get Local boundary and apply transformation of basis to solution vector */
721     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
722     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
723     /* from modified basis to original basis */
724     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
725     /* put back modified values into the global vec using INSERT_VALUES copy mode */
726     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
727     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
728   }
729   /* add solution removed in presolve */
730   if (x) {
731     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
732   }
733   /* restore rhs to its original state */
734   if (rhs) {
735     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
736   }
737   PetscFunctionReturn(0);
738 }
739 /* -------------------------------------------------------------------------- */
740 #undef __FUNCT__
741 #define __FUNCT__ "PCSetUp_BDDC"
742 /* -------------------------------------------------------------------------- */
743 /*
744    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
745                   by setting data structures and options.
746 
747    Input Parameter:
748 +  pc - the preconditioner context
749 
750    Application Interface Routine: PCSetUp()
751 
752    Notes:
753    The interface routine PCSetUp() is not usually called directly by
754    the user, but instead is called by PCApply() if necessary.
755 */
756 PetscErrorCode PCSetUp_BDDC(PC pc)
757 {
758   PetscErrorCode ierr;
759   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
760   MatStructure   flag;
761   PetscBool      computeis,computetopography,computesolvers;
762 
763   PetscFunctionBegin;
764   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */
765   /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */
766   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
767      Also, BDDC directly build the Dirichlet problem */
768   /* Get stdout for dbg */
769   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
770     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
771     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
772     if (pcbddc->current_level) {
773       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
774     }
775   }
776   /* first attempt to split work */
777   if (pc->setupcalled) {
778     computeis = PETSC_FALSE;
779     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
780     if (flag == SAME_PRECONDITIONER) {
781       computetopography = PETSC_FALSE;
782       computesolvers = PETSC_FALSE;
783     } else if (flag == SAME_NONZERO_PATTERN) {
784       computetopography = PETSC_FALSE;
785       computesolvers = PETSC_TRUE;
786     } else { /* DIFFERENT_NONZERO_PATTERN */
787       computetopography = PETSC_TRUE;
788       computesolvers = PETSC_TRUE;
789     }
790   } else {
791     computeis = PETSC_TRUE;
792     computetopography = PETSC_TRUE;
793     computesolvers = PETSC_TRUE;
794   }
795   /* Set up all the "iterative substructuring" common block without computing solvers */
796   if (computeis) {
797     /* HACK INTO PCIS */
798     PC_IS* pcis = (PC_IS*)pc->data;
799     pcis->computesolvers = PETSC_FALSE;
800     ierr = PCISSetUp(pc);CHKERRQ(ierr);
801   }
802   /* Analyze interface and set up local constraint and change of basis matrices */
803   if (computetopography) {
804     /* reset data */
805     ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
806     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
807     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
808   }
809   if (computesolvers) {
810     /* reset data */
811     ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
812     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
813     /* Create coarse and local stuffs */
814     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
815     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
816   }
817   if (pcbddc->dbg_flag && pcbddc->current_level) {
818     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
819   }
820   PetscFunctionReturn(0);
821 }
822 
823 /* -------------------------------------------------------------------------- */
824 /*
825    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
826 
827    Input Parameters:
828 .  pc - the preconditioner context
829 .  r - input vector (global)
830 
831    Output Parameter:
832 .  z - output vector (global)
833 
834    Application Interface Routine: PCApply()
835  */
836 #undef __FUNCT__
837 #define __FUNCT__ "PCApply_BDDC"
838 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
839 {
840   PC_IS             *pcis = (PC_IS*)(pc->data);
841   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
842   PetscErrorCode    ierr;
843   const PetscScalar one = 1.0;
844   const PetscScalar m_one = -1.0;
845   const PetscScalar zero = 0.0;
846 
847 /* This code is similar to that provided in nn.c for PCNN
848    NN interface preconditioner changed to BDDC
849    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
850 
851   PetscFunctionBegin;
852   if (!pcbddc->use_exact_dirichlet || pcbddc->coarse_psi_B) {
853     /* First Dirichlet solve */
854     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
855     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
856     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
857     /*
858       Assembling right hand side for BDDC operator
859       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
860       - pcis->vec1_B the interface part of the global vector z
861     */
862     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
863     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
864     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
865     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
866     ierr = VecCopy(r,z);CHKERRQ(ierr);
867     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
868     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
869     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
870   } else {
871     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
872     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
873     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
874   }
875 
876   /* Apply interface preconditioner
877      input/output vecs: pcis->vec1_B and pcis->vec1_D */
878   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
879 
880   /* Apply transpose of partition of unity operator */
881   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
882 
883   /* Second Dirichlet solve and assembling of output */
884   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
885   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
886   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
887   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
888   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr);
889   ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr);
890   if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
891   ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr);
892   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
893   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
894   PetscFunctionReturn(0);
895 }
896 /* -------------------------------------------------------------------------- */
897 
898 #undef __FUNCT__
899 #define __FUNCT__ "PCDestroy_BDDC"
900 PetscErrorCode PCDestroy_BDDC(PC pc)
901 {
902   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
903   PetscErrorCode ierr;
904 
905   PetscFunctionBegin;
906   /* free data created by PCIS */
907   ierr = PCISDestroy(pc);CHKERRQ(ierr);
908   /* free BDDC custom data  */
909   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
910   /* destroy objects related to topography */
911   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
912   /* free allocated graph structure */
913   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
914   /* free data for scaling operator */
915   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
916   /* free solvers stuff */
917   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
918   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
919   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
920   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
921   ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
922   /* free global vectors needed in presolve */
923   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
924   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
925   /* remove functions */
926   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
927   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
928   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
929   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
930   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
931   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
932   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
933   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
934   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
935   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
936   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
937   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
938   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
939   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
940   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
941   /* Free the private data structure */
942   ierr = PetscFree(pc->data);CHKERRQ(ierr);
943   PetscFunctionReturn(0);
944 }
945 /* -------------------------------------------------------------------------- */
946 
947 #undef __FUNCT__
948 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
949 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
950 {
951   FETIDPMat_ctx  mat_ctx;
952   PC_IS*         pcis;
953   PC_BDDC*       pcbddc;
954   PetscErrorCode ierr;
955 
956   PetscFunctionBegin;
957   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
958   pcis = (PC_IS*)mat_ctx->pc->data;
959   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
960 
961   /* change of basis for physical rhs if needed
962      It also changes the rhs in case of dirichlet boundaries */
963   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
964   /* store vectors for computation of fetidp final solution */
965   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
966   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
967   /* scale rhs since it should be unassembled */
968   /* TODO use counter scaling? (also below) */
969   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
970   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
971   /* Apply partition of unity */
972   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
973   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
974   if (!pcbddc->switch_static) {
975     /* compute partially subassembled Schur complement right-hand side */
976     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
977     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
978     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
979     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
980     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
981     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
982     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
983     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
984     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
985     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
986   }
987   /* BDDC rhs */
988   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
989   if (pcbddc->switch_static) {
990     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
991   }
992   /* apply BDDC */
993   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
994   /* Application of B_delta and assembling of rhs for fetidp fluxes */
995   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
996   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
997   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
998   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
999   /* restore original rhs */
1000   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
1001   PetscFunctionReturn(0);
1002 }
1003 
1004 #undef __FUNCT__
1005 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
1006 /*@
1007  PCBDDCMatFETIDPGetRHS - Get rhs for FETIDP linear system.
1008 
1009    Collective
1010 
1011    Input Parameters:
1012 +  fetidp_mat   - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators
1013 +  standard_rhs - the rhs of your linear system
1014 
1015    Output Parameters:
1016 +  fetidp_flux_rhs   - the rhs of the FETIDP linear system
1017 
1018    Level: developer
1019 
1020    Notes:
1021 
1022 .seealso: PCBDDC
1023 @*/
1024 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1025 {
1026   FETIDPMat_ctx  mat_ctx;
1027   PetscErrorCode ierr;
1028 
1029   PetscFunctionBegin;
1030   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1031   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
1032   PetscFunctionReturn(0);
1033 }
1034 /* -------------------------------------------------------------------------- */
1035 
1036 #undef __FUNCT__
1037 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
1038 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1039 {
1040   FETIDPMat_ctx  mat_ctx;
1041   PC_IS*         pcis;
1042   PC_BDDC*       pcbddc;
1043   PetscErrorCode ierr;
1044 
1045   PetscFunctionBegin;
1046   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1047   pcis = (PC_IS*)mat_ctx->pc->data;
1048   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1049 
1050   /* apply B_delta^T */
1051   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1052   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1053   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
1054   /* compute rhs for BDDC application */
1055   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1056   if (pcbddc->switch_static) {
1057     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1058   }
1059   /* apply BDDC */
1060   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
1061   /* put values into standard global vector */
1062   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1063   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1064   if (!pcbddc->switch_static) {
1065     /* compute values into the interior if solved for the partially subassembled Schur complement */
1066     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
1067     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
1068     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1069   }
1070   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1071   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1072   /* final change of basis if needed
1073      Is also sums the dirichlet part removed during RHS assembling */
1074   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
1075   PetscFunctionReturn(0);
1076 }
1077 
1078 #undef __FUNCT__
1079 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
1080 /*@
1081  PCBDDCMatFETIDPGetSolution - Get Solution for FETIDP linear system.
1082 
1083    Collective
1084 
1085    Input Parameters:
1086 +  fetidp_mat        - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators
1087 +  fetidp_flux_sol - the solution of the FETIDP linear system
1088 
1089    Output Parameters:
1090 +  standard_sol      - the solution on the global domain
1091 
1092    Level: developer
1093 
1094    Notes:
1095 
1096 .seealso: PCBDDC
1097 @*/
1098 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1099 {
1100   FETIDPMat_ctx  mat_ctx;
1101   PetscErrorCode ierr;
1102 
1103   PetscFunctionBegin;
1104   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1105   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
1106   PetscFunctionReturn(0);
1107 }
1108 /* -------------------------------------------------------------------------- */
1109 
1110 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1111 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1112 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1113 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1114 
1115 #undef __FUNCT__
1116 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
1117 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1118 {
1119 
1120   FETIDPMat_ctx  fetidpmat_ctx;
1121   Mat            newmat;
1122   FETIDPPC_ctx   fetidppc_ctx;
1123   PC             newpc;
1124   MPI_Comm       comm;
1125   PetscErrorCode ierr;
1126 
1127   PetscFunctionBegin;
1128   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1129   /* FETIDP linear matrix */
1130   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
1131   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
1132   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
1133   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1134   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
1135   ierr = MatSetUp(newmat);CHKERRQ(ierr);
1136   /* FETIDP preconditioner */
1137   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
1138   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
1139   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
1140   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
1141   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
1142   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1143   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
1144   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
1145   ierr = PCSetUp(newpc);CHKERRQ(ierr);
1146   /* return pointers for objects created */
1147   *fetidp_mat=newmat;
1148   *fetidp_pc=newpc;
1149   PetscFunctionReturn(0);
1150 }
1151 
1152 #undef __FUNCT__
1153 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
1154 /*@
1155  PCBDDCCreateFETIDPOperators - Create operators for FETIDP.
1156 
1157    Collective
1158 
1159    Input Parameters:
1160 +  pc - the BDDC preconditioning context (setup must be already called)
1161 
1162    Level: developer
1163 
1164    Notes:
1165 
1166 .seealso: PCBDDC
1167 @*/
1168 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1169 {
1170   PetscErrorCode ierr;
1171 
1172   PetscFunctionBegin;
1173   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1174   if (pc->setupcalled) {
1175     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1176   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1177   PetscFunctionReturn(0);
1178 }
1179 /* -------------------------------------------------------------------------- */
1180 /*MC
1181    PCBDDC - Balancing Domain Decomposition by Constraints.
1182 
1183    Options Database Keys:
1184 .    -pcbddc ??? -
1185 
1186    Level: intermediate
1187 
1188    Notes: The matrix used with this preconditioner must be of type MATIS
1189 
1190           Unlike more 'conventional' interface preconditioners, this iterates over ALL the
1191           degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers
1192           on the subdomains).
1193 
1194           Options for the coarse grid preconditioner can be set with -
1195           Options for the Dirichlet subproblem can be set with -
1196           Options for the Neumann subproblem can be set with -
1197 
1198    Contributed by Stefano Zampini
1199 
1200 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1201 M*/
1202 
1203 #undef __FUNCT__
1204 #define __FUNCT__ "PCCreate_BDDC"
1205 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1206 {
1207   PetscErrorCode      ierr;
1208   PC_BDDC             *pcbddc;
1209 
1210   PetscFunctionBegin;
1211   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1212   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1213   pc->data  = (void*)pcbddc;
1214 
1215   /* create PCIS data structure */
1216   ierr = PCISCreate(pc);CHKERRQ(ierr);
1217 
1218   /* BDDC customization */
1219   pcbddc->use_vertices        = PETSC_TRUE;
1220   pcbddc->use_edges           = PETSC_TRUE;
1221   pcbddc->use_faces           = PETSC_FALSE;
1222   pcbddc->use_change_of_basis = PETSC_FALSE;
1223   pcbddc->use_change_on_faces = PETSC_FALSE;
1224   pcbddc->switch_static       = PETSC_FALSE;
1225   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
1226   pcbddc->dbg_flag            = 0;
1227 
1228   pcbddc->user_primal_vertices       = 0;
1229   pcbddc->NullSpace                  = 0;
1230   pcbddc->temp_solution              = 0;
1231   pcbddc->original_rhs               = 0;
1232   pcbddc->local_mat                  = 0;
1233   pcbddc->ChangeOfBasisMatrix        = 0;
1234   pcbddc->coarse_vec                 = 0;
1235   pcbddc->coarse_rhs                 = 0;
1236   pcbddc->coarse_ksp                 = 0;
1237   pcbddc->coarse_phi_B               = 0;
1238   pcbddc->coarse_phi_D               = 0;
1239   pcbddc->coarse_psi_B               = 0;
1240   pcbddc->coarse_psi_D               = 0;
1241   pcbddc->vec1_P                     = 0;
1242   pcbddc->vec1_R                     = 0;
1243   pcbddc->vec2_R                     = 0;
1244   pcbddc->local_auxmat1              = 0;
1245   pcbddc->local_auxmat2              = 0;
1246   pcbddc->R_to_B                     = 0;
1247   pcbddc->R_to_D                     = 0;
1248   pcbddc->ksp_D                      = 0;
1249   pcbddc->ksp_R                      = 0;
1250   pcbddc->NeumannBoundaries          = 0;
1251   pcbddc->ISForDofs                  = 0;
1252   pcbddc->ConstraintMatrix           = 0;
1253   pcbddc->use_exact_dirichlet        = PETSC_TRUE;
1254   pcbddc->coarse_loc_to_glob         = 0;
1255   pcbddc->coarsening_ratio           = 8;
1256   pcbddc->current_level              = 0;
1257   pcbddc->max_levels                 = 0;
1258 
1259   /* create local graph structure */
1260   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1261 
1262   /* scaling */
1263   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1264   pcbddc->work_scaling               = 0;
1265 
1266   /* function pointers */
1267   pc->ops->apply               = PCApply_BDDC;
1268   pc->ops->applytranspose      = 0;
1269   pc->ops->setup               = PCSetUp_BDDC;
1270   pc->ops->destroy             = PCDestroy_BDDC;
1271   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1272   pc->ops->view                = 0;
1273   pc->ops->applyrichardson     = 0;
1274   pc->ops->applysymmetricleft  = 0;
1275   pc->ops->applysymmetricright = 0;
1276   pc->ops->presolve            = PCPreSolve_BDDC;
1277   pc->ops->postsolve           = PCPostSolve_BDDC;
1278 
1279   /* composing function */
1280   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1281   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
1282   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1283   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
1284   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1285   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1286   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1287   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1288   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1289   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1290   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1291   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1292   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1293   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1294   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1295   PetscFunctionReturn(0);
1296 }
1297 
1298