xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 08a5cf49253972d9eb6d66cc143f7261a0e52c9d)
1 /* TODOLIST
2 
3    ConstraintsSetup
4    - tolerances for constraints as an option (take care of single precision!)
5    - Can MAT_IGNORE_ZERO_ENTRIES be used for Constraints Matrix?
6 
7    Solvers
8    - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers)
9    - Propagate ksp prefixes for solvers to mat objects?
10    - Propagate nearnullspace info among levels
11 
12    User interface
13    - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access)
14    - DofSplitting and DM attached to pc?
15 
16    Debugging output
17    - Better management of verbosity levels of debugging output
18 
19    Build
20    - make runexe59
21 
22    Extra
23    - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
24    - add support for computing h,H and related using coordinates?
25    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
26    - Better management in PCIS code
27    - BDDC with MG framework?
28 
29    FETIDP
30    - Move FETIDP code to its own classes
31 
32    MATIS related operations contained in BDDC code
33    - Provide general case for subassembling
34    - Preallocation routines in MatISGetMPIAXAIJ
35 
36 */
37 
38 /* ----------------------------------------------------------------------------------------------------------------------------------------------
39    Implementation of BDDC preconditioner based on:
40    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
41    ---------------------------------------------------------------------------------------------------------------------------------------------- */
42 
43 #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
44 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
45 #include <petscblaslapack.h>
46 
47 /* -------------------------------------------------------------------------- */
48 #undef __FUNCT__
49 #define __FUNCT__ "PCSetFromOptions_BDDC"
50 PetscErrorCode PCSetFromOptions_BDDC(PC pc)
51 {
52   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
53   PetscErrorCode ierr;
54 
55   PetscFunctionBegin;
56   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
57   /* Verbose debugging */
58   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
59   /* Primal space cumstomization */
60   ierr = PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);CHKERRQ(ierr);
61   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
62   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
63   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
64   /* Change of basis */
65   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);
66   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);
67   if (!pcbddc->use_change_of_basis) {
68     pcbddc->use_change_on_faces = PETSC_FALSE;
69   }
70   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
71   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);
72   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
73   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
74   ierr = PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);CHKERRQ(ierr);
75   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
76   ierr = PetscOptionsTail();CHKERRQ(ierr);
77   PetscFunctionReturn(0);
78 }
79 /* -------------------------------------------------------------------------- */
80 #undef __FUNCT__
81 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
82 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
83 {
84   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
85   PetscErrorCode ierr;
86 
87   PetscFunctionBegin;
88   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
89   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
90   pcbddc->user_primal_vertices = PrimalVertices;
91   PetscFunctionReturn(0);
92 }
93 #undef __FUNCT__
94 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
95 /*@
96  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
97 
98    Not collective
99 
100    Input Parameters:
101 +  pc - the preconditioning context
102 -  PrimalVertices - index set of primal vertices in local numbering
103 
104    Level: intermediate
105 
106    Notes:
107 
108 .seealso: PCBDDC
109 @*/
110 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
111 {
112   PetscErrorCode ierr;
113 
114   PetscFunctionBegin;
115   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
116   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
117   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
118   PetscFunctionReturn(0);
119 }
120 /* -------------------------------------------------------------------------- */
121 #undef __FUNCT__
122 #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
123 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
124 {
125   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
126 
127   PetscFunctionBegin;
128   pcbddc->coarsening_ratio = k;
129   PetscFunctionReturn(0);
130 }
131 
132 #undef __FUNCT__
133 #define __FUNCT__ "PCBDDCSetCoarseningRatio"
134 /*@
135  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
136 
137    Logically collective on PC
138 
139    Input Parameters:
140 +  pc - the preconditioning context
141 -  k - coarsening ratio (H/h at the coarser level)
142 
143    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
144 
145    Level: intermediate
146 
147    Notes:
148 
149 .seealso: PCBDDC
150 @*/
151 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
152 {
153   PetscErrorCode ierr;
154 
155   PetscFunctionBegin;
156   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
157   PetscValidLogicalCollectiveInt(pc,k,2);
158   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
163 #undef __FUNCT__
164 #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
165 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
166 {
167   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
168 
169   PetscFunctionBegin;
170   pcbddc->use_exact_dirichlet_trick = flg;
171   PetscFunctionReturn(0);
172 }
173 
174 #undef __FUNCT__
175 #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
176 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
177 {
178   PetscErrorCode ierr;
179 
180   PetscFunctionBegin;
181   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
182   PetscValidLogicalCollectiveBool(pc,flg,2);
183   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
184   PetscFunctionReturn(0);
185 }
186 
187 #undef __FUNCT__
188 #define __FUNCT__ "PCBDDCSetLevel_BDDC"
189 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
190 {
191   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
192 
193   PetscFunctionBegin;
194   pcbddc->current_level = level;
195   PetscFunctionReturn(0);
196 }
197 
198 #undef __FUNCT__
199 #define __FUNCT__ "PCBDDCSetLevel"
200 PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
201 {
202   PetscErrorCode ierr;
203 
204   PetscFunctionBegin;
205   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
206   PetscValidLogicalCollectiveInt(pc,level,2);
207   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "PCBDDCSetLevels_BDDC"
213 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
214 {
215   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
216 
217   PetscFunctionBegin;
218   pcbddc->max_levels = levels;
219   PetscFunctionReturn(0);
220 }
221 
222 #undef __FUNCT__
223 #define __FUNCT__ "PCBDDCSetLevels"
224 /*@
225  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
226 
227    Logically collective on PC
228 
229    Input Parameters:
230 +  pc - the preconditioning context
231 -  levels - the maximum number of levels (max 9)
232 
233    Default value is 0, i.e. traditional one-level BDDC
234 
235    Level: intermediate
236 
237    Notes:
238 
239 .seealso: PCBDDC
240 @*/
241 PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
242 {
243   PetscErrorCode ierr;
244 
245   PetscFunctionBegin;
246   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
247   PetscValidLogicalCollectiveInt(pc,levels,2);
248   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
249   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
250   PetscFunctionReturn(0);
251 }
252 /* -------------------------------------------------------------------------- */
253 
254 #undef __FUNCT__
255 #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
256 static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
257 {
258   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
259   PetscErrorCode ierr;
260 
261   PetscFunctionBegin;
262   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
263   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
264   pcbddc->NullSpace = NullSpace;
265   PetscFunctionReturn(0);
266 }
267 
268 #undef __FUNCT__
269 #define __FUNCT__ "PCBDDCSetNullSpace"
270 /*@
271  PCBDDCSetNullSpace - Set nullspace for BDDC operator
272 
273    Logically collective on PC and MatNullSpace
274 
275    Input Parameters:
276 +  pc - the preconditioning context
277 -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
278 
279    Level: intermediate
280 
281    Notes:
282 
283 .seealso: PCBDDC
284 @*/
285 PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
286 {
287   PetscErrorCode ierr;
288 
289   PetscFunctionBegin;
290   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
291   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
292   PetscCheckSameComm(pc,1,NullSpace,2);
293   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
294   PetscFunctionReturn(0);
295 }
296 /* -------------------------------------------------------------------------- */
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
300 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
301 {
302   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
303   PetscErrorCode ierr;
304 
305   PetscFunctionBegin;
306   /* last user setting takes precendence -> destroy any other customization */
307   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
308   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
309   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
310   pcbddc->DirichletBoundaries = DirichletBoundaries;
311   pcbddc->recompute_topography = PETSC_TRUE;
312   PetscFunctionReturn(0);
313 }
314 
315 #undef __FUNCT__
316 #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
317 /*@
318  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
319 
320    Collective
321 
322    Input Parameters:
323 +  pc - the preconditioning context
324 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
325 
326    Level: intermediate
327 
328    Notes: Any process can list any global node
329 
330 .seealso: PCBDDC
331 @*/
332 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
333 {
334   PetscErrorCode ierr;
335 
336   PetscFunctionBegin;
337   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
338   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
339   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
340   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
341   PetscFunctionReturn(0);
342 }
343 /* -------------------------------------------------------------------------- */
344 
345 #undef __FUNCT__
346 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
347 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
348 {
349   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
350   PetscErrorCode ierr;
351 
352   PetscFunctionBegin;
353   /* last user setting takes precendence -> destroy any other customization */
354   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
355   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
356   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
357   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
358   pcbddc->recompute_topography = PETSC_TRUE;
359   PetscFunctionReturn(0);
360 }
361 
362 #undef __FUNCT__
363 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
364 /*@
365  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
366 
367    Collective
368 
369    Input Parameters:
370 +  pc - the preconditioning context
371 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
372 
373    Level: intermediate
374 
375    Notes:
376 
377 .seealso: PCBDDC
378 @*/
379 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
380 {
381   PetscErrorCode ierr;
382 
383   PetscFunctionBegin;
384   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
385   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
386   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
387   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
388   PetscFunctionReturn(0);
389 }
390 /* -------------------------------------------------------------------------- */
391 
392 #undef __FUNCT__
393 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
394 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
395 {
396   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
397   PetscErrorCode ierr;
398 
399   PetscFunctionBegin;
400   /* last user setting takes precendence -> destroy any other customization */
401   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
402   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
403   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
404   pcbddc->NeumannBoundaries = NeumannBoundaries;
405   pcbddc->recompute_topography = PETSC_TRUE;
406   PetscFunctionReturn(0);
407 }
408 
409 #undef __FUNCT__
410 #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
411 /*@
412  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
413 
414    Collective
415 
416    Input Parameters:
417 +  pc - the preconditioning context
418 -  NeumannBoundaries - parallel IS defining the Neumann boundaries
419 
420    Level: intermediate
421 
422    Notes: Any process can list any global node
423 
424 .seealso: PCBDDC
425 @*/
426 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
427 {
428   PetscErrorCode ierr;
429 
430   PetscFunctionBegin;
431   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
432   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
433   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
434   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
435   PetscFunctionReturn(0);
436 }
437 /* -------------------------------------------------------------------------- */
438 
439 #undef __FUNCT__
440 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
441 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
442 {
443   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
444   PetscErrorCode ierr;
445 
446   PetscFunctionBegin;
447   /* last user setting takes precendence -> destroy any other customization */
448   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
449   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
450   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
451   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
452   pcbddc->recompute_topography = PETSC_TRUE;
453   PetscFunctionReturn(0);
454 }
455 
456 #undef __FUNCT__
457 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
458 /*@
459  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
460 
461    Collective
462 
463    Input Parameters:
464 +  pc - the preconditioning context
465 -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
466 
467    Level: intermediate
468 
469    Notes:
470 
471 .seealso: PCBDDC
472 @*/
473 PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
474 {
475   PetscErrorCode ierr;
476 
477   PetscFunctionBegin;
478   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
479   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
480   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
481   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
482   PetscFunctionReturn(0);
483 }
484 /* -------------------------------------------------------------------------- */
485 
486 #undef __FUNCT__
487 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
488 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
489 {
490   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
491 
492   PetscFunctionBegin;
493   *DirichletBoundaries = pcbddc->DirichletBoundaries;
494   PetscFunctionReturn(0);
495 }
496 
497 #undef __FUNCT__
498 #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
499 /*@
500  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
501 
502    Collective
503 
504    Input Parameters:
505 .  pc - the preconditioning context
506 
507    Output Parameters:
508 .  DirichletBoundaries - index set defining the Dirichlet boundaries
509 
510    Level: intermediate
511 
512    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
513 
514 .seealso: PCBDDC
515 @*/
516 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
517 {
518   PetscErrorCode ierr;
519 
520   PetscFunctionBegin;
521   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
522   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
523   PetscFunctionReturn(0);
524 }
525 /* -------------------------------------------------------------------------- */
526 
527 #undef __FUNCT__
528 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
529 static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
530 {
531   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
532 
533   PetscFunctionBegin;
534   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
535   PetscFunctionReturn(0);
536 }
537 
538 #undef __FUNCT__
539 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
540 /*@
541  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
542 
543    Collective
544 
545    Input Parameters:
546 .  pc - the preconditioning context
547 
548    Output Parameters:
549 .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
550 
551    Level: intermediate
552 
553    Notes:
554 
555 .seealso: PCBDDC
556 @*/
557 PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
558 {
559   PetscErrorCode ierr;
560 
561   PetscFunctionBegin;
562   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
563   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
564   PetscFunctionReturn(0);
565 }
566 /* -------------------------------------------------------------------------- */
567 
568 #undef __FUNCT__
569 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
570 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
571 {
572   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
573 
574   PetscFunctionBegin;
575   *NeumannBoundaries = pcbddc->NeumannBoundaries;
576   PetscFunctionReturn(0);
577 }
578 
579 #undef __FUNCT__
580 #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
581 /*@
582  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
583 
584    Collective
585 
586    Input Parameters:
587 .  pc - the preconditioning context
588 
589    Output Parameters:
590 .  NeumannBoundaries - index set defining the Neumann boundaries
591 
592    Level: intermediate
593 
594    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
595 
596 .seealso: PCBDDC
597 @*/
598 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
599 {
600   PetscErrorCode ierr;
601 
602   PetscFunctionBegin;
603   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
604   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
605   PetscFunctionReturn(0);
606 }
607 /* -------------------------------------------------------------------------- */
608 
609 #undef __FUNCT__
610 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
611 static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
612 {
613   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
614 
615   PetscFunctionBegin;
616   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
617   PetscFunctionReturn(0);
618 }
619 
620 #undef __FUNCT__
621 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
622 /*@
623  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
624 
625    Collective
626 
627    Input Parameters:
628 .  pc - the preconditioning context
629 
630    Output Parameters:
631 .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
632 
633    Level: intermediate
634 
635    Notes:
636 
637 .seealso: PCBDDC
638 @*/
639 PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
640 {
641   PetscErrorCode ierr;
642 
643   PetscFunctionBegin;
644   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
645   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
646   PetscFunctionReturn(0);
647 }
648 /* -------------------------------------------------------------------------- */
649 
650 #undef __FUNCT__
651 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
652 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
653 {
654   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
655   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
656   PetscErrorCode ierr;
657 
658   PetscFunctionBegin;
659   /* free old CSR */
660   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
661   /* TODO: PCBDDCGraphSetAdjacency */
662   /* get CSR into graph structure */
663   if (copymode == PETSC_COPY_VALUES) {
664     ierr = PetscMalloc1((nvtxs+1),&mat_graph->xadj);CHKERRQ(ierr);
665     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
666     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
667     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
668   } else if (copymode == PETSC_OWN_POINTER) {
669     mat_graph->xadj = (PetscInt*)xadj;
670     mat_graph->adjncy = (PetscInt*)adjncy;
671   }
672   mat_graph->nvtxs_csr = nvtxs;
673   PetscFunctionReturn(0);
674 }
675 
676 #undef __FUNCT__
677 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
678 /*@
679  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
680 
681    Not collective
682 
683    Input Parameters:
684 +  pc - the preconditioning context
685 .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
686 .  xadj, adjncy - the CSR graph
687 -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
688 
689    Level: intermediate
690 
691    Notes:
692 
693 .seealso: PCBDDC,PetscCopyMode
694 @*/
695 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
696 {
697   void (*f)(void) = 0;
698   PetscErrorCode ierr;
699 
700   PetscFunctionBegin;
701   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
702   PetscValidIntPointer(xadj,3);
703   PetscValidIntPointer(xadj,4);
704   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
705     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
706   }
707   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
708   /* free arrays if PCBDDC is not the PC type */
709   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
710   if (!f && copymode == PETSC_OWN_POINTER) {
711     ierr = PetscFree(xadj);CHKERRQ(ierr);
712     ierr = PetscFree(adjncy);CHKERRQ(ierr);
713   }
714   PetscFunctionReturn(0);
715 }
716 /* -------------------------------------------------------------------------- */
717 
718 #undef __FUNCT__
719 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
720 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
721 {
722   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
723   PetscInt i;
724   PetscErrorCode ierr;
725 
726   PetscFunctionBegin;
727   /* Destroy ISes if they were already set */
728   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
729     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
730   }
731   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
732   /* last user setting takes precendence -> destroy any other customization */
733   for (i=0;i<pcbddc->n_ISForDofs;i++) {
734     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
735   }
736   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
737   pcbddc->n_ISForDofs = 0;
738   /* allocate space then set */
739   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
740   for (i=0;i<n_is;i++) {
741     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
742     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
743   }
744   pcbddc->n_ISForDofsLocal=n_is;
745   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
746   pcbddc->recompute_topography = PETSC_TRUE;
747   PetscFunctionReturn(0);
748 }
749 
750 #undef __FUNCT__
751 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
752 /*@
753  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
754 
755    Collective
756 
757    Input Parameters:
758 +  pc - the preconditioning context
759 -  n_is - number of index sets defining the fields
760 .  ISForDofs - array of IS describing the fields in local ordering
761 
762    Level: intermediate
763 
764    Notes: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to a different field.
765 
766 .seealso: PCBDDC
767 @*/
768 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
769 {
770   PetscInt       i;
771   PetscErrorCode ierr;
772 
773   PetscFunctionBegin;
774   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
775   PetscValidLogicalCollectiveInt(pc,n_is,2);
776   for (i=0;i<n_is;i++) {
777     PetscCheckSameComm(pc,1,ISForDofs[i],3);
778     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
779   }
780   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
781   PetscFunctionReturn(0);
782 }
783 /* -------------------------------------------------------------------------- */
784 
785 #undef __FUNCT__
786 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
787 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
788 {
789   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
790   PetscInt i;
791   PetscErrorCode ierr;
792 
793   PetscFunctionBegin;
794   /* Destroy ISes if they were already set */
795   for (i=0;i<pcbddc->n_ISForDofs;i++) {
796     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
797   }
798   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
799   /* last user setting takes precendence -> destroy any other customization */
800   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
801     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
802   }
803   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
804   pcbddc->n_ISForDofsLocal = 0;
805   /* allocate space then set */
806   ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
807   for (i=0;i<n_is;i++) {
808     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
809     pcbddc->ISForDofs[i]=ISForDofs[i];
810   }
811   pcbddc->n_ISForDofs=n_is;
812   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
813   pcbddc->recompute_topography = PETSC_TRUE;
814   PetscFunctionReturn(0);
815 }
816 
817 #undef __FUNCT__
818 #define __FUNCT__ "PCBDDCSetDofsSplitting"
819 /*@
820  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
821 
822    Collective
823 
824    Input Parameters:
825 +  pc - the preconditioning context
826 -  n_is - number of index sets defining the fields
827 .  ISForDofs - array of IS describing the fields in global ordering
828 
829    Level: intermediate
830 
831    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.
832 
833 .seealso: PCBDDC
834 @*/
835 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
836 {
837   PetscInt       i;
838   PetscErrorCode ierr;
839 
840   PetscFunctionBegin;
841   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
842   PetscValidLogicalCollectiveInt(pc,n_is,2);
843   for (i=0;i<n_is;i++) {
844     PetscCheckSameComm(pc,1,ISForDofs[i],3);
845     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
846   }
847   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
848   PetscFunctionReturn(0);
849 }
850 /* -------------------------------------------------------------------------- */
851 #undef __FUNCT__
852 #define __FUNCT__ "PCPreSolve_BDDC"
853 /* -------------------------------------------------------------------------- */
854 /*
855    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
856                      guess if a transformation of basis approach has been selected.
857 
858    Input Parameter:
859 +  pc - the preconditioner contex
860 
861    Application Interface Routine: PCPreSolve()
862 
863    Notes:
864    The interface routine PCPreSolve() is not usually called directly by
865    the user, but instead is called by KSPSolve().
866 */
867 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
868 {
869   PetscErrorCode ierr;
870   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
871   PC_IS          *pcis = (PC_IS*)(pc->data);
872   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
873   Mat            temp_mat;
874   IS             dirIS;
875   Vec            used_vec;
876   PetscBool      guess_nonzero;
877 
878   PetscFunctionBegin;
879   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
880   if (ksp) {
881     PetscBool iscg;
882     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
883     if (!iscg) {
884       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
885     }
886   }
887   /* Creates parallel work vectors used in presolve */
888   if (!pcbddc->original_rhs) {
889     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
890   }
891   if (!pcbddc->temp_solution) {
892     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
893   }
894   if (x) {
895     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
896     used_vec = x;
897   } else {
898     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
899     used_vec = pcbddc->temp_solution;
900     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
901   }
902   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
903   if (ksp) {
904     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
905     if (!guess_nonzero) {
906       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
907     }
908   }
909 
910   /* store the original rhs */
911   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
912 
913   /* Take into account zeroed rows -> change rhs and store solution removed */
914   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
915   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
916   if (rhs && dirIS) {
917     PetscInt    dirsize,i,*is_indices;
918     PetscScalar *array_x,*array_diagonal;
919 
920     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
921     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
922     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
923     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
924     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
925     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
926     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
927     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
928     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
929     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
930     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
931     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
932     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
933     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
934     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
935     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
936 
937     /* remove the computed solution from the rhs */
938     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
939     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
940     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
941   }
942 
943   /* store partially computed solution and set initial guess */
944   if (x) {
945     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
946     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
947     if (pcbddc->use_exact_dirichlet_trick) {
948       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
949       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
950       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
951       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
952       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
953       if (ksp) {
954         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
955       }
956     }
957   }
958 
959   /* prepare MatMult and rhs for solver */
960   if (pcbddc->use_change_of_basis) {
961     /* swap pointers for local matrices */
962     temp_mat = matis->A;
963     matis->A = pcbddc->local_mat;
964     pcbddc->local_mat = temp_mat;
965     if (rhs) {
966       /* Get local rhs and apply transformation of basis */
967       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
968       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
969       /* from original basis to modified basis */
970       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
971       /* put back modified values into the global vec using INSERT_VALUES copy mode */
972       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
973       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
974     }
975   }
976 
977   /* remove nullspace if present */
978   if (ksp && pcbddc->NullSpace) {
979     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
980     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
981   }
982   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
983   PetscFunctionReturn(0);
984 }
985 /* -------------------------------------------------------------------------- */
986 #undef __FUNCT__
987 #define __FUNCT__ "PCPostSolve_BDDC"
988 /* -------------------------------------------------------------------------- */
989 /*
990    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
991                      approach has been selected. Also, restores rhs to its original state.
992 
993    Input Parameter:
994 +  pc - the preconditioner contex
995 
996    Application Interface Routine: PCPostSolve()
997 
998    Notes:
999    The interface routine PCPostSolve() is not usually called directly by
1000    the user, but instead is called by KSPSolve().
1001 */
1002 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1003 {
1004   PetscErrorCode ierr;
1005   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1006   PC_IS          *pcis   = (PC_IS*)(pc->data);
1007   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
1008   Mat            temp_mat;
1009 
1010   PetscFunctionBegin;
1011   if (pcbddc->use_change_of_basis) {
1012     /* swap pointers for local matrices */
1013     temp_mat = matis->A;
1014     matis->A = pcbddc->local_mat;
1015     pcbddc->local_mat = temp_mat;
1016   }
1017   if (pcbddc->use_change_of_basis && x) {
1018     /* Get Local boundary and apply transformation of basis to solution vector */
1019     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1020     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1021     /* from modified basis to original basis */
1022     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
1023     /* put back modified values into the global vec using INSERT_VALUES copy mode */
1024     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1025     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1026   }
1027   /* add solution removed in presolve */
1028   if (x) {
1029     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
1030   }
1031   /* restore rhs to its original state */
1032   if (rhs) {
1033     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1034   }
1035   PetscFunctionReturn(0);
1036 }
1037 /* -------------------------------------------------------------------------- */
1038 #undef __FUNCT__
1039 #define __FUNCT__ "PCSetUp_BDDC"
1040 /* -------------------------------------------------------------------------- */
1041 /*
1042    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1043                   by setting data structures and options.
1044 
1045    Input Parameter:
1046 +  pc - the preconditioner context
1047 
1048    Application Interface Routine: PCSetUp()
1049 
1050    Notes:
1051    The interface routine PCSetUp() is not usually called directly by
1052    the user, but instead is called by PCApply() if necessary.
1053 */
1054 PetscErrorCode PCSetUp_BDDC(PC pc)
1055 {
1056   PetscErrorCode   ierr;
1057   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
1058   MatNullSpace     nearnullspace;
1059   PetscBool        computeis,computetopography,computesolvers;
1060   PetscBool        new_nearnullspace_provided;
1061 
1062   PetscFunctionBegin;
1063   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1064   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1065      Also, BDDC directly build the Dirichlet problem */
1066 
1067   /* split work */
1068   if (pc->setupcalled) {
1069     computeis = PETSC_FALSE;
1070     if (pc->flag == SAME_NONZERO_PATTERN) {
1071       computetopography = PETSC_FALSE;
1072       computesolvers = PETSC_TRUE;
1073     } else { /* DIFFERENT_NONZERO_PATTERN */
1074       computetopography = PETSC_TRUE;
1075       computesolvers = PETSC_TRUE;
1076     }
1077   } else {
1078     computeis = PETSC_TRUE;
1079     computetopography = PETSC_TRUE;
1080     computesolvers = PETSC_TRUE;
1081   }
1082   if (pcbddc->recompute_topography) {
1083     computetopography = PETSC_TRUE;
1084   }
1085 
1086   /* Get stdout for dbg */
1087   if (pcbddc->dbg_flag) {
1088     if (!pcbddc->dbg_viewer) {
1089       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1090       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1091     }
1092     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1093   }
1094 
1095   /* Set up all the "iterative substructuring" common block without computing solvers */
1096   if (computeis) {
1097     /* HACK INTO PCIS */
1098     PC_IS* pcis = (PC_IS*)pc->data;
1099     pcis->computesolvers = PETSC_FALSE;
1100     ierr = PCISSetUp(pc);CHKERRQ(ierr);
1101     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr);
1102   }
1103 
1104   /* Analyze interface */
1105   if (computetopography) {
1106     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1107   }
1108 
1109   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1110   new_nearnullspace_provided = PETSC_FALSE;
1111   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1112   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1113     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1114       new_nearnullspace_provided = PETSC_TRUE;
1115     } else {
1116       /* determine if the two nullspaces are different (should be lightweight) */
1117       if (nearnullspace != pcbddc->onearnullspace) {
1118         new_nearnullspace_provided = PETSC_TRUE;
1119       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1120         PetscInt         i;
1121         const Vec        *nearnullvecs;
1122         PetscObjectState state;
1123         PetscInt         nnsp_size;
1124         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1125         for (i=0;i<nnsp_size;i++) {
1126           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1127           if (pcbddc->onearnullvecs_state[i] != state) {
1128             new_nearnullspace_provided = PETSC_TRUE;
1129             break;
1130           }
1131         }
1132       }
1133     }
1134   } else {
1135     if (!nearnullspace) { /* both nearnullspaces are null */
1136       new_nearnullspace_provided = PETSC_FALSE;
1137     } else { /* nearnullspace attached later */
1138       new_nearnullspace_provided = PETSC_TRUE;
1139     }
1140   }
1141 
1142   /* Setup constraints and related work vectors */
1143   /* reset primal space flags */
1144   pcbddc->new_primal_space = PETSC_FALSE;
1145   pcbddc->new_primal_space_local = PETSC_FALSE;
1146   if (computetopography || new_nearnullspace_provided) {
1147     /* It also sets the primal space flags */
1148     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1149     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1150     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
1151   }
1152 
1153   if (computesolvers || pcbddc->new_primal_space) {
1154     /* reset data */
1155     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1156     /* Create coarse and local stuffs */
1157     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1158     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
1159   }
1160 
1161   if (pcbddc->dbg_flag) {
1162     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1163   }
1164   PetscFunctionReturn(0);
1165 }
1166 
1167 /* -------------------------------------------------------------------------- */
1168 /*
1169    PCApply_BDDC - Applies the BDDC operator to a vector.
1170 
1171    Input Parameters:
1172 .  pc - the preconditioner context
1173 .  r - input vector (global)
1174 
1175    Output Parameter:
1176 .  z - output vector (global)
1177 
1178    Application Interface Routine: PCApply()
1179  */
1180 #undef __FUNCT__
1181 #define __FUNCT__ "PCApply_BDDC"
1182 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1183 {
1184   PC_IS             *pcis = (PC_IS*)(pc->data);
1185   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1186   PetscErrorCode    ierr;
1187   const PetscScalar one = 1.0;
1188   const PetscScalar m_one = -1.0;
1189   const PetscScalar zero = 0.0;
1190 
1191 /* This code is similar to that provided in nn.c for PCNN
1192    NN interface preconditioner changed to BDDC
1193    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
1194 
1195   PetscFunctionBegin;
1196   if (!pcbddc->use_exact_dirichlet_trick) {
1197     /* First Dirichlet solve */
1198     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1199     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1200     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1201     /*
1202       Assembling right hand side for BDDC operator
1203       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1204       - pcis->vec1_B the interface part of the global vector z
1205     */
1206     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1207     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1208     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1209     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1210     ierr = VecCopy(r,z);CHKERRQ(ierr);
1211     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1212     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1213     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1214   } else {
1215     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1216     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1217     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1218   }
1219 
1220   /* Apply interface preconditioner
1221      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1222   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
1223 
1224   /* Apply transpose of partition of unity operator */
1225   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1226 
1227   /* Second Dirichlet solve and assembling of output */
1228   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1229   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1230   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1231   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1232   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1233   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1234   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1235   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
1236   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1237   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1238   PetscFunctionReturn(0);
1239 }
1240 
1241 /* -------------------------------------------------------------------------- */
1242 /*
1243    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
1244 
1245    Input Parameters:
1246 .  pc - the preconditioner context
1247 .  r - input vector (global)
1248 
1249    Output Parameter:
1250 .  z - output vector (global)
1251 
1252    Application Interface Routine: PCApplyTranspose()
1253  */
1254 #undef __FUNCT__
1255 #define __FUNCT__ "PCApplyTranspose_BDDC"
1256 PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
1257 {
1258   PC_IS             *pcis = (PC_IS*)(pc->data);
1259   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1260   PetscErrorCode    ierr;
1261   const PetscScalar one = 1.0;
1262   const PetscScalar m_one = -1.0;
1263   const PetscScalar zero = 0.0;
1264 
1265   PetscFunctionBegin;
1266   if (!pcbddc->use_exact_dirichlet_trick) {
1267     /* First Dirichlet solve */
1268     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1269     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1270     ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1271     /*
1272       Assembling right hand side for BDDC operator
1273       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1274       - pcis->vec1_B the interface part of the global vector z
1275     */
1276     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1277     ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1278     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1279     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1280     ierr = VecCopy(r,z);CHKERRQ(ierr);
1281     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1282     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1283     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1284   } else {
1285     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1286     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1287     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1288   }
1289 
1290   /* Apply interface preconditioner
1291      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1292   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
1293 
1294   /* Apply transpose of partition of unity operator */
1295   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1296 
1297   /* Second Dirichlet solve and assembling of output */
1298   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1299   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1300   ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1301   if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1302   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1303   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1304   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1305   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
1306   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1307   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1308   PetscFunctionReturn(0);
1309 }
1310 /* -------------------------------------------------------------------------- */
1311 
1312 #undef __FUNCT__
1313 #define __FUNCT__ "PCDestroy_BDDC"
1314 PetscErrorCode PCDestroy_BDDC(PC pc)
1315 {
1316   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1317   PetscErrorCode ierr;
1318 
1319   PetscFunctionBegin;
1320   /* free data created by PCIS */
1321   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1322   /* free BDDC custom data  */
1323   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1324   /* destroy objects related to topography */
1325   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1326   /* free allocated graph structure */
1327   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1328   /* free data for scaling operator */
1329   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1330   /* free solvers stuff */
1331   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1332   /* free global vectors needed in presolve */
1333   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
1334   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1335   ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr);
1336   /* remove functions */
1337   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1338   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
1339   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1340   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
1341   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1342   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1343   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1344   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1345   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1346   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1347   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1348   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1349   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1350   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1351   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1352   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1353   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1354   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1355   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1356   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1357   /* Free the private data structure */
1358   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1359   PetscFunctionReturn(0);
1360 }
1361 /* -------------------------------------------------------------------------- */
1362 
1363 #undef __FUNCT__
1364 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
1365 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1366 {
1367   FETIDPMat_ctx  mat_ctx;
1368   PC_IS*         pcis;
1369   PC_BDDC*       pcbddc;
1370   PetscErrorCode ierr;
1371 
1372   PetscFunctionBegin;
1373   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1374   pcis = (PC_IS*)mat_ctx->pc->data;
1375   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1376 
1377   /* change of basis for physical rhs if needed
1378      It also changes the rhs in case of dirichlet boundaries */
1379   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
1380   /* store vectors for computation of fetidp final solution */
1381   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1382   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1383   /* scale rhs since it should be unassembled */
1384   /* TODO use counter scaling? (also below) */
1385   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1386   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1387   /* Apply partition of unity */
1388   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1389   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1390   if (!pcbddc->switch_static) {
1391     /* compute partially subassembled Schur complement right-hand side */
1392     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1393     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
1394     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1395     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
1396     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1397     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1398     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1399     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1400     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1401     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1402   }
1403   /* BDDC rhs */
1404   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
1405   if (pcbddc->switch_static) {
1406     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1407   }
1408   /* apply BDDC */
1409   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
1410   /* Application of B_delta and assembling of rhs for fetidp fluxes */
1411   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
1412   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
1413   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1414   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1415   /* restore original rhs */
1416   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
1417   PetscFunctionReturn(0);
1418 }
1419 
1420 #undef __FUNCT__
1421 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
1422 /*@
1423  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
1424 
1425    Collective
1426 
1427    Input Parameters:
1428 +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
1429 .  standard_rhs - the right-hand side for your linear system
1430 
1431    Output Parameters:
1432 -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
1433 
1434    Level: developer
1435 
1436    Notes:
1437 
1438 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1439 @*/
1440 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1441 {
1442   FETIDPMat_ctx  mat_ctx;
1443   PetscErrorCode ierr;
1444 
1445   PetscFunctionBegin;
1446   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1447   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
1448   PetscFunctionReturn(0);
1449 }
1450 /* -------------------------------------------------------------------------- */
1451 
1452 #undef __FUNCT__
1453 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
1454 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1455 {
1456   FETIDPMat_ctx  mat_ctx;
1457   PC_IS*         pcis;
1458   PC_BDDC*       pcbddc;
1459   PetscErrorCode ierr;
1460 
1461   PetscFunctionBegin;
1462   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1463   pcis = (PC_IS*)mat_ctx->pc->data;
1464   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1465 
1466   /* apply B_delta^T */
1467   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1468   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1469   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
1470   /* compute rhs for BDDC application */
1471   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1472   if (pcbddc->switch_static) {
1473     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1474   }
1475   /* apply BDDC */
1476   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
1477   /* put values into standard global vector */
1478   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1479   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1480   if (!pcbddc->switch_static) {
1481     /* compute values into the interior if solved for the partially subassembled Schur complement */
1482     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
1483     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
1484     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1485   }
1486   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1487   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1488   /* final change of basis if needed
1489      Is also sums the dirichlet part removed during RHS assembling */
1490   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
1491   PetscFunctionReturn(0);
1492 }
1493 
1494 #undef __FUNCT__
1495 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
1496 /*@
1497  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
1498 
1499    Collective
1500 
1501    Input Parameters:
1502 +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
1503 .  fetidp_flux_sol - the solution of the FETIDP linear system
1504 
1505    Output Parameters:
1506 -  standard_sol      - the solution defined on the physical domain
1507 
1508    Level: developer
1509 
1510    Notes:
1511 
1512 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1513 @*/
1514 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1515 {
1516   FETIDPMat_ctx  mat_ctx;
1517   PetscErrorCode ierr;
1518 
1519   PetscFunctionBegin;
1520   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1521   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
1522   PetscFunctionReturn(0);
1523 }
1524 /* -------------------------------------------------------------------------- */
1525 
1526 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1527 extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1528 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1529 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1530 extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1531 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1532 
1533 #undef __FUNCT__
1534 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
1535 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1536 {
1537 
1538   FETIDPMat_ctx  fetidpmat_ctx;
1539   Mat            newmat;
1540   FETIDPPC_ctx   fetidppc_ctx;
1541   PC             newpc;
1542   MPI_Comm       comm;
1543   PetscErrorCode ierr;
1544 
1545   PetscFunctionBegin;
1546   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1547   /* FETIDP linear matrix */
1548   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
1549   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
1550   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
1551   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1552   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
1553   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
1554   ierr = MatSetUp(newmat);CHKERRQ(ierr);
1555   /* FETIDP preconditioner */
1556   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
1557   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
1558   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
1559   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
1560   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
1561   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1562   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
1563   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
1564   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
1565   ierr = PCSetUp(newpc);CHKERRQ(ierr);
1566   /* return pointers for objects created */
1567   *fetidp_mat=newmat;
1568   *fetidp_pc=newpc;
1569   PetscFunctionReturn(0);
1570 }
1571 
1572 #undef __FUNCT__
1573 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
1574 /*@
1575  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
1576 
1577    Collective
1578 
1579    Input Parameters:
1580 +  pc - the BDDC preconditioning context already setup
1581 
1582    Output Parameters:
1583 .  fetidp_mat - shell FETIDP matrix object
1584 .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
1585 
1586    Options Database Keys:
1587 -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
1588 
1589    Level: developer
1590 
1591    Notes:
1592      Currently the only operation provided for FETIDP matrix is MatMult
1593 
1594 .seealso: PCBDDC
1595 @*/
1596 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1597 {
1598   PetscErrorCode ierr;
1599 
1600   PetscFunctionBegin;
1601   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1602   if (pc->setupcalled) {
1603     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1604   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1605   PetscFunctionReturn(0);
1606 }
1607 /* -------------------------------------------------------------------------- */
1608 /*MC
1609    PCBDDC - Balancing Domain Decomposition by Constraints.
1610 
1611    An implementation of the BDDC preconditioner based on
1612 
1613 .vb
1614    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
1615    [2] A. Klawonn and O. B. Widlund. "Dual-Primal FETI Methods for Linear Elasticity", http://cs.nyu.edu/csweb/Research/TechReports/TR2004-855/TR2004-855.pdf
1616    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
1617 .ve
1618 
1619    The matrix to be preconditioned (Pmat) must be of type MATIS.
1620 
1621    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
1622 
1623    It also works with unsymmetric and indefinite problems.
1624 
1625    Unlike 'conventional' interface preconditioners, PCBDDC iterates over all degrees of freedom, not just those on the interface. This allows the use of approximate solvers on the subdomains.
1626 
1627    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
1628 
1629    Boundary nodes are split in vertices, edges and faces using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph
1630 
1631    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace.
1632 
1633    Change of basis is performed similarly to [2] when requested. When more the one constraint is present on a single connected component (i.e. an edge or a face), a robust method based on local QR factorizations is used.
1634 
1635    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
1636 
1637    Options Database Keys:
1638 
1639 .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
1640 .    -pc_bddc_use_edges <1> - use or not edges in primal space
1641 .    -pc_bddc_use_faces <0> - use or not faces in primal space
1642 .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
1643 .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
1644 .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
1645 .    -pc_bddc_levels <0> - maximum number of levels for multilevel
1646 .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
1647 -    -pc_bddc_check_level <0> - set verbosity level of debugging output
1648 
1649    Options for Dirichlet, Neumann or coarse solver can be set with
1650 .vb
1651       -pc_bddc_dirichlet_
1652       -pc_bddc_neumann_
1653       -pc_bddc_coarse_
1654 .ve
1655    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
1656 
1657    When using a multilevel approach, solvers' options at the N-th level can be specified as
1658 .vb
1659       -pc_bddc_dirichlet_lN_
1660       -pc_bddc_neumann_lN_
1661       -pc_bddc_coarse_lN_
1662 .ve
1663    Note that level number ranges from the finest 0 to the coarsest N.
1664 
1665    Level: intermediate
1666 
1667    Developer notes:
1668 
1669    New deluxe scaling operator will be available soon.
1670 
1671    Contributed by Stefano Zampini
1672 
1673 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1674 M*/
1675 
1676 #undef __FUNCT__
1677 #define __FUNCT__ "PCCreate_BDDC"
1678 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1679 {
1680   PetscErrorCode      ierr;
1681   PC_BDDC             *pcbddc;
1682 
1683   PetscFunctionBegin;
1684   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1685   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
1686   pc->data  = (void*)pcbddc;
1687 
1688   /* create PCIS data structure */
1689   ierr = PCISCreate(pc);CHKERRQ(ierr);
1690 
1691   /* BDDC customization */
1692   pcbddc->use_local_adj       = PETSC_TRUE;
1693   pcbddc->use_vertices        = PETSC_TRUE;
1694   pcbddc->use_edges           = PETSC_TRUE;
1695   pcbddc->use_faces           = PETSC_FALSE;
1696   pcbddc->use_change_of_basis = PETSC_FALSE;
1697   pcbddc->use_change_on_faces = PETSC_FALSE;
1698   pcbddc->switch_static       = PETSC_FALSE;
1699   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
1700   pcbddc->dbg_flag            = 0;
1701 
1702   pcbddc->BtoNmap                    = 0;
1703   pcbddc->local_primal_size          = 0;
1704   pcbddc->n_vertices                 = 0;
1705   pcbddc->n_actual_vertices          = 0;
1706   pcbddc->n_constraints              = 0;
1707   pcbddc->primal_indices_local_idxs  = 0;
1708   pcbddc->recompute_topography       = PETSC_FALSE;
1709   pcbddc->coarse_size                = -1;
1710   pcbddc->new_primal_space           = PETSC_FALSE;
1711   pcbddc->new_primal_space_local     = PETSC_FALSE;
1712   pcbddc->global_primal_indices      = 0;
1713   pcbddc->onearnullspace             = 0;
1714   pcbddc->onearnullvecs_state        = 0;
1715   pcbddc->user_primal_vertices       = 0;
1716   pcbddc->NullSpace                  = 0;
1717   pcbddc->temp_solution              = 0;
1718   pcbddc->original_rhs               = 0;
1719   pcbddc->local_mat                  = 0;
1720   pcbddc->ChangeOfBasisMatrix        = 0;
1721   pcbddc->coarse_vec                 = 0;
1722   pcbddc->coarse_rhs                 = 0;
1723   pcbddc->coarse_ksp                 = 0;
1724   pcbddc->coarse_phi_B               = 0;
1725   pcbddc->coarse_phi_D               = 0;
1726   pcbddc->coarse_psi_B               = 0;
1727   pcbddc->coarse_psi_D               = 0;
1728   pcbddc->vec1_P                     = 0;
1729   pcbddc->vec1_R                     = 0;
1730   pcbddc->vec2_R                     = 0;
1731   pcbddc->local_auxmat1              = 0;
1732   pcbddc->local_auxmat2              = 0;
1733   pcbddc->R_to_B                     = 0;
1734   pcbddc->R_to_D                     = 0;
1735   pcbddc->ksp_D                      = 0;
1736   pcbddc->ksp_R                      = 0;
1737   pcbddc->NeumannBoundaries          = 0;
1738   pcbddc->NeumannBoundariesLocal     = 0;
1739   pcbddc->DirichletBoundaries        = 0;
1740   pcbddc->DirichletBoundariesLocal   = 0;
1741   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
1742   pcbddc->n_ISForDofs                = 0;
1743   pcbddc->n_ISForDofsLocal           = 0;
1744   pcbddc->ISForDofs                  = 0;
1745   pcbddc->ISForDofsLocal             = 0;
1746   pcbddc->ConstraintMatrix           = 0;
1747   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
1748   pcbddc->coarse_loc_to_glob         = 0;
1749   pcbddc->coarsening_ratio           = 8;
1750   pcbddc->current_level              = 0;
1751   pcbddc->max_levels                 = 0;
1752   pcbddc->use_coarse_estimates       = PETSC_FALSE;
1753   pcbddc->coarse_subassembling       = 0;
1754   pcbddc->coarse_subassembling_init  = 0;
1755 
1756   /* create local graph structure */
1757   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1758 
1759   /* scaling */
1760   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1761   pcbddc->deluxe_ctx                 = 0;
1762   pcbddc->work_scaling               = 0;
1763 
1764   /* function pointers */
1765   pc->ops->apply               = PCApply_BDDC;
1766   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
1767   pc->ops->setup               = PCSetUp_BDDC;
1768   pc->ops->destroy             = PCDestroy_BDDC;
1769   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1770   pc->ops->view                = 0;
1771   pc->ops->applyrichardson     = 0;
1772   pc->ops->applysymmetricleft  = 0;
1773   pc->ops->applysymmetricright = 0;
1774   pc->ops->presolve            = PCPreSolve_BDDC;
1775   pc->ops->postsolve           = PCPostSolve_BDDC;
1776 
1777   /* composing function */
1778   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1779   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
1780   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1781   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
1782   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1783   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1784   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1785   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1786   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1787   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1788   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1789   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1790   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1791   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1792   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1793   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
1794   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1795   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1796   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1797   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1798   PetscFunctionReturn(0);
1799 }
1800 
1801