xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 63602bcaa893421da8af814d137379896d1618ad)
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__ "PCBDDCSetDofsSplittingLocal_BDDC"
720 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
721 {
722   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
723   PetscInt i;
724   PetscErrorCode ierr;
725 
726   PetscFunctionBegin;
727   /* Destroy ISes if they were already set */
728   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
729     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
730   }
731   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
732   /* last user setting takes precendence -> destroy any other customization */
733   for (i=0;i<pcbddc->n_ISForDofs;i++) {
734     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
735   }
736   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
737   pcbddc->n_ISForDofs = 0;
738   /* allocate space then set */
739   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
740   for (i=0;i<n_is;i++) {
741     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
742     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
743   }
744   pcbddc->n_ISForDofsLocal=n_is;
745   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
746   PetscFunctionReturn(0);
747 }
748 
749 #undef __FUNCT__
750 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
751 /*@
752  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
753 
754    Collective
755 
756    Input Parameters:
757 +  pc - the preconditioning context
758 -  n_is - number of index sets defining the fields
759 .  ISForDofs - array of IS describing the fields in local ordering
760 
761    Level: intermediate
762 
763    Notes: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to a different field.
764 
765 .seealso: PCBDDC
766 @*/
767 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
768 {
769   PetscInt       i;
770   PetscErrorCode ierr;
771 
772   PetscFunctionBegin;
773   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
774   PetscValidLogicalCollectiveInt(pc,n_is,2);
775   for (i=0;i<n_is;i++) {
776     PetscCheckSameComm(pc,1,ISForDofs[i],3);
777     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
778   }
779   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
780   PetscFunctionReturn(0);
781 }
782 /* -------------------------------------------------------------------------- */
783 
784 #undef __FUNCT__
785 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
786 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
787 {
788   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
789   PetscInt i;
790   PetscErrorCode ierr;
791 
792   PetscFunctionBegin;
793   /* Destroy ISes if they were already set */
794   for (i=0;i<pcbddc->n_ISForDofs;i++) {
795     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
796   }
797   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
798   /* last user setting takes precendence -> destroy any other customization */
799   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
800     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
801   }
802   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
803   pcbddc->n_ISForDofsLocal = 0;
804   /* allocate space then set */
805   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
806   for (i=0;i<n_is;i++) {
807     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
808     pcbddc->ISForDofs[i]=ISForDofs[i];
809   }
810   pcbddc->n_ISForDofs=n_is;
811   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
812   PetscFunctionReturn(0);
813 }
814 
815 #undef __FUNCT__
816 #define __FUNCT__ "PCBDDCSetDofsSplitting"
817 /*@
818  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
819 
820    Collective
821 
822    Input Parameters:
823 +  pc - the preconditioning context
824 -  n_is - number of index sets defining the fields
825 .  ISForDofs - array of IS describing the fields in global ordering
826 
827    Level: intermediate
828 
829    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.
830 
831 .seealso: PCBDDC
832 @*/
833 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
834 {
835   PetscInt       i;
836   PetscErrorCode ierr;
837 
838   PetscFunctionBegin;
839   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
840   PetscValidLogicalCollectiveInt(pc,n_is,2);
841   for (i=0;i<n_is;i++) {
842     PetscCheckSameComm(pc,1,ISForDofs[i],3);
843     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
844   }
845   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
846   PetscFunctionReturn(0);
847 }
848 /* -------------------------------------------------------------------------- */
849 #undef __FUNCT__
850 #define __FUNCT__ "PCPreSolve_BDDC"
851 /* -------------------------------------------------------------------------- */
852 /*
853    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
854                      guess if a transformation of basis approach has been selected.
855 
856    Input Parameter:
857 +  pc - the preconditioner contex
858 
859    Application Interface Routine: PCPreSolve()
860 
861    Notes:
862    The interface routine PCPreSolve() is not usually called directly by
863    the user, but instead is called by KSPSolve().
864 */
865 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
866 {
867   PetscErrorCode ierr;
868   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
869   PC_IS          *pcis = (PC_IS*)(pc->data);
870   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
871   Mat            temp_mat;
872   IS             dirIS;
873   Vec            used_vec;
874   PetscBool      guess_nonzero;
875 
876   PetscFunctionBegin;
877   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
878   if (ksp) {
879     PetscBool iscg;
880     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
881     if (!iscg) {
882       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
883     }
884   }
885   /* Creates parallel work vectors used in presolve */
886   if (!pcbddc->original_rhs) {
887     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
888   }
889   if (!pcbddc->temp_solution) {
890     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
891   }
892   if (x) {
893     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
894     used_vec = x;
895   } else {
896     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
897     used_vec = pcbddc->temp_solution;
898     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
899   }
900   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
901   if (ksp) {
902     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
903     if (!guess_nonzero) {
904       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
905     }
906   }
907 
908   /* store the original rhs */
909   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
910 
911   /* Take into account zeroed rows -> change rhs and store solution removed */
912   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
913   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
914   if (rhs && dirIS) {
915     PetscInt    dirsize,i,*is_indices;
916     PetscScalar *array_x,*array_diagonal;
917 
918     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
919     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
920     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
921     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
922     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
923     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
924     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
925     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
926     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
927     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
928     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
929     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
930     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
931     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
932     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
933     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
934 
935     /* remove the computed solution from the rhs */
936     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
937     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
938     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
939   }
940 
941   /* store partially computed solution and set initial guess */
942   if (x) {
943     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
944     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
945     if (pcbddc->use_exact_dirichlet_trick) {
946       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
947       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
948       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
949       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
950       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
951       if (ksp) {
952         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
953       }
954     }
955   }
956 
957   /* prepare MatMult and rhs for solver */
958   if (pcbddc->use_change_of_basis) {
959     /* swap pointers for local matrices */
960     temp_mat = matis->A;
961     matis->A = pcbddc->local_mat;
962     pcbddc->local_mat = temp_mat;
963     if (rhs) {
964       /* Get local rhs and apply transformation of basis */
965       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
966       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
967       /* from original basis to modified basis */
968       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
969       /* put back modified values into the global vec using INSERT_VALUES copy mode */
970       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
971       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
972     }
973   }
974 
975   /* remove nullspace if present */
976   if (ksp && pcbddc->NullSpace) {
977     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
978     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
979   }
980   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
981   PetscFunctionReturn(0);
982 }
983 /* -------------------------------------------------------------------------- */
984 #undef __FUNCT__
985 #define __FUNCT__ "PCPostSolve_BDDC"
986 /* -------------------------------------------------------------------------- */
987 /*
988    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
989                      approach has been selected. Also, restores rhs to its original state.
990 
991    Input Parameter:
992 +  pc - the preconditioner contex
993 
994    Application Interface Routine: PCPostSolve()
995 
996    Notes:
997    The interface routine PCPostSolve() is not usually called directly by
998    the user, but instead is called by KSPSolve().
999 */
1000 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1001 {
1002   PetscErrorCode ierr;
1003   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1004   PC_IS          *pcis   = (PC_IS*)(pc->data);
1005   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
1006   Mat            temp_mat;
1007 
1008   PetscFunctionBegin;
1009   if (pcbddc->use_change_of_basis) {
1010     /* swap pointers for local matrices */
1011     temp_mat = matis->A;
1012     matis->A = pcbddc->local_mat;
1013     pcbddc->local_mat = temp_mat;
1014   }
1015   if (pcbddc->use_change_of_basis && x) {
1016     /* Get Local boundary and apply transformation of basis to solution vector */
1017     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1018     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1019     /* from modified basis to original basis */
1020     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
1021     /* put back modified values into the global vec using INSERT_VALUES copy mode */
1022     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1023     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1024   }
1025   /* add solution removed in presolve */
1026   if (x) {
1027     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
1028   }
1029   /* restore rhs to its original state */
1030   if (rhs) {
1031     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1032   }
1033   PetscFunctionReturn(0);
1034 }
1035 /* -------------------------------------------------------------------------- */
1036 #undef __FUNCT__
1037 #define __FUNCT__ "PCSetUp_BDDC"
1038 /* -------------------------------------------------------------------------- */
1039 /*
1040    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1041                   by setting data structures and options.
1042 
1043    Input Parameter:
1044 +  pc - the preconditioner context
1045 
1046    Application Interface Routine: PCSetUp()
1047 
1048    Notes:
1049    The interface routine PCSetUp() is not usually called directly by
1050    the user, but instead is called by PCApply() if necessary.
1051 */
1052 PetscErrorCode PCSetUp_BDDC(PC pc)
1053 {
1054   PetscErrorCode ierr;
1055   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
1056   MatStructure   flag;
1057   PetscBool      computeis,computetopography,computesolvers;
1058 
1059   PetscFunctionBegin;
1060   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */
1061   /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */
1062   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1063      Also, BDDC directly build the Dirichlet problem */
1064   /* Get stdout for dbg */
1065   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
1066     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
1067     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1068     if (pcbddc->current_level) {
1069       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
1070     }
1071   }
1072   /* first attempt to split work */
1073   if (pc->setupcalled) {
1074     computeis = PETSC_FALSE;
1075     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
1076     if (flag == SAME_PRECONDITIONER) {
1077       computetopography = PETSC_FALSE;
1078       computesolvers = PETSC_FALSE;
1079     } else if (flag == SAME_NONZERO_PATTERN) {
1080       computetopography = PETSC_FALSE;
1081       computesolvers = PETSC_TRUE;
1082     } else { /* DIFFERENT_NONZERO_PATTERN */
1083       computetopography = PETSC_TRUE;
1084       computesolvers = PETSC_TRUE;
1085     }
1086   } else {
1087     computeis = PETSC_TRUE;
1088     computetopography = PETSC_TRUE;
1089     computesolvers = PETSC_TRUE;
1090   }
1091   /* Set up all the "iterative substructuring" common block without computing solvers */
1092   if (computeis) {
1093     /* HACK INTO PCIS */
1094     PC_IS* pcis = (PC_IS*)pc->data;
1095     pcis->computesolvers = PETSC_FALSE;
1096     ierr = PCISSetUp(pc);CHKERRQ(ierr);
1097   }
1098   /* Analyze interface and set up local constraint and change of basis matrices */
1099   if (computetopography) {
1100     /* reset data */
1101     ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1102     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1103     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1104   }
1105   if (computesolvers) {
1106     /* reset data */
1107     ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1108     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1109     /* Create coarse and local stuffs */
1110     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1111     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
1112   }
1113   if (pcbddc->dbg_flag && pcbddc->current_level) {
1114     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
1115   }
1116   PetscFunctionReturn(0);
1117 }
1118 
1119 /* -------------------------------------------------------------------------- */
1120 /*
1121    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
1122 
1123    Input Parameters:
1124 .  pc - the preconditioner context
1125 .  r - input vector (global)
1126 
1127    Output Parameter:
1128 .  z - output vector (global)
1129 
1130    Application Interface Routine: PCApply()
1131  */
1132 #undef __FUNCT__
1133 #define __FUNCT__ "PCApply_BDDC"
1134 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1135 {
1136   PC_IS             *pcis = (PC_IS*)(pc->data);
1137   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1138   PetscErrorCode    ierr;
1139   const PetscScalar one = 1.0;
1140   const PetscScalar m_one = -1.0;
1141   const PetscScalar zero = 0.0;
1142 
1143 /* This code is similar to that provided in nn.c for PCNN
1144    NN interface preconditioner changed to BDDC
1145    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
1146 
1147   PetscFunctionBegin;
1148   if (!pcbddc->use_exact_dirichlet_trick) {
1149     /* First Dirichlet solve */
1150     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1151     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1152     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1153     /*
1154       Assembling right hand side for BDDC operator
1155       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1156       - pcis->vec1_B the interface part of the global vector z
1157     */
1158     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1159     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1160     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1161     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1162     ierr = VecCopy(r,z);CHKERRQ(ierr);
1163     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1164     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1165     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1166   } else {
1167     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1168     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1169     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1170   }
1171 
1172   /* Apply interface preconditioner
1173      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1174   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
1175 
1176   /* Apply transpose of partition of unity operator */
1177   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1178 
1179   /* Second Dirichlet solve and assembling of output */
1180   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1181   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1182   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1183   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1184   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr);
1185   ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr);
1186   if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1187   ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr);
1188   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1189   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1190   PetscFunctionReturn(0);
1191 }
1192 /* -------------------------------------------------------------------------- */
1193 
1194 #undef __FUNCT__
1195 #define __FUNCT__ "PCDestroy_BDDC"
1196 PetscErrorCode PCDestroy_BDDC(PC pc)
1197 {
1198   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1199   PetscErrorCode ierr;
1200 
1201   PetscFunctionBegin;
1202   /* free data created by PCIS */
1203   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1204   /* free BDDC custom data  */
1205   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1206   /* destroy objects related to topography */
1207   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1208   /* free allocated graph structure */
1209   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1210   /* free data for scaling operator */
1211   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1212   /* free solvers stuff */
1213   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1214   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
1215   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
1216   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
1217   ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1218   /* free global vectors needed in presolve */
1219   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
1220   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1221   /* remove functions */
1222   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1223   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
1224   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1225   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
1226   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1227   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1228   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1229   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1230   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1231   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1232   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1233   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1234   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1235   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1236   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1237   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1238   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1239   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1240   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1241   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1242   /* Free the private data structure */
1243   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1244   PetscFunctionReturn(0);
1245 }
1246 /* -------------------------------------------------------------------------- */
1247 
1248 #undef __FUNCT__
1249 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
1250 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1251 {
1252   FETIDPMat_ctx  mat_ctx;
1253   PC_IS*         pcis;
1254   PC_BDDC*       pcbddc;
1255   PetscErrorCode ierr;
1256 
1257   PetscFunctionBegin;
1258   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1259   pcis = (PC_IS*)mat_ctx->pc->data;
1260   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1261 
1262   /* change of basis for physical rhs if needed
1263      It also changes the rhs in case of dirichlet boundaries */
1264   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
1265   /* store vectors for computation of fetidp final solution */
1266   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1267   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1268   /* scale rhs since it should be unassembled */
1269   /* TODO use counter scaling? (also below) */
1270   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1271   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1272   /* Apply partition of unity */
1273   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1274   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1275   if (!pcbddc->switch_static) {
1276     /* compute partially subassembled Schur complement right-hand side */
1277     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1278     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
1279     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1280     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
1281     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1282     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1283     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1284     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1285     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1286     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1287   }
1288   /* BDDC rhs */
1289   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
1290   if (pcbddc->switch_static) {
1291     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1292   }
1293   /* apply BDDC */
1294   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
1295   /* Application of B_delta and assembling of rhs for fetidp fluxes */
1296   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
1297   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
1298   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1299   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1300   /* restore original rhs */
1301   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
1302   PetscFunctionReturn(0);
1303 }
1304 
1305 #undef __FUNCT__
1306 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
1307 /*@
1308  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
1309 
1310    Collective
1311 
1312    Input Parameters:
1313 +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
1314 .  standard_rhs - the right-hand side for your linear system
1315 
1316    Output Parameters:
1317 -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
1318 
1319    Level: developer
1320 
1321    Notes:
1322 
1323 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1324 @*/
1325 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1326 {
1327   FETIDPMat_ctx  mat_ctx;
1328   PetscErrorCode ierr;
1329 
1330   PetscFunctionBegin;
1331   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1332   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
1333   PetscFunctionReturn(0);
1334 }
1335 /* -------------------------------------------------------------------------- */
1336 
1337 #undef __FUNCT__
1338 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
1339 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1340 {
1341   FETIDPMat_ctx  mat_ctx;
1342   PC_IS*         pcis;
1343   PC_BDDC*       pcbddc;
1344   PetscErrorCode ierr;
1345 
1346   PetscFunctionBegin;
1347   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1348   pcis = (PC_IS*)mat_ctx->pc->data;
1349   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1350 
1351   /* apply B_delta^T */
1352   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1353   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1354   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
1355   /* compute rhs for BDDC application */
1356   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1357   if (pcbddc->switch_static) {
1358     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1359   }
1360   /* apply BDDC */
1361   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
1362   /* put values into standard global vector */
1363   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1364   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1365   if (!pcbddc->switch_static) {
1366     /* compute values into the interior if solved for the partially subassembled Schur complement */
1367     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
1368     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
1369     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1370   }
1371   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1372   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1373   /* final change of basis if needed
1374      Is also sums the dirichlet part removed during RHS assembling */
1375   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
1376   PetscFunctionReturn(0);
1377 }
1378 
1379 #undef __FUNCT__
1380 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
1381 /*@
1382  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
1383 
1384    Collective
1385 
1386    Input Parameters:
1387 +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
1388 .  fetidp_flux_sol - the solution of the FETIDP linear system
1389 
1390    Output Parameters:
1391 -  standard_sol      - the solution defined on the physical domain
1392 
1393    Level: developer
1394 
1395    Notes:
1396 
1397 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1398 @*/
1399 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1400 {
1401   FETIDPMat_ctx  mat_ctx;
1402   PetscErrorCode ierr;
1403 
1404   PetscFunctionBegin;
1405   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1406   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
1407   PetscFunctionReturn(0);
1408 }
1409 /* -------------------------------------------------------------------------- */
1410 
1411 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1412 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1413 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1414 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1415 
1416 #undef __FUNCT__
1417 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
1418 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1419 {
1420 
1421   FETIDPMat_ctx  fetidpmat_ctx;
1422   Mat            newmat;
1423   FETIDPPC_ctx   fetidppc_ctx;
1424   PC             newpc;
1425   MPI_Comm       comm;
1426   PetscErrorCode ierr;
1427 
1428   PetscFunctionBegin;
1429   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1430   /* FETIDP linear matrix */
1431   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
1432   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
1433   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
1434   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1435   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
1436   ierr = MatSetUp(newmat);CHKERRQ(ierr);
1437   /* FETIDP preconditioner */
1438   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
1439   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
1440   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
1441   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
1442   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
1443   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1444   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
1445   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
1446   ierr = PCSetUp(newpc);CHKERRQ(ierr);
1447   /* return pointers for objects created */
1448   *fetidp_mat=newmat;
1449   *fetidp_pc=newpc;
1450   PetscFunctionReturn(0);
1451 }
1452 
1453 #undef __FUNCT__
1454 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
1455 /*@
1456  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
1457 
1458    Collective
1459 
1460    Input Parameters:
1461 +  pc - the BDDC preconditioning context already setup
1462 
1463    Output Parameters:
1464 .  fetidp_mat - shell FETIDP matrix object
1465 .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
1466 
1467    Options Database Keys:
1468 -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
1469 
1470    Level: developer
1471 
1472    Notes:
1473      Currently the only operation provided for FETIDP matrix is MatMult
1474 
1475 .seealso: PCBDDC
1476 @*/
1477 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1478 {
1479   PetscErrorCode ierr;
1480 
1481   PetscFunctionBegin;
1482   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1483   if (pc->setupcalled) {
1484     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1485   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1486   PetscFunctionReturn(0);
1487 }
1488 /* -------------------------------------------------------------------------- */
1489 /*MC
1490    PCBDDC - Balancing Domain Decomposition by Constraints.
1491 
1492    An implementation of the BDDC preconditioner based on
1493 
1494 .vb
1495    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
1496    [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
1497    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
1498 .ve
1499 
1500    The matrix to be preconditioned (Pmat) must be of type MATIS.
1501 
1502    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
1503 
1504    It also works with unsymmetric and indefinite problems.
1505 
1506    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.
1507 
1508    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
1509 
1510    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
1511 
1512    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace.
1513 
1514    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.
1515 
1516    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
1517 
1518    Options Database Keys:
1519 
1520 .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
1521 .    -pc_bddc_use_edges <1> - use or not edges in primal space
1522 .    -pc_bddc_use_faces <0> - use or not faces in primal space
1523 .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
1524 .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
1525 .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
1526 .    -pc_bddc_levels <0> - maximum number of levels for multilevel
1527 .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
1528 -    -pc_bddc_check_level <0> - set verbosity level of debugging output
1529 
1530    Options for Dirichlet, Neumann or coarse solver can be set with
1531 .vb
1532       -pc_bddc_dirichlet_
1533       -pc_bddc_neumann_
1534       -pc_bddc_coarse_
1535 .ve
1536    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
1537 
1538    When using a multilevel approach, solvers' options at the N-th level can be specified as
1539 .vb
1540       -pc_bddc_dirichlet_N_
1541       -pc_bddc_neumann_N_
1542       -pc_bddc_coarse_N_
1543 .ve
1544    Note that level number ranges from the finest 0 to the coarsest N
1545 
1546    Level: intermediate
1547 
1548    Developer notes:
1549      Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose
1550 
1551      New deluxe scaling operator will be available soon.
1552 
1553    Contributed by Stefano Zampini
1554 
1555 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1556 M*/
1557 
1558 #undef __FUNCT__
1559 #define __FUNCT__ "PCCreate_BDDC"
1560 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1561 {
1562   PetscErrorCode      ierr;
1563   PC_BDDC             *pcbddc;
1564 
1565   PetscFunctionBegin;
1566   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1567   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1568   pc->data  = (void*)pcbddc;
1569 
1570   /* create PCIS data structure */
1571   ierr = PCISCreate(pc);CHKERRQ(ierr);
1572 
1573   /* BDDC customization */
1574   pcbddc->use_vertices        = PETSC_TRUE;
1575   pcbddc->use_edges           = PETSC_TRUE;
1576   pcbddc->use_faces           = PETSC_FALSE;
1577   pcbddc->use_change_of_basis = PETSC_FALSE;
1578   pcbddc->use_change_on_faces = PETSC_FALSE;
1579   pcbddc->switch_static       = PETSC_FALSE;
1580   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
1581   pcbddc->dbg_flag            = 0;
1582 
1583   pcbddc->user_primal_vertices       = 0;
1584   pcbddc->NullSpace                  = 0;
1585   pcbddc->temp_solution              = 0;
1586   pcbddc->original_rhs               = 0;
1587   pcbddc->local_mat                  = 0;
1588   pcbddc->ChangeOfBasisMatrix        = 0;
1589   pcbddc->coarse_vec                 = 0;
1590   pcbddc->coarse_rhs                 = 0;
1591   pcbddc->coarse_ksp                 = 0;
1592   pcbddc->coarse_phi_B               = 0;
1593   pcbddc->coarse_phi_D               = 0;
1594   pcbddc->coarse_psi_B               = 0;
1595   pcbddc->coarse_psi_D               = 0;
1596   pcbddc->vec1_P                     = 0;
1597   pcbddc->vec1_R                     = 0;
1598   pcbddc->vec2_R                     = 0;
1599   pcbddc->local_auxmat1              = 0;
1600   pcbddc->local_auxmat2              = 0;
1601   pcbddc->R_to_B                     = 0;
1602   pcbddc->R_to_D                     = 0;
1603   pcbddc->ksp_D                      = 0;
1604   pcbddc->ksp_R                      = 0;
1605   pcbddc->NeumannBoundaries          = 0;
1606   pcbddc->NeumannBoundariesLocal     = 0;
1607   pcbddc->DirichletBoundaries        = 0;
1608   pcbddc->DirichletBoundariesLocal   = 0;
1609   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
1610   pcbddc->n_ISForDofs                = 0;
1611   pcbddc->n_ISForDofsLocal           = 0;
1612   pcbddc->ISForDofs                  = 0;
1613   pcbddc->ISForDofsLocal             = 0;
1614   pcbddc->ConstraintMatrix           = 0;
1615   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
1616   pcbddc->coarse_loc_to_glob         = 0;
1617   pcbddc->coarsening_ratio           = 8;
1618   pcbddc->current_level              = 0;
1619   pcbddc->max_levels                 = 0;
1620 
1621   /* create local graph structure */
1622   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1623 
1624   /* scaling */
1625   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1626   pcbddc->work_scaling               = 0;
1627 
1628   /* function pointers */
1629   pc->ops->apply               = PCApply_BDDC;
1630   pc->ops->applytranspose      = 0;
1631   pc->ops->setup               = PCSetUp_BDDC;
1632   pc->ops->destroy             = PCDestroy_BDDC;
1633   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1634   pc->ops->view                = 0;
1635   pc->ops->applyrichardson     = 0;
1636   pc->ops->applysymmetricleft  = 0;
1637   pc->ops->applysymmetricright = 0;
1638   pc->ops->presolve            = PCPreSolve_BDDC;
1639   pc->ops->postsolve           = PCPostSolve_BDDC;
1640 
1641   /* composing function */
1642   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1643   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
1644   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1645   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
1646   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1647   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1648   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1649   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1650   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1651   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1652   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1653   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1654   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1655   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1656   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1657   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
1658   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1659   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1660   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1661   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1662   PetscFunctionReturn(0);
1663 }
1664 
1665