xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 785d12438418abfe1dc590cf0e72c5e5ca69a49b)
1 /* TODOLIST
2 
3    ConstraintsSetup
4    - assure same constraints between neighbours by sorting vals by global index before SVD!
5    - tolerances for constraints as an option (take care of single precision!)
6    - Allow different constraints customizations among different linear solves (requires also reset/destroy of ksp_R and coarse_ksp)
7    - MAT_IGNORE_ZERO_ENTRIES for Constraints Matrix
8 
9    Solvers
10    - Try to reduce the work when reusing the solvers
11    - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers)
12    - reuse already allocated coarse matrix if possible
13    - Propagate ksp prefixes for solvers to mat objects?
14    - Propagate nearnullspace info among levels
15 
16    User interface
17    - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access)
18    - Provide PCApplyTranpose_BDDC
19    - DofSplitting and DM attached to pc?
20 
21    Debugging output
22    - Better management of verbosity levels of debugging output
23    - Crashes on some architecture -> call SynchronizedAllow before every SynchronizedPrintf
24 
25    Build
26    - make runexe59
27 
28    Extra
29    - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
30    - Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver?
31    - add support for computing h,H and related using coordinates?
32    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
33    - Better management in PCIS code
34    - BDDC with MG framework?
35 
36    FETIDP
37    - Move FETIDP code to its own classes
38 
39    MATIS related operations contained in BDDC code
40    - Provide general case for subassembling
41    - Preallocation routines in MatConvert_IS_AIJ
42 
43 */
44 
45 /* ----------------------------------------------------------------------------------------------------------------------------------------------
46    Implementation of BDDC preconditioner based on:
47    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
48    ---------------------------------------------------------------------------------------------------------------------------------------------- */
49 
50 #include "bddc.h" /*I "petscpc.h" I*/  /* includes for fortran wrappers */
51 #include "bddcprivate.h"
52 #include <petscblaslapack.h>
53 
54 /* -------------------------------------------------------------------------- */
55 #undef __FUNCT__
56 #define __FUNCT__ "PCSetFromOptions_BDDC"
57 PetscErrorCode PCSetFromOptions_BDDC(PC pc)
58 {
59   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
60   PetscErrorCode ierr;
61 
62   PetscFunctionBegin;
63   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
64   /* Verbose debugging */
65   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
66   /* Primal space cumstomization */
67   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
68   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
69   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
70   /* Change of basis */
71   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);
72   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);
73   if (!pcbddc->use_change_of_basis) {
74     pcbddc->use_change_on_faces = PETSC_FALSE;
75   }
76   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
77   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);
78   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
79   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
80   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
81   ierr = PetscOptionsTail();CHKERRQ(ierr);
82   PetscFunctionReturn(0);
83 }
84 /* -------------------------------------------------------------------------- */
85 #undef __FUNCT__
86 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
87 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
88 {
89   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
90   PetscErrorCode ierr;
91 
92   PetscFunctionBegin;
93   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
94   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
95   pcbddc->user_primal_vertices = PrimalVertices;
96   PetscFunctionReturn(0);
97 }
98 #undef __FUNCT__
99 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
100 /*@
101  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
102 
103    Not collective
104 
105    Input Parameters:
106 +  pc - the preconditioning context
107 -  PrimalVertices - index set of primal vertices in local numbering
108 
109    Level: intermediate
110 
111    Notes:
112 
113 .seealso: PCBDDC
114 @*/
115 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
116 {
117   PetscErrorCode ierr;
118 
119   PetscFunctionBegin;
120   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
121   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
122   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
123   PetscFunctionReturn(0);
124 }
125 /* -------------------------------------------------------------------------- */
126 #undef __FUNCT__
127 #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
128 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
129 {
130   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
131 
132   PetscFunctionBegin;
133   pcbddc->coarsening_ratio = k;
134   PetscFunctionReturn(0);
135 }
136 
137 #undef __FUNCT__
138 #define __FUNCT__ "PCBDDCSetCoarseningRatio"
139 /*@
140  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
141 
142    Logically collective on PC
143 
144    Input Parameters:
145 +  pc - the preconditioning context
146 -  k - coarsening ratio (H/h at the coarser level)
147 
148    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
149 
150    Level: intermediate
151 
152    Notes:
153 
154 .seealso: PCBDDC
155 @*/
156 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
157 {
158   PetscErrorCode ierr;
159 
160   PetscFunctionBegin;
161   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
162   PetscValidLogicalCollectiveInt(pc,k,2);
163   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
164   PetscFunctionReturn(0);
165 }
166 
167 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
168 #undef __FUNCT__
169 #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
170 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
171 {
172   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
173 
174   PetscFunctionBegin;
175   pcbddc->use_exact_dirichlet_trick = flg;
176   PetscFunctionReturn(0);
177 }
178 
179 #undef __FUNCT__
180 #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
181 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
182 {
183   PetscErrorCode ierr;
184 
185   PetscFunctionBegin;
186   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
187   PetscValidLogicalCollectiveBool(pc,flg,2);
188   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
189   PetscFunctionReturn(0);
190 }
191 
192 #undef __FUNCT__
193 #define __FUNCT__ "PCBDDCSetLevel_BDDC"
194 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
195 {
196   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
197 
198   PetscFunctionBegin;
199   pcbddc->current_level = level;
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "PCBDDCSetLevel"
205 PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
206 {
207   PetscErrorCode ierr;
208 
209   PetscFunctionBegin;
210   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
211   PetscValidLogicalCollectiveInt(pc,level,2);
212   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
213   PetscFunctionReturn(0);
214 }
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "PCBDDCSetLevels_BDDC"
218 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
219 {
220   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
221 
222   PetscFunctionBegin;
223   pcbddc->max_levels = levels;
224   PetscFunctionReturn(0);
225 }
226 
227 #undef __FUNCT__
228 #define __FUNCT__ "PCBDDCSetLevels"
229 /*@
230  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
231 
232    Logically collective on PC
233 
234    Input Parameters:
235 +  pc - the preconditioning context
236 -  levels - the maximum number of levels (max 9)
237 
238    Default value is 0, i.e. traditional one-level BDDC
239 
240    Level: intermediate
241 
242    Notes:
243 
244 .seealso: PCBDDC
245 @*/
246 PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
247 {
248   PetscErrorCode ierr;
249 
250   PetscFunctionBegin;
251   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
252   PetscValidLogicalCollectiveInt(pc,levels,2);
253   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
254   PetscFunctionReturn(0);
255 }
256 /* -------------------------------------------------------------------------- */
257 
258 #undef __FUNCT__
259 #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
260 static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
261 {
262   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
263   PetscErrorCode ierr;
264 
265   PetscFunctionBegin;
266   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
267   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
268   pcbddc->NullSpace = NullSpace;
269   PetscFunctionReturn(0);
270 }
271 
272 #undef __FUNCT__
273 #define __FUNCT__ "PCBDDCSetNullSpace"
274 /*@
275  PCBDDCSetNullSpace - Set nullspace for BDDC operator
276 
277    Logically collective on PC and MatNullSpace
278 
279    Input Parameters:
280 +  pc - the preconditioning context
281 -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
282 
283    Level: intermediate
284 
285    Notes:
286 
287 .seealso: PCBDDC
288 @*/
289 PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
290 {
291   PetscErrorCode ierr;
292 
293   PetscFunctionBegin;
294   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
295   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
296   PetscCheckSameComm(pc,1,NullSpace,2);
297   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
298   PetscFunctionReturn(0);
299 }
300 /* -------------------------------------------------------------------------- */
301 
302 #undef __FUNCT__
303 #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
304 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
305 {
306   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
307   PetscErrorCode ierr;
308 
309   PetscFunctionBegin;
310   /* last user setting takes precendence -> destroy any other customization */
311   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
312   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
313   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
314   pcbddc->DirichletBoundaries = DirichletBoundaries;
315   PetscFunctionReturn(0);
316 }
317 
318 #undef __FUNCT__
319 #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
320 /*@
321  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
322 
323    Collective
324 
325    Input Parameters:
326 +  pc - the preconditioning context
327 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
328 
329    Level: intermediate
330 
331    Notes: Any process can list any global node
332 
333 .seealso: PCBDDC
334 @*/
335 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
336 {
337   PetscErrorCode ierr;
338 
339   PetscFunctionBegin;
340   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
341   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
342   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
343   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
344   PetscFunctionReturn(0);
345 }
346 /* -------------------------------------------------------------------------- */
347 
348 #undef __FUNCT__
349 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
350 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
351 {
352   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
353   PetscErrorCode ierr;
354 
355   PetscFunctionBegin;
356   /* last user setting takes precendence -> destroy any other customization */
357   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
358   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
359   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
360   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
361   PetscFunctionReturn(0);
362 }
363 
364 #undef __FUNCT__
365 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
366 /*@
367  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
368 
369    Collective
370 
371    Input Parameters:
372 +  pc - the preconditioning context
373 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
374 
375    Level: intermediate
376 
377    Notes:
378 
379 .seealso: PCBDDC
380 @*/
381 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
382 {
383   PetscErrorCode ierr;
384 
385   PetscFunctionBegin;
386   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
387   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
388   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
389   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
390   PetscFunctionReturn(0);
391 }
392 /* -------------------------------------------------------------------------- */
393 
394 #undef __FUNCT__
395 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
396 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
397 {
398   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
399   PetscErrorCode ierr;
400 
401   PetscFunctionBegin;
402   /* last user setting takes precendence -> destroy any other customization */
403   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
404   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
405   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
406   pcbddc->NeumannBoundaries = NeumannBoundaries;
407   PetscFunctionReturn(0);
408 }
409 
410 #undef __FUNCT__
411 #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
412 /*@
413  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
414 
415    Collective
416 
417    Input Parameters:
418 +  pc - the preconditioning context
419 -  NeumannBoundaries - parallel IS defining the Neumann boundaries
420 
421    Level: intermediate
422 
423    Notes: Any process can list any global node
424 
425 .seealso: PCBDDC
426 @*/
427 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
428 {
429   PetscErrorCode ierr;
430 
431   PetscFunctionBegin;
432   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
433   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
434   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
435   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
436   PetscFunctionReturn(0);
437 }
438 /* -------------------------------------------------------------------------- */
439 
440 #undef __FUNCT__
441 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
442 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
443 {
444   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
445   PetscErrorCode ierr;
446 
447   PetscFunctionBegin;
448   /* last user setting takes precendence -> destroy any other customization */
449   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
450   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
451   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
452   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
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 = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
665     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&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__ "PCBDDCSetDofsSplitting_BDDC"
720 static PetscErrorCode PCBDDCSetDofsSplitting_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_ISForDofs;i++) {
729     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
730   }
731   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
732   /* allocate space then set */
733   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
734   for (i=0;i<n_is;i++) {
735     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
736     pcbddc->ISForDofs[i]=ISForDofs[i];
737   }
738   pcbddc->n_ISForDofs=n_is;
739   pcbddc->user_provided_isfordofs = PETSC_TRUE;
740   PetscFunctionReturn(0);
741 }
742 
743 #undef __FUNCT__
744 #define __FUNCT__ "PCBDDCSetDofsSplitting"
745 /*@
746  PCBDDCSetDofsSplitting - Set index sets defining fields of the local Neumann matrix
747 
748    Not collective
749 
750    Input Parameters:
751 +  pc - the preconditioning context
752 -  n_is - number of index sets defining the fields
753 .  ISForDofs - array of IS describing the fields
754 
755    Level: intermediate
756 
757    Notes:
758 
759 .seealso: PCBDDC
760 @*/
761 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
762 {
763   PetscInt       i;
764   PetscErrorCode ierr;
765 
766   PetscFunctionBegin;
767   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
768   for (i=0;i<n_is;i++) {
769     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2);
770   }
771   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
772   PetscFunctionReturn(0);
773 }
774 /* -------------------------------------------------------------------------- */
775 #undef __FUNCT__
776 #define __FUNCT__ "PCPreSolve_BDDC"
777 /* -------------------------------------------------------------------------- */
778 /*
779    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
780                      guess if a transformation of basis approach has been selected.
781 
782    Input Parameter:
783 +  pc - the preconditioner contex
784 
785    Application Interface Routine: PCPreSolve()
786 
787    Notes:
788    The interface routine PCPreSolve() is not usually called directly by
789    the user, but instead is called by KSPSolve().
790 */
791 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
792 {
793   PetscErrorCode ierr;
794   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
795   PC_IS          *pcis = (PC_IS*)(pc->data);
796   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
797   Mat            temp_mat;
798   IS             dirIS;
799   Vec            used_vec;
800   PetscBool      guess_nonzero;
801 
802   PetscFunctionBegin;
803   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
804   if (ksp) {
805     PetscBool iscg;
806     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
807     if (!iscg) {
808       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
809     }
810   }
811   /* Creates parallel work vectors used in presolve */
812   if (!pcbddc->original_rhs) {
813     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
814   }
815   if (!pcbddc->temp_solution) {
816     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
817   }
818   if (x) {
819     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
820     used_vec = x;
821   } else {
822     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
823     used_vec = pcbddc->temp_solution;
824     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
825   }
826   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
827   if (ksp) {
828     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
829     if (!guess_nonzero) {
830       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
831     }
832   }
833 
834   /* store the original rhs */
835   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
836 
837   /* Take into account zeroed rows -> change rhs and store solution removed */
838   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
839   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
840   if (rhs && dirIS) {
841     PetscInt    dirsize,i,*is_indices;
842     PetscScalar *array_x,*array_diagonal;
843 
844     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
845     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
846     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
847     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
848     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
849     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
850     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
851     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
852     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
853     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
854     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
855     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
856     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
857     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
858     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
859     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
860 
861     /* remove the computed solution from the rhs */
862     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
863     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
864     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
865   }
866 
867   /* store partially computed solution and set initial guess */
868   if (x) {
869     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
870     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
871     if (pcbddc->use_exact_dirichlet_trick) {
872       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
873       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
874       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
875       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
876       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
877       if (ksp) {
878         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
879       }
880     }
881   }
882 
883   /* prepare MatMult and rhs for solver */
884   if (pcbddc->use_change_of_basis) {
885     /* swap pointers for local matrices */
886     temp_mat = matis->A;
887     matis->A = pcbddc->local_mat;
888     pcbddc->local_mat = temp_mat;
889     if (rhs) {
890       /* Get local rhs and apply transformation of basis */
891       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
892       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
893       /* from original basis to modified basis */
894       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
895       /* put back modified values into the global vec using INSERT_VALUES copy mode */
896       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
897       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
898     }
899   }
900 
901   /* remove nullspace if present */
902   if (ksp && pcbddc->NullSpace) {
903     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
904     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
905   }
906   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
907   PetscFunctionReturn(0);
908 }
909 /* -------------------------------------------------------------------------- */
910 #undef __FUNCT__
911 #define __FUNCT__ "PCPostSolve_BDDC"
912 /* -------------------------------------------------------------------------- */
913 /*
914    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
915                      approach has been selected. Also, restores rhs to its original state.
916 
917    Input Parameter:
918 +  pc - the preconditioner contex
919 
920    Application Interface Routine: PCPostSolve()
921 
922    Notes:
923    The interface routine PCPostSolve() is not usually called directly by
924    the user, but instead is called by KSPSolve().
925 */
926 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
927 {
928   PetscErrorCode ierr;
929   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
930   PC_IS          *pcis   = (PC_IS*)(pc->data);
931   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
932   Mat            temp_mat;
933 
934   PetscFunctionBegin;
935   if (pcbddc->use_change_of_basis) {
936     /* swap pointers for local matrices */
937     temp_mat = matis->A;
938     matis->A = pcbddc->local_mat;
939     pcbddc->local_mat = temp_mat;
940   }
941   if (pcbddc->use_change_of_basis && x) {
942     /* Get Local boundary and apply transformation of basis to solution vector */
943     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
944     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
945     /* from modified basis to original basis */
946     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
947     /* put back modified values into the global vec using INSERT_VALUES copy mode */
948     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
949     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
950   }
951   /* add solution removed in presolve */
952   if (x) {
953     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
954   }
955   /* restore rhs to its original state */
956   if (rhs) {
957     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
958   }
959   PetscFunctionReturn(0);
960 }
961 /* -------------------------------------------------------------------------- */
962 #undef __FUNCT__
963 #define __FUNCT__ "PCSetUp_BDDC"
964 /* -------------------------------------------------------------------------- */
965 /*
966    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
967                   by setting data structures and options.
968 
969    Input Parameter:
970 +  pc - the preconditioner context
971 
972    Application Interface Routine: PCSetUp()
973 
974    Notes:
975    The interface routine PCSetUp() is not usually called directly by
976    the user, but instead is called by PCApply() if necessary.
977 */
978 PetscErrorCode PCSetUp_BDDC(PC pc)
979 {
980   PetscErrorCode ierr;
981   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
982   MatStructure   flag;
983   PetscBool      computeis,computetopography,computesolvers;
984 
985   PetscFunctionBegin;
986   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */
987   /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */
988   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
989      Also, BDDC directly build the Dirichlet problem */
990   /* Get stdout for dbg */
991   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
992     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
993     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
994     if (pcbddc->current_level) {
995       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
996     }
997   }
998   /* first attempt to split work */
999   if (pc->setupcalled) {
1000     computeis = PETSC_FALSE;
1001     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
1002     if (flag == SAME_PRECONDITIONER) {
1003       computetopography = PETSC_FALSE;
1004       computesolvers = PETSC_FALSE;
1005     } else if (flag == SAME_NONZERO_PATTERN) {
1006       computetopography = PETSC_FALSE;
1007       computesolvers = PETSC_TRUE;
1008     } else { /* DIFFERENT_NONZERO_PATTERN */
1009       computetopography = PETSC_TRUE;
1010       computesolvers = PETSC_TRUE;
1011     }
1012   } else {
1013     computeis = PETSC_TRUE;
1014     computetopography = PETSC_TRUE;
1015     computesolvers = PETSC_TRUE;
1016   }
1017   /* Set up all the "iterative substructuring" common block without computing solvers */
1018   if (computeis) {
1019     /* HACK INTO PCIS */
1020     PC_IS* pcis = (PC_IS*)pc->data;
1021     pcis->computesolvers = PETSC_FALSE;
1022     ierr = PCISSetUp(pc);CHKERRQ(ierr);
1023   }
1024   /* Analyze interface and set up local constraint and change of basis matrices */
1025   if (computetopography) {
1026     /* reset data */
1027     ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1028     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1029     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1030   }
1031   if (computesolvers) {
1032     /* reset data */
1033     ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1034     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1035     /* Create coarse and local stuffs */
1036     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1037     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
1038   }
1039   if (pcbddc->dbg_flag && pcbddc->current_level) {
1040     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
1041   }
1042   PetscFunctionReturn(0);
1043 }
1044 
1045 /* -------------------------------------------------------------------------- */
1046 /*
1047    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
1048 
1049    Input Parameters:
1050 .  pc - the preconditioner context
1051 .  r - input vector (global)
1052 
1053    Output Parameter:
1054 .  z - output vector (global)
1055 
1056    Application Interface Routine: PCApply()
1057  */
1058 #undef __FUNCT__
1059 #define __FUNCT__ "PCApply_BDDC"
1060 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1061 {
1062   PC_IS             *pcis = (PC_IS*)(pc->data);
1063   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1064   PetscErrorCode    ierr;
1065   const PetscScalar one = 1.0;
1066   const PetscScalar m_one = -1.0;
1067   const PetscScalar zero = 0.0;
1068 
1069 /* This code is similar to that provided in nn.c for PCNN
1070    NN interface preconditioner changed to BDDC
1071    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
1072 
1073   PetscFunctionBegin;
1074   if (!pcbddc->use_exact_dirichlet_trick) {
1075     /* First Dirichlet solve */
1076     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1077     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1078     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1079     /*
1080       Assembling right hand side for BDDC operator
1081       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1082       - pcis->vec1_B the interface part of the global vector z
1083     */
1084     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1085     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1086     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1087     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1088     ierr = VecCopy(r,z);CHKERRQ(ierr);
1089     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1090     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1091     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1092   } else {
1093     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1094     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1095     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1096   }
1097 
1098   /* Apply interface preconditioner
1099      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1100   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
1101 
1102   /* Apply transpose of partition of unity operator */
1103   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1104 
1105   /* Second Dirichlet solve and assembling of output */
1106   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1107   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1108   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1109   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1110   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr);
1111   ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr);
1112   if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1113   ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr);
1114   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1115   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1116   PetscFunctionReturn(0);
1117 }
1118 /* -------------------------------------------------------------------------- */
1119 
1120 #undef __FUNCT__
1121 #define __FUNCT__ "PCDestroy_BDDC"
1122 PetscErrorCode PCDestroy_BDDC(PC pc)
1123 {
1124   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1125   PetscErrorCode ierr;
1126 
1127   PetscFunctionBegin;
1128   /* free data created by PCIS */
1129   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1130   /* free BDDC custom data  */
1131   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1132   /* destroy objects related to topography */
1133   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1134   /* free allocated graph structure */
1135   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1136   /* free data for scaling operator */
1137   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1138   /* free solvers stuff */
1139   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1140   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
1141   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
1142   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
1143   ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1144   /* free global vectors needed in presolve */
1145   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
1146   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1147   /* remove functions */
1148   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1149   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
1150   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1151   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
1152   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1153   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1154   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1155   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1156   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1157   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1158   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1159   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1160   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1161   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1162   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1163   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1164   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1165   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1166   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1167   /* Free the private data structure */
1168   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1169   PetscFunctionReturn(0);
1170 }
1171 /* -------------------------------------------------------------------------- */
1172 
1173 #undef __FUNCT__
1174 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
1175 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1176 {
1177   FETIDPMat_ctx  mat_ctx;
1178   PC_IS*         pcis;
1179   PC_BDDC*       pcbddc;
1180   PetscErrorCode ierr;
1181 
1182   PetscFunctionBegin;
1183   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1184   pcis = (PC_IS*)mat_ctx->pc->data;
1185   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1186 
1187   /* change of basis for physical rhs if needed
1188      It also changes the rhs in case of dirichlet boundaries */
1189   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
1190   /* store vectors for computation of fetidp final solution */
1191   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1192   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1193   /* scale rhs since it should be unassembled */
1194   /* TODO use counter scaling? (also below) */
1195   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1196   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1197   /* Apply partition of unity */
1198   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1199   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1200   if (!pcbddc->switch_static) {
1201     /* compute partially subassembled Schur complement right-hand side */
1202     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1203     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
1204     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1205     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
1206     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1207     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1208     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1209     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1210     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1211     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1212   }
1213   /* BDDC rhs */
1214   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
1215   if (pcbddc->switch_static) {
1216     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1217   }
1218   /* apply BDDC */
1219   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
1220   /* Application of B_delta and assembling of rhs for fetidp fluxes */
1221   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
1222   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
1223   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1224   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1225   /* restore original rhs */
1226   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
1227   PetscFunctionReturn(0);
1228 }
1229 
1230 #undef __FUNCT__
1231 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
1232 /*@
1233  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
1234 
1235    Collective
1236 
1237    Input Parameters:
1238 +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
1239 .  standard_rhs - the right-hand side for your linear system
1240 
1241    Output Parameters:
1242 -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
1243 
1244    Level: developer
1245 
1246    Notes:
1247 
1248 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1249 @*/
1250 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1251 {
1252   FETIDPMat_ctx  mat_ctx;
1253   PetscErrorCode ierr;
1254 
1255   PetscFunctionBegin;
1256   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1257   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
1258   PetscFunctionReturn(0);
1259 }
1260 /* -------------------------------------------------------------------------- */
1261 
1262 #undef __FUNCT__
1263 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
1264 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1265 {
1266   FETIDPMat_ctx  mat_ctx;
1267   PC_IS*         pcis;
1268   PC_BDDC*       pcbddc;
1269   PetscErrorCode ierr;
1270 
1271   PetscFunctionBegin;
1272   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1273   pcis = (PC_IS*)mat_ctx->pc->data;
1274   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1275 
1276   /* apply B_delta^T */
1277   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1278   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1279   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
1280   /* compute rhs for BDDC application */
1281   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1282   if (pcbddc->switch_static) {
1283     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1284   }
1285   /* apply BDDC */
1286   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
1287   /* put values into standard global vector */
1288   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1289   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1290   if (!pcbddc->switch_static) {
1291     /* compute values into the interior if solved for the partially subassembled Schur complement */
1292     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
1293     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
1294     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1295   }
1296   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1297   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1298   /* final change of basis if needed
1299      Is also sums the dirichlet part removed during RHS assembling */
1300   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
1301   PetscFunctionReturn(0);
1302 }
1303 
1304 #undef __FUNCT__
1305 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
1306 /*@
1307  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
1308 
1309    Collective
1310 
1311    Input Parameters:
1312 +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
1313 .  fetidp_flux_sol - the solution of the FETIDP linear system
1314 
1315    Output Parameters:
1316 -  standard_sol      - the solution defined on the physical domain
1317 
1318    Level: developer
1319 
1320    Notes:
1321 
1322 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1323 @*/
1324 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1325 {
1326   FETIDPMat_ctx  mat_ctx;
1327   PetscErrorCode ierr;
1328 
1329   PetscFunctionBegin;
1330   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1331   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
1332   PetscFunctionReturn(0);
1333 }
1334 /* -------------------------------------------------------------------------- */
1335 
1336 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1337 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1338 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1339 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1340 
1341 #undef __FUNCT__
1342 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
1343 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1344 {
1345 
1346   FETIDPMat_ctx  fetidpmat_ctx;
1347   Mat            newmat;
1348   FETIDPPC_ctx   fetidppc_ctx;
1349   PC             newpc;
1350   MPI_Comm       comm;
1351   PetscErrorCode ierr;
1352 
1353   PetscFunctionBegin;
1354   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1355   /* FETIDP linear matrix */
1356   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
1357   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
1358   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
1359   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1360   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
1361   ierr = MatSetUp(newmat);CHKERRQ(ierr);
1362   /* FETIDP preconditioner */
1363   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
1364   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
1365   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
1366   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
1367   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
1368   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1369   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
1370   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
1371   ierr = PCSetUp(newpc);CHKERRQ(ierr);
1372   /* return pointers for objects created */
1373   *fetidp_mat=newmat;
1374   *fetidp_pc=newpc;
1375   PetscFunctionReturn(0);
1376 }
1377 
1378 #undef __FUNCT__
1379 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
1380 /*@
1381  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
1382 
1383    Collective
1384 
1385    Input Parameters:
1386 +  pc - the BDDC preconditioning context already setup
1387 
1388    Output Parameters:
1389 .  fetidp_mat - shell FETIDP matrix object
1390 .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
1391 
1392    Options Database Keys:
1393 -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
1394 
1395    Level: developer
1396 
1397    Notes:
1398      Currently the only operation provided for FETIDP matrix is MatMult
1399 
1400 .seealso: PCBDDC
1401 @*/
1402 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1403 {
1404   PetscErrorCode ierr;
1405 
1406   PetscFunctionBegin;
1407   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1408   if (pc->setupcalled) {
1409     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1410   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1411   PetscFunctionReturn(0);
1412 }
1413 /* -------------------------------------------------------------------------- */
1414 /*MC
1415    PCBDDC - Balancing Domain Decomposition by Constraints.
1416 
1417    An implementation of the BDDC preconditioner based on
1418 
1419 .vb
1420    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
1421    [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
1422    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
1423 .ve
1424 
1425    The matrix to be preconditioned (Pmat) must be of type MATIS.
1426 
1427    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
1428 
1429    It also works with unsymmetric and indefinite problems.
1430 
1431    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.
1432 
1433    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
1434 
1435    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
1436 
1437    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace.
1438 
1439    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.
1440 
1441    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
1442 
1443    Options Database Keys:
1444 
1445 .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
1446 .    -pc_bddc_use_edges <1> - use or not edges in primal space
1447 .    -pc_bddc_use_faces <0> - use or not faces in primal space
1448 .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
1449 .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
1450 .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
1451 .    -pc_bddc_levels <0> - maximum number of levels for multilevel
1452 .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
1453 -    -pc_bddc_check_level <0> - set verbosity level of debugging output
1454 
1455    Options for Dirichlet, Neumann or coarse solver can be set with
1456 .vb
1457       -pc_bddc_dirichlet_
1458       -pc_bddc_neumann_
1459       -pc_bddc_coarse_
1460 .ve
1461    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
1462 
1463    When using a multilevel approach, solvers' options at the N-th level can be specified as
1464 .vb
1465       -pc_bddc_dirichlet_N_
1466       -pc_bddc_neumann_N_
1467       -pc_bddc_coarse_N_
1468 .ve
1469    Note that level number ranges from the finest 0 to the coarsest N
1470 
1471    Level: intermediate
1472 
1473    Developer notes:
1474      Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose
1475 
1476      New deluxe scaling operator will be available soon.
1477 
1478    Contributed by Stefano Zampini
1479 
1480 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1481 M*/
1482 
1483 #undef __FUNCT__
1484 #define __FUNCT__ "PCCreate_BDDC"
1485 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1486 {
1487   PetscErrorCode      ierr;
1488   PC_BDDC             *pcbddc;
1489 
1490   PetscFunctionBegin;
1491   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1492   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1493   pc->data  = (void*)pcbddc;
1494 
1495   /* create PCIS data structure */
1496   ierr = PCISCreate(pc);CHKERRQ(ierr);
1497 
1498   /* BDDC customization */
1499   pcbddc->use_vertices        = PETSC_TRUE;
1500   pcbddc->use_edges           = PETSC_TRUE;
1501   pcbddc->use_faces           = PETSC_FALSE;
1502   pcbddc->use_change_of_basis = PETSC_FALSE;
1503   pcbddc->use_change_on_faces = PETSC_FALSE;
1504   pcbddc->switch_static       = PETSC_FALSE;
1505   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
1506   pcbddc->dbg_flag            = 0;
1507 
1508   pcbddc->user_primal_vertices       = 0;
1509   pcbddc->NullSpace                  = 0;
1510   pcbddc->temp_solution              = 0;
1511   pcbddc->original_rhs               = 0;
1512   pcbddc->local_mat                  = 0;
1513   pcbddc->ChangeOfBasisMatrix        = 0;
1514   pcbddc->coarse_vec                 = 0;
1515   pcbddc->coarse_rhs                 = 0;
1516   pcbddc->coarse_ksp                 = 0;
1517   pcbddc->coarse_phi_B               = 0;
1518   pcbddc->coarse_phi_D               = 0;
1519   pcbddc->coarse_psi_B               = 0;
1520   pcbddc->coarse_psi_D               = 0;
1521   pcbddc->vec1_P                     = 0;
1522   pcbddc->vec1_R                     = 0;
1523   pcbddc->vec2_R                     = 0;
1524   pcbddc->local_auxmat1              = 0;
1525   pcbddc->local_auxmat2              = 0;
1526   pcbddc->R_to_B                     = 0;
1527   pcbddc->R_to_D                     = 0;
1528   pcbddc->ksp_D                      = 0;
1529   pcbddc->ksp_R                      = 0;
1530   pcbddc->NeumannBoundaries          = 0;
1531   pcbddc->NeumannBoundariesLocal     = 0;
1532   pcbddc->DirichletBoundaries        = 0;
1533   pcbddc->DirichletBoundariesLocal   = 0;
1534   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
1535   pcbddc->n_ISForDofs                = 0;
1536   pcbddc->ISForDofs                  = 0;
1537   pcbddc->ConstraintMatrix           = 0;
1538   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
1539   pcbddc->coarse_loc_to_glob         = 0;
1540   pcbddc->coarsening_ratio           = 8;
1541   pcbddc->current_level              = 0;
1542   pcbddc->max_levels                 = 0;
1543 
1544   /* create local graph structure */
1545   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1546 
1547   /* scaling */
1548   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1549   pcbddc->work_scaling               = 0;
1550 
1551   /* function pointers */
1552   pc->ops->apply               = PCApply_BDDC;
1553   pc->ops->applytranspose      = 0;
1554   pc->ops->setup               = PCSetUp_BDDC;
1555   pc->ops->destroy             = PCDestroy_BDDC;
1556   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1557   pc->ops->view                = 0;
1558   pc->ops->applyrichardson     = 0;
1559   pc->ops->applysymmetricleft  = 0;
1560   pc->ops->applysymmetricright = 0;
1561   pc->ops->presolve            = PCPreSolve_BDDC;
1562   pc->ops->postsolve           = PCPostSolve_BDDC;
1563 
1564   /* composing function */
1565   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1566   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
1567   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1568   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
1569   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1570   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1571   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1572   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1573   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1574   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1575   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1576   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1577   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1578   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1579   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1580   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1581   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1582   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1583   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1584   PetscFunctionReturn(0);
1585 }
1586 
1587