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