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