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