xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 727cdba63cef3c72c66ef27a9d74c90cecc5d6d6)
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    - MAT_IGNORE_ZERO_ENTRIES for Constraints Matrix
7 
8    Solvers
9    - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers)
10    - reuse already allocated coarse matrix if possible
11    - Propagate ksp prefixes for solvers to mat objects?
12    - Propagate nearnullspace info among levels
13 
14    User interface
15    - Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet)
16    - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access)
17    - Provide PCApplyTranpose_BDDC
18    - DofSplitting and DM attached to pc?
19 
20    Debugging output
21    - Better management of verbosity levels of debugging output
22    - Crashes on some architecture -> call SynchronizedAllow before every SynchronizedPrintf
23 
24    Build
25    - make runexe59
26 
27    Extra
28    - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
29    - Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver?
30    - add support for computing h,H and related using coordinates?
31    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
32    - Better management in PCIS code
33    - BDDC with MG framework?
34 
35    FETIDP
36    - Move FETIDP code to its own classes
37 
38    MATIS related operations contained in BDDC code
39    - Add MAT_REUSE in MatConvert_IS_AIJ
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   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
311   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
312   pcbddc->DirichletBoundaries=DirichletBoundaries;
313   pcbddc->recompute_topography = PETSC_TRUE;
314   PetscFunctionReturn(0);
315 }
316 
317 #undef __FUNCT__
318 #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
319 /*@
320  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
321 
322    Not collective
323 
324    Input Parameters:
325 +  pc - the preconditioning context
326 -  DirichletBoundaries - sequential IS defining the subdomain part of Dirichlet boundaries (in local ordering)
327 
328    Level: intermediate
329 
330    Notes:
331 
332 .seealso: PCBDDC
333 @*/
334 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
335 {
336   PetscErrorCode ierr;
337 
338   PetscFunctionBegin;
339   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
340   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
341   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
342   PetscFunctionReturn(0);
343 }
344 /* -------------------------------------------------------------------------- */
345 
346 #undef __FUNCT__
347 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
348 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
349 {
350   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
351   PetscErrorCode ierr;
352 
353   PetscFunctionBegin;
354   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
355   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
356   pcbddc->NeumannBoundaries=NeumannBoundaries;
357   pcbddc->recompute_topography = PETSC_TRUE;
358   PetscFunctionReturn(0);
359 }
360 
361 #undef __FUNCT__
362 #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
363 /*@
364  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
365 
366    Not collective
367 
368    Input Parameters:
369 +  pc - the preconditioning context
370 -  NeumannBoundaries - sequential 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 PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
379 {
380   PetscErrorCode ierr;
381 
382   PetscFunctionBegin;
383   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
384   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
385   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
386   PetscFunctionReturn(0);
387 }
388 /* -------------------------------------------------------------------------- */
389 
390 #undef __FUNCT__
391 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
392 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
393 {
394   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
395 
396   PetscFunctionBegin;
397   *DirichletBoundaries = pcbddc->DirichletBoundaries;
398   PetscFunctionReturn(0);
399 }
400 
401 #undef __FUNCT__
402 #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
403 /*@
404  PCBDDCGetDirichletBoundaries - Get IS for local Dirichlet boundaries
405 
406    Not collective
407 
408    Input Parameters:
409 .  pc - the preconditioning context
410 
411    Output Parameters:
412 .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
413 
414    Level: intermediate
415 
416    Notes:
417 
418 .seealso: PCBDDC
419 @*/
420 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
421 {
422   PetscErrorCode ierr;
423 
424   PetscFunctionBegin;
425   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
426   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
427   PetscFunctionReturn(0);
428 }
429 /* -------------------------------------------------------------------------- */
430 
431 #undef __FUNCT__
432 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
433 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
434 {
435   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
436 
437   PetscFunctionBegin;
438   *NeumannBoundaries = pcbddc->NeumannBoundaries;
439   PetscFunctionReturn(0);
440 }
441 
442 #undef __FUNCT__
443 #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
444 /*@
445  PCBDDCGetNeumannBoundaries - Get IS for local Neumann boundaries
446 
447    Not collective
448 
449    Input Parameters:
450 .  pc - the preconditioning context
451 
452    Output Parameters:
453 .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
454 
455    Level: intermediate
456 
457    Notes:
458 
459 .seealso: PCBDDC
460 @*/
461 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
462 {
463   PetscErrorCode ierr;
464 
465   PetscFunctionBegin;
466   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
467   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
468   PetscFunctionReturn(0);
469 }
470 /* -------------------------------------------------------------------------- */
471 
472 #undef __FUNCT__
473 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
474 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
475 {
476   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
477   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
478   PetscErrorCode ierr;
479 
480   PetscFunctionBegin;
481   /* free old CSR */
482   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
483   /* TODO: PCBDDCGraphSetAdjacency */
484   /* get CSR into graph structure */
485   if (copymode == PETSC_COPY_VALUES) {
486     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
487     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
488     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
489     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
490   } else if (copymode == PETSC_OWN_POINTER) {
491     mat_graph->xadj = (PetscInt*)xadj;
492     mat_graph->adjncy = (PetscInt*)adjncy;
493   }
494   mat_graph->nvtxs_csr = nvtxs;
495   PetscFunctionReturn(0);
496 }
497 
498 #undef __FUNCT__
499 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
500 /*@
501  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
502 
503    Not collective
504 
505    Input Parameters:
506 +  pc - the preconditioning context
507 .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
508 .  xadj, adjncy - the CSR graph
509 -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
510 
511    Level: intermediate
512 
513    Notes:
514 
515 .seealso: PCBDDC,PetscCopyMode
516 @*/
517 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
518 {
519   void (*f)(void) = 0;
520   PetscErrorCode ierr;
521 
522   PetscFunctionBegin;
523   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
524   PetscValidIntPointer(xadj,3);
525   PetscValidIntPointer(xadj,4);
526   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
527     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
528   }
529   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
530   /* free arrays if PCBDDC is not the PC type */
531   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
532   if (!f && copymode == PETSC_OWN_POINTER) {
533     ierr = PetscFree(xadj);CHKERRQ(ierr);
534     ierr = PetscFree(adjncy);CHKERRQ(ierr);
535   }
536   PetscFunctionReturn(0);
537 }
538 /* -------------------------------------------------------------------------- */
539 
540 #undef __FUNCT__
541 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
542 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
543 {
544   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
545   PetscInt i;
546   PetscErrorCode ierr;
547 
548   PetscFunctionBegin;
549   /* Destroy ISes if they were already set */
550   for (i=0;i<pcbddc->n_ISForDofs;i++) {
551     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
552   }
553   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
554   /* allocate space then set */
555   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
556   for (i=0;i<n_is;i++) {
557     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
558     pcbddc->ISForDofs[i]=ISForDofs[i];
559   }
560   pcbddc->n_ISForDofs=n_is;
561   pcbddc->user_provided_isfordofs = PETSC_TRUE;
562   PetscFunctionReturn(0);
563 }
564 
565 #undef __FUNCT__
566 #define __FUNCT__ "PCBDDCSetDofsSplitting"
567 /*@
568  PCBDDCSetDofsSplitting - Set index sets defining fields of the local Neumann matrix
569 
570    Not collective
571 
572    Input Parameters:
573 +  pc - the preconditioning context
574 -  n_is - number of index sets defining the fields
575 .  ISForDofs - array of IS describing the fields
576 
577    Level: intermediate
578 
579    Notes:
580 
581 .seealso: PCBDDC
582 @*/
583 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
584 {
585   PetscInt       i;
586   PetscErrorCode ierr;
587 
588   PetscFunctionBegin;
589   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
590   for (i=0;i<n_is;i++) {
591     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2);
592   }
593   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
594   PetscFunctionReturn(0);
595 }
596 /* -------------------------------------------------------------------------- */
597 #undef __FUNCT__
598 #define __FUNCT__ "PCPreSolve_BDDC"
599 /* -------------------------------------------------------------------------- */
600 /*
601    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
602                      guess if a transformation of basis approach has been selected.
603 
604    Input Parameter:
605 +  pc - the preconditioner contex
606 
607    Application Interface Routine: PCPreSolve()
608 
609    Notes:
610    The interface routine PCPreSolve() is not usually called directly by
611    the user, but instead is called by KSPSolve().
612 */
613 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
614 {
615   PetscErrorCode ierr;
616   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
617   PC_IS          *pcis = (PC_IS*)(pc->data);
618   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
619   Mat            temp_mat;
620   IS             dirIS;
621   PetscInt       dirsize,i,*is_indices;
622   PetscScalar    *array_x,*array_diagonal;
623   Vec            used_vec;
624   PetscBool      guess_nonzero,flg,bddc_has_dirichlet_boundaries;
625 
626   PetscFunctionBegin;
627   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
628   if (ksp) {
629     PetscBool iscg;
630     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
631     if (!iscg) {
632       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
633     }
634   }
635   /* Creates parallel work vectors used in presolve */
636   if (!pcbddc->original_rhs) {
637     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
638   }
639   if (!pcbddc->temp_solution) {
640     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
641   }
642   if (x) {
643     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
644     used_vec = x;
645   } else {
646     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
647     used_vec = pcbddc->temp_solution;
648     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
649   }
650   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
651   if (ksp) {
652     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
653     if (!guess_nonzero) {
654       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
655     }
656   }
657 
658   /* TODO: remove when Dirichlet boundaries will be shared */
659   ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr);
660   flg = PETSC_FALSE;
661   if (dirIS) flg = PETSC_TRUE;
662   ierr = MPI_Allreduce(&flg,&bddc_has_dirichlet_boundaries,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
663 
664   /* store the original rhs */
665   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
666 
667   /* Take into account zeroed rows -> change rhs and store solution removed */
668   if (rhs && bddc_has_dirichlet_boundaries) {
669     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
670     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
671     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
672     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
673     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
674     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
675     if (dirIS) {
676       ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr);
677       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
678       ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
679       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
680       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
681       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
682       ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
683       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
684     }
685     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
686     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
687 
688     /* remove the computed solution from the rhs */
689     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
690     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
691     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
692   }
693 
694   /* store partially computed solution and set initial guess */
695   if (x) {
696     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
697     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
698     if (pcbddc->use_exact_dirichlet_trick) {
699       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
700       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
701       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
702       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
703       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
704       if (ksp) {
705         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
706       }
707     }
708   }
709 
710   /* prepare MatMult and rhs for solver */
711   if (pcbddc->use_change_of_basis) {
712     /* swap pointers for local matrices */
713     temp_mat = matis->A;
714     matis->A = pcbddc->local_mat;
715     pcbddc->local_mat = temp_mat;
716     if (rhs) {
717       /* Get local rhs and apply transformation of basis */
718       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
719       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
720       /* from original basis to modified basis */
721       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
722       /* put back modified values into the global vec using INSERT_VALUES copy mode */
723       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
724       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
725     }
726   }
727 
728   /* remove nullspace if present */
729   if (ksp && pcbddc->NullSpace) {
730     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
731     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
732   }
733   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
734   PetscFunctionReturn(0);
735 }
736 /* -------------------------------------------------------------------------- */
737 #undef __FUNCT__
738 #define __FUNCT__ "PCPostSolve_BDDC"
739 /* -------------------------------------------------------------------------- */
740 /*
741    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
742                      approach has been selected. Also, restores rhs to its original state.
743 
744    Input Parameter:
745 +  pc - the preconditioner contex
746 
747    Application Interface Routine: PCPostSolve()
748 
749    Notes:
750    The interface routine PCPostSolve() is not usually called directly by
751    the user, but instead is called by KSPSolve().
752 */
753 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
754 {
755   PetscErrorCode ierr;
756   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
757   PC_IS          *pcis   = (PC_IS*)(pc->data);
758   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
759   Mat            temp_mat;
760 
761   PetscFunctionBegin;
762   if (pcbddc->use_change_of_basis) {
763     /* swap pointers for local matrices */
764     temp_mat = matis->A;
765     matis->A = pcbddc->local_mat;
766     pcbddc->local_mat = temp_mat;
767   }
768   if (pcbddc->use_change_of_basis && x) {
769     /* Get Local boundary and apply transformation of basis to solution vector */
770     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
771     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
772     /* from modified basis to original basis */
773     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
774     /* put back modified values into the global vec using INSERT_VALUES copy mode */
775     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
776     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
777   }
778   /* add solution removed in presolve */
779   if (x) {
780     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
781   }
782   /* restore rhs to its original state */
783   if (rhs) {
784     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
785   }
786   PetscFunctionReturn(0);
787 }
788 /* -------------------------------------------------------------------------- */
789 #undef __FUNCT__
790 #define __FUNCT__ "PCSetUp_BDDC"
791 /* -------------------------------------------------------------------------- */
792 /*
793    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
794                   by setting data structures and options.
795 
796    Input Parameter:
797 +  pc - the preconditioner context
798 
799    Application Interface Routine: PCSetUp()
800 
801    Notes:
802    The interface routine PCSetUp() is not usually called directly by
803    the user, but instead is called by PCApply() if necessary.
804 */
805 PetscErrorCode PCSetUp_BDDC(PC pc)
806 {
807   PetscErrorCode   ierr;
808   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
809   MatNullSpace     nearnullspace;
810   const Vec        *nearnullvecs,*onearnullvecs;
811   MatStructure     flag;
812   PetscObjectState state;
813   PetscInt         nnsp_size,onnsp_size;
814   PetscBool        computeis,computetopography,computesolvers;
815   PetscBool        new_nearnullspace_provided,nnsp_has_cnst,onnsp_has_cnst;
816 
817   PetscFunctionBegin;
818   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
819   /* PCIS does not support MatStructure flags different from SAME_PRECONDITIONER */
820   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
821      Also, BDDC directly build the Dirichlet problem */
822 
823   /* split work */
824   if (pc->setupcalled) {
825     computeis = PETSC_FALSE;
826     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
827     if (flag == SAME_PRECONDITIONER) {
828       PetscFunctionReturn(0);
829     } else if (flag == SAME_NONZERO_PATTERN) {
830       computetopography = PETSC_FALSE;
831       computesolvers = PETSC_TRUE;
832     } else { /* DIFFERENT_NONZERO_PATTERN */
833       computetopography = PETSC_TRUE;
834       computesolvers = PETSC_TRUE;
835     }
836   } else {
837     computeis = PETSC_TRUE;
838     computetopography = PETSC_TRUE;
839     computesolvers = PETSC_TRUE;
840   }
841   if (pcbddc->recompute_topography) {
842     computetopography = PETSC_TRUE;
843   }
844 
845   /* Get stdout for dbg */
846   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
847     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
848     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
849     if (pcbddc->current_level) {
850       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
851     }
852   }
853 
854   /* Set up all the "iterative substructuring" common block without computing solvers */
855   if (computeis) {
856     /* HACK INTO PCIS */
857     PC_IS* pcis = (PC_IS*)pc->data;
858     pcis->computesolvers = PETSC_FALSE;
859     ierr = PCISSetUp(pc);CHKERRQ(ierr);
860   }
861 
862   /* Analyze interface */
863   if (computetopography) {
864     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
865   }
866 
867   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
868   new_nearnullspace_provided = PETSC_FALSE;
869   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
870   if (pcbddc->onearnullspace) { /* already used nearnullspace */
871     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
872       new_nearnullspace_provided = PETSC_TRUE;
873     } else {
874       /* determine if the two nullspaces are different (should be lightweight) */
875       if (nearnullspace != pcbddc->onearnullspace) {
876         new_nearnullspace_provided = PETSC_TRUE;
877       } else { /* maybe the user has changed the content of the nearnullspace */
878         ierr = MatNullSpaceGetVecs(nearnullspace,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
879         ierr = MatNullSpaceGetVecs(pcbddc->onearnullspace,&onnsp_has_cnst,&onnsp_size,&onearnullvecs);CHKERRQ(ierr);
880         if ( (nnsp_has_cnst != onnsp_has_cnst) || (nnsp_size != onnsp_size) ) {
881           new_nearnullspace_provided = PETSC_TRUE;
882         } else { /* nullspaces have the same size, so check vectors or their ObjectStateId */
883           PetscInt i;
884           for (i=0;i<nnsp_size;i++) {
885             ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
886             if (nearnullvecs[i] != onearnullvecs[i] || pcbddc->onearnullvecs_state[i] != state) {
887               new_nearnullspace_provided = PETSC_TRUE;
888               break;
889             }
890           }
891         }
892       }
893     }
894   } else {
895     if (!nearnullspace) { /* both nearnullspaces are null */
896       new_nearnullspace_provided = PETSC_FALSE;
897     } else { /* nearnullspace attached later */
898       new_nearnullspace_provided = PETSC_TRUE;
899     }
900   }
901 
902   /* Setup constraints and related work vectors */
903   /* reset primal space flags */
904   pcbddc->new_primal_space = PETSC_FALSE;
905   pcbddc->new_primal_space_local = PETSC_FALSE;
906   if (computetopography || new_nearnullspace_provided) {
907     /* It also sets the primal space flags */
908     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
909     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
910     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
911   }
912 
913   if (computesolvers || pcbddc->new_primal_space) {
914     /* reset data */
915     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
916     /* Create coarse and local stuffs */
917     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
918     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
919   }
920   if (pcbddc->dbg_flag && pcbddc->current_level) {
921     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
922   }
923   PetscFunctionReturn(0);
924 }
925 
926 /* -------------------------------------------------------------------------- */
927 /*
928    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
929 
930    Input Parameters:
931 .  pc - the preconditioner context
932 .  r - input vector (global)
933 
934    Output Parameter:
935 .  z - output vector (global)
936 
937    Application Interface Routine: PCApply()
938  */
939 #undef __FUNCT__
940 #define __FUNCT__ "PCApply_BDDC"
941 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
942 {
943   PC_IS             *pcis = (PC_IS*)(pc->data);
944   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
945   PetscErrorCode    ierr;
946   const PetscScalar one = 1.0;
947   const PetscScalar m_one = -1.0;
948   const PetscScalar zero = 0.0;
949 
950 /* This code is similar to that provided in nn.c for PCNN
951    NN interface preconditioner changed to BDDC
952    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
953 
954   PetscFunctionBegin;
955   if (!pcbddc->use_exact_dirichlet_trick) {
956     /* First Dirichlet solve */
957     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
958     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
959     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
960     /*
961       Assembling right hand side for BDDC operator
962       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
963       - pcis->vec1_B the interface part of the global vector z
964     */
965     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
966     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
967     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
968     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
969     ierr = VecCopy(r,z);CHKERRQ(ierr);
970     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
971     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
972     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
973   } else {
974     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
975     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
976     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
977   }
978 
979   /* Apply interface preconditioner
980      input/output vecs: pcis->vec1_B and pcis->vec1_D */
981   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
982 
983   /* Apply transpose of partition of unity operator */
984   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
985 
986   /* Second Dirichlet solve and assembling of output */
987   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
988   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
989   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
990   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
991   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
992   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
993   if (pcbddc->switch_static) { ierr = VecAXPY (pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
994   ierr = VecAXPY (pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
995   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
996   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
997   PetscFunctionReturn(0);
998 }
999 /* -------------------------------------------------------------------------- */
1000 
1001 #undef __FUNCT__
1002 #define __FUNCT__ "PCDestroy_BDDC"
1003 PetscErrorCode PCDestroy_BDDC(PC pc)
1004 {
1005   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1006   PetscErrorCode ierr;
1007 
1008   PetscFunctionBegin;
1009   /* free data created by PCIS */
1010   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1011   /* free BDDC custom data  */
1012   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1013   /* destroy objects related to topography */
1014   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1015   /* free allocated graph structure */
1016   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1017   /* free data for scaling operator */
1018   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1019   /* free solvers stuff */
1020   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1021   /* free global vectors needed in presolve */
1022   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
1023   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1024   /* remove functions */
1025   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1026   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
1027   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1028   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
1029   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1030   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1031   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1032   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1033   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1034   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1035   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1036   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1037   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1038   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1039   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1040   /* Free the private data structure */
1041   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1042   PetscFunctionReturn(0);
1043 }
1044 /* -------------------------------------------------------------------------- */
1045 
1046 #undef __FUNCT__
1047 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
1048 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1049 {
1050   FETIDPMat_ctx  mat_ctx;
1051   PC_IS*         pcis;
1052   PC_BDDC*       pcbddc;
1053   PetscErrorCode ierr;
1054 
1055   PetscFunctionBegin;
1056   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1057   pcis = (PC_IS*)mat_ctx->pc->data;
1058   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1059 
1060   /* change of basis for physical rhs if needed
1061      It also changes the rhs in case of dirichlet boundaries */
1062   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
1063   /* store vectors for computation of fetidp final solution */
1064   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1065   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1066   /* scale rhs since it should be unassembled */
1067   /* TODO use counter scaling? (also below) */
1068   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1069   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1070   /* Apply partition of unity */
1071   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1072   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1073   if (!pcbddc->switch_static) {
1074     /* compute partially subassembled Schur complement right-hand side */
1075     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1076     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
1077     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1078     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
1079     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1080     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1081     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1082     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1083     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1084     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1085   }
1086   /* BDDC rhs */
1087   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
1088   if (pcbddc->switch_static) {
1089     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1090   }
1091   /* apply BDDC */
1092   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
1093   /* Application of B_delta and assembling of rhs for fetidp fluxes */
1094   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
1095   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
1096   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1097   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1098   /* restore original rhs */
1099   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
1100   PetscFunctionReturn(0);
1101 }
1102 
1103 #undef __FUNCT__
1104 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
1105 /*@
1106  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
1107 
1108    Collective
1109 
1110    Input Parameters:
1111 +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
1112 .  standard_rhs - the right-hand side for your linear system
1113 
1114    Output Parameters:
1115 -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
1116 
1117    Level: developer
1118 
1119    Notes:
1120 
1121 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1122 @*/
1123 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1124 {
1125   FETIDPMat_ctx  mat_ctx;
1126   PetscErrorCode ierr;
1127 
1128   PetscFunctionBegin;
1129   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1130   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
1131   PetscFunctionReturn(0);
1132 }
1133 /* -------------------------------------------------------------------------- */
1134 
1135 #undef __FUNCT__
1136 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
1137 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1138 {
1139   FETIDPMat_ctx  mat_ctx;
1140   PC_IS*         pcis;
1141   PC_BDDC*       pcbddc;
1142   PetscErrorCode ierr;
1143 
1144   PetscFunctionBegin;
1145   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1146   pcis = (PC_IS*)mat_ctx->pc->data;
1147   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1148 
1149   /* apply B_delta^T */
1150   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1151   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1152   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
1153   /* compute rhs for BDDC application */
1154   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1155   if (pcbddc->switch_static) {
1156     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1157   }
1158   /* apply BDDC */
1159   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
1160   /* put values into standard global vector */
1161   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1162   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1163   if (!pcbddc->switch_static) {
1164     /* compute values into the interior if solved for the partially subassembled Schur complement */
1165     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
1166     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
1167     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1168   }
1169   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1170   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1171   /* final change of basis if needed
1172      Is also sums the dirichlet part removed during RHS assembling */
1173   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
1174   PetscFunctionReturn(0);
1175 }
1176 
1177 #undef __FUNCT__
1178 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
1179 /*@
1180  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
1181 
1182    Collective
1183 
1184    Input Parameters:
1185 +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
1186 .  fetidp_flux_sol - the solution of the FETIDP linear system
1187 
1188    Output Parameters:
1189 -  standard_sol      - the solution defined on the physical domain
1190 
1191    Level: developer
1192 
1193    Notes:
1194 
1195 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1196 @*/
1197 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1198 {
1199   FETIDPMat_ctx  mat_ctx;
1200   PetscErrorCode ierr;
1201 
1202   PetscFunctionBegin;
1203   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1204   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
1205   PetscFunctionReturn(0);
1206 }
1207 /* -------------------------------------------------------------------------- */
1208 
1209 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1210 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1211 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1212 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1213 
1214 #undef __FUNCT__
1215 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
1216 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1217 {
1218 
1219   FETIDPMat_ctx  fetidpmat_ctx;
1220   Mat            newmat;
1221   FETIDPPC_ctx   fetidppc_ctx;
1222   PC             newpc;
1223   MPI_Comm       comm;
1224   PetscErrorCode ierr;
1225 
1226   PetscFunctionBegin;
1227   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1228   /* FETIDP linear matrix */
1229   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
1230   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
1231   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
1232   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1233   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
1234   ierr = MatSetUp(newmat);CHKERRQ(ierr);
1235   /* FETIDP preconditioner */
1236   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
1237   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
1238   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
1239   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
1240   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
1241   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1242   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
1243   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
1244   ierr = PCSetUp(newpc);CHKERRQ(ierr);
1245   /* return pointers for objects created */
1246   *fetidp_mat=newmat;
1247   *fetidp_pc=newpc;
1248   PetscFunctionReturn(0);
1249 }
1250 
1251 #undef __FUNCT__
1252 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
1253 /*@
1254  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
1255 
1256    Collective
1257 
1258    Input Parameters:
1259 +  pc - the BDDC preconditioning context already setup
1260 
1261    Output Parameters:
1262 .  fetidp_mat - shell FETIDP matrix object
1263 .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
1264 
1265    Options Database Keys:
1266 -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
1267 
1268    Level: developer
1269 
1270    Notes:
1271      Currently the only operation provided for FETIDP matrix is MatMult
1272 
1273 .seealso: PCBDDC
1274 @*/
1275 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1276 {
1277   PetscErrorCode ierr;
1278 
1279   PetscFunctionBegin;
1280   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1281   if (pc->setupcalled) {
1282     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1283   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1284   PetscFunctionReturn(0);
1285 }
1286 /* -------------------------------------------------------------------------- */
1287 /*MC
1288    PCBDDC - Balancing Domain Decomposition by Constraints.
1289 
1290    An implementation of the BDDC preconditioner based on
1291 
1292 .vb
1293    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
1294    [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
1295    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
1296 .ve
1297 
1298    The matrix to be preconditioned (Pmat) must be of type MATIS.
1299 
1300    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
1301 
1302    It also works with unsymmetric and indefinite problems.
1303 
1304    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.
1305 
1306    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
1307 
1308    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
1309 
1310    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace.
1311 
1312    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.
1313 
1314    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
1315 
1316    Options Database Keys:
1317 
1318 .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
1319 .    -pc_bddc_use_edges <1> - use or not edges in primal space
1320 .    -pc_bddc_use_faces <0> - use or not faces in primal space
1321 .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
1322 .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
1323 .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
1324 .    -pc_bddc_levels <0> - maximum number of levels for multilevel
1325 .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
1326 -    -pc_bddc_check_level <0> - set verbosity level of debugging output
1327 
1328    Options for Dirichlet, Neumann or coarse solver can be set with
1329 .vb
1330       -pc_bddc_dirichlet_
1331       -pc_bddc_neumann_
1332       -pc_bddc_coarse_
1333 .ve
1334    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
1335 
1336    When using a multilevel approach, solvers' options at the N-th level can be specified as
1337 .vb
1338       -pc_bddc_dirichlet_N_
1339       -pc_bddc_neumann_N_
1340       -pc_bddc_coarse_N_
1341 .ve
1342    Note that level number ranges from the finest 0 to the coarsest N
1343 
1344    Level: intermediate
1345 
1346    Developer notes:
1347      Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose
1348 
1349      New deluxe scaling operator will be available soon.
1350 
1351    Contributed by Stefano Zampini
1352 
1353 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1354 M*/
1355 
1356 #undef __FUNCT__
1357 #define __FUNCT__ "PCCreate_BDDC"
1358 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1359 {
1360   PetscErrorCode      ierr;
1361   PC_BDDC             *pcbddc;
1362 
1363   PetscFunctionBegin;
1364   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1365   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1366   pc->data  = (void*)pcbddc;
1367 
1368   /* create PCIS data structure */
1369   ierr = PCISCreate(pc);CHKERRQ(ierr);
1370 
1371   /* BDDC customization */
1372   pcbddc->use_vertices        = PETSC_TRUE;
1373   pcbddc->use_edges           = PETSC_TRUE;
1374   pcbddc->use_faces           = PETSC_FALSE;
1375   pcbddc->use_change_of_basis = PETSC_FALSE;
1376   pcbddc->use_change_on_faces = PETSC_FALSE;
1377   pcbddc->switch_static       = PETSC_FALSE;
1378   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
1379   pcbddc->dbg_flag            = 0;
1380 
1381   pcbddc->local_primal_size          = 0;
1382   pcbddc->primal_indices_local_idxs  = 0;
1383   pcbddc->recompute_topography       = PETSC_FALSE;
1384   pcbddc->coarse_size                = 0;
1385   pcbddc->new_primal_space           = PETSC_FALSE;
1386   pcbddc->new_primal_space_local     = PETSC_FALSE;
1387   pcbddc->global_primal_indices      = 0;
1388   pcbddc->onearnullspace             = 0;
1389   pcbddc->onearnullvecs_state        = 0;
1390   pcbddc->user_primal_vertices       = 0;
1391   pcbddc->NullSpace                  = 0;
1392   pcbddc->temp_solution              = 0;
1393   pcbddc->original_rhs               = 0;
1394   pcbddc->local_mat                  = 0;
1395   pcbddc->ChangeOfBasisMatrix        = 0;
1396   pcbddc->coarse_vec                 = 0;
1397   pcbddc->coarse_rhs                 = 0;
1398   pcbddc->coarse_ksp                 = 0;
1399   pcbddc->coarse_phi_B               = 0;
1400   pcbddc->coarse_phi_D               = 0;
1401   pcbddc->coarse_psi_B               = 0;
1402   pcbddc->coarse_psi_D               = 0;
1403   pcbddc->vec1_P                     = 0;
1404   pcbddc->vec1_R                     = 0;
1405   pcbddc->vec2_R                     = 0;
1406   pcbddc->local_auxmat1              = 0;
1407   pcbddc->local_auxmat2              = 0;
1408   pcbddc->R_to_B                     = 0;
1409   pcbddc->R_to_D                     = 0;
1410   pcbddc->ksp_D                      = 0;
1411   pcbddc->ksp_R                      = 0;
1412   pcbddc->NeumannBoundaries          = 0;
1413   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
1414   pcbddc->n_ISForDofs                = 0;
1415   pcbddc->ISForDofs                  = 0;
1416   pcbddc->ConstraintMatrix           = 0;
1417   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
1418   pcbddc->coarse_loc_to_glob         = 0;
1419   pcbddc->coarsening_ratio           = 8;
1420   pcbddc->current_level              = 0;
1421   pcbddc->max_levels                 = 0;
1422 
1423   /* create local graph structure */
1424   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1425 
1426   /* scaling */
1427   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1428   pcbddc->work_scaling               = 0;
1429 
1430   /* function pointers */
1431   pc->ops->apply               = PCApply_BDDC;
1432   pc->ops->applytranspose      = 0;
1433   pc->ops->setup               = PCSetUp_BDDC;
1434   pc->ops->destroy             = PCDestroy_BDDC;
1435   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1436   pc->ops->view                = 0;
1437   pc->ops->applyrichardson     = 0;
1438   pc->ops->applysymmetricleft  = 0;
1439   pc->ops->applysymmetricright = 0;
1440   pc->ops->presolve            = PCPreSolve_BDDC;
1441   pc->ops->postsolve           = PCPostSolve_BDDC;
1442 
1443   /* composing function */
1444   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1445   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
1446   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1447   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
1448   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1449   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1450   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1451   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1452   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1453   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1454   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1455   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1456   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1457   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1458   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1459   PetscFunctionReturn(0);
1460 }
1461 
1462