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