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