xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 906d46d482912ff65fe95caf6add84f15a3f7e20)
1 /* TODOLIST
2 
3    ConstraintsSetup
4    - tolerances for constraints as an option (take care of single precision!)
5 
6    Solvers
7    - Add support for cholesky for coarse solver (similar to local solvers)
8    - Propagate ksp prefixes for solvers to mat objects?
9    - Propagate nearnullspace info among levels
10 
11    User interface
12    - ** DofSplitting and DM attached to pc?
13 
14    Debugging output
15    - * Better management of verbosity levels of debugging output
16 
17    Build
18    - make runexe59
19 
20    Extra
21    - ** GetRid of PCBDDCApplySchur, use MatSchur instead
22    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
23    - add support for computing h,H and related using coordinates?
24    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
25    - Better management in PCIS code
26    - BDDC with MG framework?
27 
28    FETIDP
29    - Move FETIDP code to its own classes
30 
31    MATIS related operations contained in BDDC code
32    - Provide general case for subassembling
33    - *** Preallocation routines in MatISGetMPIAXAIJ
34 
35 */
36 
37 /* ----------------------------------------------------------------------------------------------------------------------------------------------
38    Implementation of BDDC preconditioner based on:
39    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
40    ---------------------------------------------------------------------------------------------------------------------------------------------- */
41 
42 #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
43 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
44 #include <petscblaslapack.h>
45 
46 /* -------------------------------------------------------------------------- */
47 #undef __FUNCT__
48 #define __FUNCT__ "PCSetFromOptions_BDDC"
49 PetscErrorCode PCSetFromOptions_BDDC(PC pc)
50 {
51   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
52   PetscErrorCode ierr;
53 
54   PetscFunctionBegin;
55   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
56   /* Verbose debugging */
57   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
58   /* Primal space cumstomization */
59   ierr = PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);CHKERRQ(ierr);
60   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
61   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
62   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
63   /* Change of basis */
64   ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not internal change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr);
65   ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not internal change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr);
66   if (!pcbddc->use_change_of_basis) {
67     pcbddc->use_change_on_faces = PETSC_FALSE;
68   }
69   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
70   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);
71   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
72   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
73   ierr = PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);CHKERRQ(ierr);
74   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
75   ierr = PetscOptionsTail();CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 /* -------------------------------------------------------------------------- */
79 #undef __FUNCT__
80 #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
81 static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
82 {
83   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
84   PetscErrorCode ierr;
85 
86   PetscFunctionBegin;
87   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
88   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
89   pcbddc->user_ChangeOfBasisMatrix = change;
90   PetscFunctionReturn(0);
91 }
92 #undef __FUNCT__
93 #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
94 /*@
95  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
96 
97    Collective on PC
98 
99    Input Parameters:
100 +  pc - the preconditioning context
101 -  change - the change of basis matrix
102 
103    Level: intermediate
104 
105    Notes:
106 
107 .seealso: PCBDDC
108 @*/
109 PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
110 {
111   PetscErrorCode ierr;
112 
113   PetscFunctionBegin;
114   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
115   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
116   PetscCheckSameComm(pc,1,change,2);
117   if (pc->mat) {
118     PetscInt rows_c,cols_c,rows,cols;
119     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
120     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
121     if (rows_c != rows) {
122       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
123     }
124     if (cols_c != cols) {
125       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
126     }
127     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
128     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
129     if (rows_c != rows) {
130       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
131     }
132     if (cols_c != cols) {
133       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
134     }
135   }
136   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
137   PetscFunctionReturn(0);
138 }
139 /* -------------------------------------------------------------------------- */
140 #undef __FUNCT__
141 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
142 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
143 {
144   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
145   PetscErrorCode ierr;
146 
147   PetscFunctionBegin;
148   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
149   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
150   pcbddc->user_primal_vertices = PrimalVertices;
151   PetscFunctionReturn(0);
152 }
153 #undef __FUNCT__
154 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
155 /*@
156  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
157 
158    Not collective
159 
160    Input Parameters:
161 +  pc - the preconditioning context
162 -  PrimalVertices - index set of primal vertices in local numbering
163 
164    Level: intermediate
165 
166    Notes:
167 
168 .seealso: PCBDDC
169 @*/
170 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
171 {
172   PetscErrorCode ierr;
173 
174   PetscFunctionBegin;
175   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
176   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
177   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
178   PetscFunctionReturn(0);
179 }
180 /* -------------------------------------------------------------------------- */
181 #undef __FUNCT__
182 #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
183 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
184 {
185   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
186 
187   PetscFunctionBegin;
188   pcbddc->coarsening_ratio = k;
189   PetscFunctionReturn(0);
190 }
191 
192 #undef __FUNCT__
193 #define __FUNCT__ "PCBDDCSetCoarseningRatio"
194 /*@
195  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
196 
197    Logically collective on PC
198 
199    Input Parameters:
200 +  pc - the preconditioning context
201 -  k - coarsening ratio (H/h at the coarser level)
202 
203    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
204 
205    Level: intermediate
206 
207    Notes:
208 
209 .seealso: PCBDDC
210 @*/
211 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
212 {
213   PetscErrorCode ierr;
214 
215   PetscFunctionBegin;
216   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
217   PetscValidLogicalCollectiveInt(pc,k,2);
218   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
219   PetscFunctionReturn(0);
220 }
221 
222 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
223 #undef __FUNCT__
224 #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
225 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
226 {
227   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
228 
229   PetscFunctionBegin;
230   pcbddc->use_exact_dirichlet_trick = flg;
231   PetscFunctionReturn(0);
232 }
233 
234 #undef __FUNCT__
235 #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
236 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
237 {
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
242   PetscValidLogicalCollectiveBool(pc,flg,2);
243   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
244   PetscFunctionReturn(0);
245 }
246 
247 #undef __FUNCT__
248 #define __FUNCT__ "PCBDDCSetLevel_BDDC"
249 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
250 {
251   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
252 
253   PetscFunctionBegin;
254   pcbddc->current_level = level;
255   PetscFunctionReturn(0);
256 }
257 
258 #undef __FUNCT__
259 #define __FUNCT__ "PCBDDCSetLevel"
260 PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
261 {
262   PetscErrorCode ierr;
263 
264   PetscFunctionBegin;
265   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
266   PetscValidLogicalCollectiveInt(pc,level,2);
267   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
268   PetscFunctionReturn(0);
269 }
270 
271 #undef __FUNCT__
272 #define __FUNCT__ "PCBDDCSetLevels_BDDC"
273 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
274 {
275   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
276 
277   PetscFunctionBegin;
278   pcbddc->max_levels = levels;
279   PetscFunctionReturn(0);
280 }
281 
282 #undef __FUNCT__
283 #define __FUNCT__ "PCBDDCSetLevels"
284 /*@
285  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
286 
287    Logically collective on PC
288 
289    Input Parameters:
290 +  pc - the preconditioning context
291 -  levels - the maximum number of levels (max 9)
292 
293    Default value is 0, i.e. traditional one-level BDDC
294 
295    Level: intermediate
296 
297    Notes:
298 
299 .seealso: PCBDDC
300 @*/
301 PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
302 {
303   PetscErrorCode ierr;
304 
305   PetscFunctionBegin;
306   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
307   PetscValidLogicalCollectiveInt(pc,levels,2);
308   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
309   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
310   PetscFunctionReturn(0);
311 }
312 /* -------------------------------------------------------------------------- */
313 
314 #undef __FUNCT__
315 #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
316 static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
317 {
318   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
319   PetscErrorCode ierr;
320 
321   PetscFunctionBegin;
322   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
323   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
324   pcbddc->NullSpace = NullSpace;
325   PetscFunctionReturn(0);
326 }
327 
328 #undef __FUNCT__
329 #define __FUNCT__ "PCBDDCSetNullSpace"
330 /*@
331  PCBDDCSetNullSpace - Set nullspace for BDDC operator
332 
333    Logically collective on PC and MatNullSpace
334 
335    Input Parameters:
336 +  pc - the preconditioning context
337 -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
338 
339    Level: intermediate
340 
341    Notes:
342 
343 .seealso: PCBDDC
344 @*/
345 PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
346 {
347   PetscErrorCode ierr;
348 
349   PetscFunctionBegin;
350   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
351   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
352   PetscCheckSameComm(pc,1,NullSpace,2);
353   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
354   PetscFunctionReturn(0);
355 }
356 /* -------------------------------------------------------------------------- */
357 
358 #undef __FUNCT__
359 #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
360 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
361 {
362   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
363   PetscErrorCode ierr;
364 
365   PetscFunctionBegin;
366   /* last user setting takes precendence -> destroy any other customization */
367   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
368   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
369   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
370   pcbddc->DirichletBoundaries = DirichletBoundaries;
371   pcbddc->recompute_topography = PETSC_TRUE;
372   PetscFunctionReturn(0);
373 }
374 
375 #undef __FUNCT__
376 #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
377 /*@
378  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
379 
380    Collective
381 
382    Input Parameters:
383 +  pc - the preconditioning context
384 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
385 
386    Level: intermediate
387 
388    Notes: Any process can list any global node
389 
390 .seealso: PCBDDC
391 @*/
392 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
393 {
394   PetscErrorCode ierr;
395 
396   PetscFunctionBegin;
397   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
398   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
399   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
400   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
401   PetscFunctionReturn(0);
402 }
403 /* -------------------------------------------------------------------------- */
404 
405 #undef __FUNCT__
406 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
407 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
408 {
409   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
410   PetscErrorCode ierr;
411 
412   PetscFunctionBegin;
413   /* last user setting takes precendence -> destroy any other customization */
414   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
415   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
416   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
417   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
418   pcbddc->recompute_topography = PETSC_TRUE;
419   PetscFunctionReturn(0);
420 }
421 
422 #undef __FUNCT__
423 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
424 /*@
425  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
426 
427    Collective
428 
429    Input Parameters:
430 +  pc - the preconditioning context
431 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
432 
433    Level: intermediate
434 
435    Notes:
436 
437 .seealso: PCBDDC
438 @*/
439 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
440 {
441   PetscErrorCode ierr;
442 
443   PetscFunctionBegin;
444   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
445   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
446   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
447   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
448   PetscFunctionReturn(0);
449 }
450 /* -------------------------------------------------------------------------- */
451 
452 #undef __FUNCT__
453 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
454 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
455 {
456   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
457   PetscErrorCode ierr;
458 
459   PetscFunctionBegin;
460   /* last user setting takes precendence -> destroy any other customization */
461   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
462   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
463   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
464   pcbddc->NeumannBoundaries = NeumannBoundaries;
465   pcbddc->recompute_topography = PETSC_TRUE;
466   PetscFunctionReturn(0);
467 }
468 
469 #undef __FUNCT__
470 #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
471 /*@
472  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
473 
474    Collective
475 
476    Input Parameters:
477 +  pc - the preconditioning context
478 -  NeumannBoundaries - parallel IS defining the Neumann boundaries
479 
480    Level: intermediate
481 
482    Notes: Any process can list any global node
483 
484 .seealso: PCBDDC
485 @*/
486 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
487 {
488   PetscErrorCode ierr;
489 
490   PetscFunctionBegin;
491   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
492   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
493   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
494   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
495   PetscFunctionReturn(0);
496 }
497 /* -------------------------------------------------------------------------- */
498 
499 #undef __FUNCT__
500 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
501 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
502 {
503   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
504   PetscErrorCode ierr;
505 
506   PetscFunctionBegin;
507   /* last user setting takes precendence -> destroy any other customization */
508   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
509   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
510   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
511   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
512   pcbddc->recompute_topography = PETSC_TRUE;
513   PetscFunctionReturn(0);
514 }
515 
516 #undef __FUNCT__
517 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
518 /*@
519  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
520 
521    Collective
522 
523    Input Parameters:
524 +  pc - the preconditioning context
525 -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
526 
527    Level: intermediate
528 
529    Notes:
530 
531 .seealso: PCBDDC
532 @*/
533 PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
534 {
535   PetscErrorCode ierr;
536 
537   PetscFunctionBegin;
538   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
539   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
540   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
541   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
542   PetscFunctionReturn(0);
543 }
544 /* -------------------------------------------------------------------------- */
545 
546 #undef __FUNCT__
547 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
548 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
549 {
550   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
551 
552   PetscFunctionBegin;
553   *DirichletBoundaries = pcbddc->DirichletBoundaries;
554   PetscFunctionReturn(0);
555 }
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
559 /*@
560  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
561 
562    Collective
563 
564    Input Parameters:
565 .  pc - the preconditioning context
566 
567    Output Parameters:
568 .  DirichletBoundaries - index set defining the Dirichlet boundaries
569 
570    Level: intermediate
571 
572    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
573 
574 .seealso: PCBDDC
575 @*/
576 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
577 {
578   PetscErrorCode ierr;
579 
580   PetscFunctionBegin;
581   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
582   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
583   PetscFunctionReturn(0);
584 }
585 /* -------------------------------------------------------------------------- */
586 
587 #undef __FUNCT__
588 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
589 static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
590 {
591   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
592 
593   PetscFunctionBegin;
594   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
595   PetscFunctionReturn(0);
596 }
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
600 /*@
601  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
602 
603    Collective
604 
605    Input Parameters:
606 .  pc - the preconditioning context
607 
608    Output Parameters:
609 .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
610 
611    Level: intermediate
612 
613    Notes:
614 
615 .seealso: PCBDDC
616 @*/
617 PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
618 {
619   PetscErrorCode ierr;
620 
621   PetscFunctionBegin;
622   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
623   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
624   PetscFunctionReturn(0);
625 }
626 /* -------------------------------------------------------------------------- */
627 
628 #undef __FUNCT__
629 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
630 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
631 {
632   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
633 
634   PetscFunctionBegin;
635   *NeumannBoundaries = pcbddc->NeumannBoundaries;
636   PetscFunctionReturn(0);
637 }
638 
639 #undef __FUNCT__
640 #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
641 /*@
642  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
643 
644    Collective
645 
646    Input Parameters:
647 .  pc - the preconditioning context
648 
649    Output Parameters:
650 .  NeumannBoundaries - index set defining the Neumann boundaries
651 
652    Level: intermediate
653 
654    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
655 
656 .seealso: PCBDDC
657 @*/
658 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
659 {
660   PetscErrorCode ierr;
661 
662   PetscFunctionBegin;
663   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
664   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
665   PetscFunctionReturn(0);
666 }
667 /* -------------------------------------------------------------------------- */
668 
669 #undef __FUNCT__
670 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
671 static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
672 {
673   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
674 
675   PetscFunctionBegin;
676   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
677   PetscFunctionReturn(0);
678 }
679 
680 #undef __FUNCT__
681 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
682 /*@
683  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
684 
685    Collective
686 
687    Input Parameters:
688 .  pc - the preconditioning context
689 
690    Output Parameters:
691 .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
692 
693    Level: intermediate
694 
695    Notes:
696 
697 .seealso: PCBDDC
698 @*/
699 PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
700 {
701   PetscErrorCode ierr;
702 
703   PetscFunctionBegin;
704   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
705   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
706   PetscFunctionReturn(0);
707 }
708 /* -------------------------------------------------------------------------- */
709 
710 #undef __FUNCT__
711 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
712 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
713 {
714   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
715   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
716   PetscErrorCode ierr;
717 
718   PetscFunctionBegin;
719   /* free old CSR */
720   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
721   /* TODO: PCBDDCGraphSetAdjacency */
722   /* get CSR into graph structure */
723   if (copymode == PETSC_COPY_VALUES) {
724     ierr = PetscMalloc1((nvtxs+1),&mat_graph->xadj);CHKERRQ(ierr);
725     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
726     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
727     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
728   } else if (copymode == PETSC_OWN_POINTER) {
729     mat_graph->xadj = (PetscInt*)xadj;
730     mat_graph->adjncy = (PetscInt*)adjncy;
731   }
732   mat_graph->nvtxs_csr = nvtxs;
733   PetscFunctionReturn(0);
734 }
735 
736 #undef __FUNCT__
737 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
738 /*@
739  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
740 
741    Not collective
742 
743    Input Parameters:
744 +  pc - the preconditioning context
745 .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
746 .  xadj, adjncy - the CSR graph
747 -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
748 
749    Level: intermediate
750 
751    Notes:
752 
753 .seealso: PCBDDC,PetscCopyMode
754 @*/
755 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
756 {
757   void (*f)(void) = 0;
758   PetscErrorCode ierr;
759 
760   PetscFunctionBegin;
761   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
762   PetscValidIntPointer(xadj,3);
763   PetscValidIntPointer(xadj,4);
764   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
765     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
766   }
767   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
768   /* free arrays if PCBDDC is not the PC type */
769   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
770   if (!f && copymode == PETSC_OWN_POINTER) {
771     ierr = PetscFree(xadj);CHKERRQ(ierr);
772     ierr = PetscFree(adjncy);CHKERRQ(ierr);
773   }
774   PetscFunctionReturn(0);
775 }
776 /* -------------------------------------------------------------------------- */
777 
778 #undef __FUNCT__
779 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
780 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
781 {
782   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
783   PetscInt i;
784   PetscErrorCode ierr;
785 
786   PetscFunctionBegin;
787   /* Destroy ISes if they were already set */
788   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
789     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
790   }
791   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
792   /* last user setting takes precendence -> destroy any other customization */
793   for (i=0;i<pcbddc->n_ISForDofs;i++) {
794     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
795   }
796   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
797   pcbddc->n_ISForDofs = 0;
798   /* allocate space then set */
799   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
800   for (i=0;i<n_is;i++) {
801     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
802     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
803   }
804   pcbddc->n_ISForDofsLocal=n_is;
805   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
806   pcbddc->recompute_topography = PETSC_TRUE;
807   PetscFunctionReturn(0);
808 }
809 
810 #undef __FUNCT__
811 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
812 /*@
813  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
814 
815    Collective
816 
817    Input Parameters:
818 +  pc - the preconditioning context
819 -  n_is - number of index sets defining the fields
820 .  ISForDofs - array of IS describing the fields in local ordering
821 
822    Level: intermediate
823 
824    Notes: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to a different field.
825 
826 .seealso: PCBDDC
827 @*/
828 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
829 {
830   PetscInt       i;
831   PetscErrorCode ierr;
832 
833   PetscFunctionBegin;
834   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
835   PetscValidLogicalCollectiveInt(pc,n_is,2);
836   for (i=0;i<n_is;i++) {
837     PetscCheckSameComm(pc,1,ISForDofs[i],3);
838     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
839   }
840   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
841   PetscFunctionReturn(0);
842 }
843 /* -------------------------------------------------------------------------- */
844 
845 #undef __FUNCT__
846 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
847 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
848 {
849   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
850   PetscInt i;
851   PetscErrorCode ierr;
852 
853   PetscFunctionBegin;
854   /* Destroy ISes if they were already set */
855   for (i=0;i<pcbddc->n_ISForDofs;i++) {
856     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
857   }
858   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
859   /* last user setting takes precendence -> destroy any other customization */
860   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
861     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
862   }
863   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
864   pcbddc->n_ISForDofsLocal = 0;
865   /* allocate space then set */
866   ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
867   for (i=0;i<n_is;i++) {
868     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
869     pcbddc->ISForDofs[i]=ISForDofs[i];
870   }
871   pcbddc->n_ISForDofs=n_is;
872   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
873   pcbddc->recompute_topography = PETSC_TRUE;
874   PetscFunctionReturn(0);
875 }
876 
877 #undef __FUNCT__
878 #define __FUNCT__ "PCBDDCSetDofsSplitting"
879 /*@
880  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
881 
882    Collective
883 
884    Input Parameters:
885 +  pc - the preconditioning context
886 -  n_is - number of index sets defining the fields
887 .  ISForDofs - array of IS describing the fields in global ordering
888 
889    Level: intermediate
890 
891    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.
892 
893 .seealso: PCBDDC
894 @*/
895 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
896 {
897   PetscInt       i;
898   PetscErrorCode ierr;
899 
900   PetscFunctionBegin;
901   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
902   PetscValidLogicalCollectiveInt(pc,n_is,2);
903   for (i=0;i<n_is;i++) {
904     PetscCheckSameComm(pc,1,ISForDofs[i],3);
905     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
906   }
907   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
908   PetscFunctionReturn(0);
909 }
910 
911 /* -------------------------------------------------------------------------- */
912 #undef __FUNCT__
913 #define __FUNCT__ "PCPreSolve_BDDC"
914 /* -------------------------------------------------------------------------- */
915 /*
916    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
917                      guess if a transformation of basis approach has been selected.
918 
919    Input Parameter:
920 +  pc - the preconditioner contex
921 
922    Application Interface Routine: PCPreSolve()
923 
924    Notes:
925    The interface routine PCPreSolve() is not usually called directly by
926    the user, but instead is called by KSPSolve().
927 */
928 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
929 {
930   PetscErrorCode ierr;
931   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
932   PC_IS          *pcis = (PC_IS*)(pc->data);
933   IS             dirIS;
934   Vec            used_vec;
935 
936   PetscFunctionBegin;
937   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
938   if (ksp) {
939     PetscBool iscg;
940     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
941     if (!iscg) {
942       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
943     }
944   }
945   /* Creates parallel work vectors used in presolve */
946   if (!pcbddc->original_rhs) {
947     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
948   }
949   if (!pcbddc->temp_solution) {
950     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
951   }
952   if (x) {
953     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
954     used_vec = x;
955   } else {
956     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
957     used_vec = pcbddc->temp_solution;
958     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
959   }
960   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
961   if (ksp) {
962     PetscBool guess_nonzero;
963     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
964     if (!guess_nonzero) {
965       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
966     }
967   }
968 
969   /* store the original rhs */
970   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
971 
972   /* Take into account zeroed rows -> change rhs and store solution removed */
973   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
974   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
975   if (rhs && dirIS) {
976     Mat_IS      *matis = (Mat_IS*)pc->pmat->data;
977     PetscInt    dirsize,i,*is_indices;
978     PetscScalar *array_x,*array_diagonal;
979 
980     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
981     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
982     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
983     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
984     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
985     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
986     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
987     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
988     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
989     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
990     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
991     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
992     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
993     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
994     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
995     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
996 
997     /* remove the computed solution from the rhs */
998     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
999     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
1000     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
1001   }
1002 
1003   /* store partially computed solution and set initial guess */
1004   if (x) {
1005     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
1006     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
1007     if (pcbddc->use_exact_dirichlet_trick) {
1008       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1009       ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1010       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1011       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1012       ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1013       if (ksp) {
1014         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1015       }
1016     }
1017   }
1018 
1019   if (pcbddc->ChangeOfBasisMatrix) {
1020     PCBDDCChange_ctx change_ctx;
1021 
1022     /* get change ctx */
1023     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1024 
1025     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1026     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1027     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1028     change_ctx->original_mat = pc->mat;
1029 
1030     /* change iteration matrix */
1031     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1032     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1033     pc->mat = pcbddc->new_global_mat;
1034 
1035     /* change rhs */
1036     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1037     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
1038   }
1039 
1040   /* remove nullspace if present */
1041   if (ksp && pcbddc->NullSpace) {
1042     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
1043     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1044   }
1045   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1046   PetscFunctionReturn(0);
1047 }
1048 
1049 /* -------------------------------------------------------------------------- */
1050 #undef __FUNCT__
1051 #define __FUNCT__ "PCPostSolve_BDDC"
1052 /* -------------------------------------------------------------------------- */
1053 /*
1054    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1055                      approach has been selected. Also, restores rhs to its original state.
1056 
1057    Input Parameter:
1058 +  pc - the preconditioner contex
1059 
1060    Application Interface Routine: PCPostSolve()
1061 
1062    Notes:
1063    The interface routine PCPostSolve() is not usually called directly by
1064    the user, but instead is called by KSPSolve().
1065 */
1066 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1067 {
1068   PetscErrorCode ierr;
1069   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1070 
1071   PetscFunctionBegin;
1072   if (pcbddc->ChangeOfBasisMatrix) {
1073     PCBDDCChange_ctx change_ctx;
1074 
1075     /* get change ctx */
1076     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1077 
1078     /* restore iteration matrix */
1079     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1080     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1081     pc->mat = change_ctx->original_mat;
1082 
1083     /* get solution in original basis */
1084     if (x) {
1085       PC_IS *pcis = (PC_IS*)(pc->data);
1086       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1087       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
1088     }
1089   }
1090 
1091   /* add solution removed in presolve */
1092   if (x) {
1093     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
1094   }
1095 
1096   /* restore rhs to its original state */
1097   if (rhs) {
1098     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1099   }
1100   PetscFunctionReturn(0);
1101 }
1102 /* -------------------------------------------------------------------------- */
1103 #undef __FUNCT__
1104 #define __FUNCT__ "PCSetUp_BDDC"
1105 /* -------------------------------------------------------------------------- */
1106 /*
1107    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1108                   by setting data structures and options.
1109 
1110    Input Parameter:
1111 +  pc - the preconditioner context
1112 
1113    Application Interface Routine: PCSetUp()
1114 
1115    Notes:
1116    The interface routine PCSetUp() is not usually called directly by
1117    the user, but instead is called by PCApply() if necessary.
1118 */
1119 PetscErrorCode PCSetUp_BDDC(PC pc)
1120 {
1121   PetscErrorCode   ierr;
1122   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
1123   MatNullSpace     nearnullspace;
1124   PetscBool        computeis,computetopography,computesolvers;
1125   PetscBool        new_nearnullspace_provided;
1126 
1127   PetscFunctionBegin;
1128   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1129   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1130      Also, BDDC directly build the Dirichlet problem */
1131   /* split work */
1132   if (pc->setupcalled) {
1133     computeis = PETSC_FALSE;
1134     if (pc->flag == SAME_NONZERO_PATTERN) {
1135       computetopography = PETSC_FALSE;
1136       computesolvers = PETSC_TRUE;
1137     } else { /* DIFFERENT_NONZERO_PATTERN */
1138       computetopography = PETSC_TRUE;
1139       computesolvers = PETSC_TRUE;
1140     }
1141   } else {
1142     computeis = PETSC_TRUE;
1143     computetopography = PETSC_TRUE;
1144     computesolvers = PETSC_TRUE;
1145   }
1146   if (pcbddc->recompute_topography) {
1147     computetopography = PETSC_TRUE;
1148   }
1149 
1150   /* Get stdout for dbg */
1151   if (pcbddc->dbg_flag) {
1152     if (!pcbddc->dbg_viewer) {
1153       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1154       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1155     }
1156     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1157   }
1158 
1159   /* Set up all the "iterative substructuring" common block without computing solvers */
1160   if (computeis) {
1161     /* HACK INTO PCIS */
1162     PC_IS* pcis = (PC_IS*)pc->data;
1163     pcis->computesolvers = PETSC_FALSE;
1164     ierr = PCISSetUp(pc);CHKERRQ(ierr);
1165     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr);
1166   }
1167 
1168   /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1169   if (pcbddc->user_ChangeOfBasisMatrix) {
1170     pcbddc->use_change_of_basis = PETSC_FALSE;
1171   }
1172 
1173   /* Analyze interface */
1174   if (computetopography) {
1175     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1176   }
1177 
1178   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1179   new_nearnullspace_provided = PETSC_FALSE;
1180   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1181   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1182     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1183       new_nearnullspace_provided = PETSC_TRUE;
1184     } else {
1185       /* determine if the two nullspaces are different (should be lightweight) */
1186       if (nearnullspace != pcbddc->onearnullspace) {
1187         new_nearnullspace_provided = PETSC_TRUE;
1188       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1189         PetscInt         i;
1190         const Vec        *nearnullvecs;
1191         PetscObjectState state;
1192         PetscInt         nnsp_size;
1193         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1194         for (i=0;i<nnsp_size;i++) {
1195           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1196           if (pcbddc->onearnullvecs_state[i] != state) {
1197             new_nearnullspace_provided = PETSC_TRUE;
1198             break;
1199           }
1200         }
1201       }
1202     }
1203   } else {
1204     if (!nearnullspace) { /* both nearnullspaces are null */
1205       new_nearnullspace_provided = PETSC_FALSE;
1206     } else { /* nearnullspace attached later */
1207       new_nearnullspace_provided = PETSC_TRUE;
1208     }
1209   }
1210 
1211   /* Setup constraints and related work vectors */
1212   /* reset primal space flags */
1213   pcbddc->new_primal_space = PETSC_FALSE;
1214   pcbddc->new_primal_space_local = PETSC_FALSE;
1215   if (computetopography || new_nearnullspace_provided) {
1216     /* It also sets the primal space flags */
1217     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1218     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1219     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
1220   }
1221 
1222   if (computesolvers || pcbddc->new_primal_space) {
1223     /* reset data */
1224     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1225     /* Create coarse and local stuffs */
1226     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1227     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
1228   }
1229 
1230   if (pcbddc->dbg_flag) {
1231     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1232   }
1233   PetscFunctionReturn(0);
1234 }
1235 
1236 /* -------------------------------------------------------------------------- */
1237 /*
1238    PCApply_BDDC - Applies the BDDC operator to a vector.
1239 
1240    Input Parameters:
1241 .  pc - the preconditioner context
1242 .  r - input vector (global)
1243 
1244    Output Parameter:
1245 .  z - output vector (global)
1246 
1247    Application Interface Routine: PCApply()
1248  */
1249 #undef __FUNCT__
1250 #define __FUNCT__ "PCApply_BDDC"
1251 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1252 {
1253   PC_IS             *pcis = (PC_IS*)(pc->data);
1254   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1255   PetscErrorCode    ierr;
1256   const PetscScalar one = 1.0;
1257   const PetscScalar m_one = -1.0;
1258   const PetscScalar zero = 0.0;
1259 
1260 /* This code is similar to that provided in nn.c for PCNN
1261    NN interface preconditioner changed to BDDC
1262    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
1263 
1264   PetscFunctionBegin;
1265   if (!pcbddc->use_exact_dirichlet_trick) {
1266     /* First Dirichlet solve */
1267     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1268     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1269     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1270     /*
1271       Assembling right hand side for BDDC operator
1272       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1273       - pcis->vec1_B the interface part of the global vector z
1274     */
1275     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1276     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1277     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1278     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1279     ierr = VecCopy(r,z);CHKERRQ(ierr);
1280     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1281     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1282     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1283   } else {
1284     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1285     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1286     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1287   }
1288 
1289   /* Apply interface preconditioner
1290      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1291   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
1292 
1293   /* Apply transpose of partition of unity operator */
1294   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1295 
1296   /* Second Dirichlet solve and assembling of output */
1297   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1298   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1299   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1300   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1301   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1302   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1303   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1304   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
1305   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1306   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1307   PetscFunctionReturn(0);
1308 }
1309 
1310 /* -------------------------------------------------------------------------- */
1311 /*
1312    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
1313 
1314    Input Parameters:
1315 .  pc - the preconditioner context
1316 .  r - input vector (global)
1317 
1318    Output Parameter:
1319 .  z - output vector (global)
1320 
1321    Application Interface Routine: PCApplyTranspose()
1322  */
1323 #undef __FUNCT__
1324 #define __FUNCT__ "PCApplyTranspose_BDDC"
1325 PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
1326 {
1327   PC_IS             *pcis = (PC_IS*)(pc->data);
1328   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1329   PetscErrorCode    ierr;
1330   const PetscScalar one = 1.0;
1331   const PetscScalar m_one = -1.0;
1332   const PetscScalar zero = 0.0;
1333 
1334   PetscFunctionBegin;
1335   if (!pcbddc->use_exact_dirichlet_trick) {
1336     /* First Dirichlet solve */
1337     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1338     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1339     ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1340     /*
1341       Assembling right hand side for BDDC operator
1342       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1343       - pcis->vec1_B the interface part of the global vector z
1344     */
1345     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1346     ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1347     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1348     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1349     ierr = VecCopy(r,z);CHKERRQ(ierr);
1350     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1351     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1352     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1353   } else {
1354     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1355     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1356     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1357   }
1358 
1359   /* Apply interface preconditioner
1360      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1361   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
1362 
1363   /* Apply transpose of partition of unity operator */
1364   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1365 
1366   /* Second Dirichlet solve and assembling of output */
1367   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1368   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1369   ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1370   if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1371   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1372   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1373   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1374   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
1375   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1376   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1377   PetscFunctionReturn(0);
1378 }
1379 /* -------------------------------------------------------------------------- */
1380 
1381 #undef __FUNCT__
1382 #define __FUNCT__ "PCDestroy_BDDC"
1383 PetscErrorCode PCDestroy_BDDC(PC pc)
1384 {
1385   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1386   PetscErrorCode ierr;
1387 
1388   PetscFunctionBegin;
1389   /* free data created by PCIS */
1390   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1391   /* free BDDC custom data  */
1392   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1393   /* destroy objects related to topography */
1394   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1395   /* free allocated graph structure */
1396   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1397   /* free data for scaling operator */
1398   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1399   /* free solvers stuff */
1400   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1401   /* free global vectors needed in presolve */
1402   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
1403   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1404   /* free stuff for change of basis hooks */
1405   if (pcbddc->new_global_mat) {
1406     PCBDDCChange_ctx change_ctx;
1407     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1408     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1409     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1410     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1411     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1412   }
1413   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
1414   /* remove map from local boundary to local numbering */
1415   ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr);
1416   /* remove functions */
1417   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1418   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1419   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
1420   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1421   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
1422   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1423   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1424   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1425   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1426   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1427   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1428   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1429   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1430   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1431   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1432   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1433   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1434   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1435   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1436   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1437   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1438   /* Free the private data structure */
1439   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1440   PetscFunctionReturn(0);
1441 }
1442 /* -------------------------------------------------------------------------- */
1443 
1444 #undef __FUNCT__
1445 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
1446 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1447 {
1448   FETIDPMat_ctx  mat_ctx;
1449   PC_IS*         pcis;
1450   PC_BDDC*       pcbddc;
1451   PetscErrorCode ierr;
1452 
1453   PetscFunctionBegin;
1454   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1455   pcis = (PC_IS*)mat_ctx->pc->data;
1456   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1457 
1458   /* change of basis for physical rhs if needed
1459      It also changes the rhs in case of dirichlet boundaries */
1460   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
1461   /* store vectors for computation of fetidp final solution */
1462   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1463   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1464   /* scale rhs since it should be unassembled */
1465   /* TODO use counter scaling? (also below) */
1466   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1467   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1468   /* Apply partition of unity */
1469   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1470   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1471   if (!pcbddc->switch_static) {
1472     /* compute partially subassembled Schur complement right-hand side */
1473     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1474     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
1475     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1476     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
1477     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1478     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1479     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1480     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1481     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1482     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1483   }
1484   /* BDDC rhs */
1485   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
1486   if (pcbddc->switch_static) {
1487     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1488   }
1489   /* apply BDDC */
1490   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
1491   /* Application of B_delta and assembling of rhs for fetidp fluxes */
1492   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
1493   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
1494   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1495   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1496   /* restore original rhs */
1497   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
1498   PetscFunctionReturn(0);
1499 }
1500 
1501 #undef __FUNCT__
1502 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
1503 /*@
1504  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
1505 
1506    Collective
1507 
1508    Input Parameters:
1509 +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
1510 .  standard_rhs - the right-hand side for your linear system
1511 
1512    Output Parameters:
1513 -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
1514 
1515    Level: developer
1516 
1517    Notes:
1518 
1519 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1520 @*/
1521 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1522 {
1523   FETIDPMat_ctx  mat_ctx;
1524   PetscErrorCode ierr;
1525 
1526   PetscFunctionBegin;
1527   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1528   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
1529   PetscFunctionReturn(0);
1530 }
1531 /* -------------------------------------------------------------------------- */
1532 
1533 #undef __FUNCT__
1534 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
1535 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1536 {
1537   FETIDPMat_ctx  mat_ctx;
1538   PC_IS*         pcis;
1539   PC_BDDC*       pcbddc;
1540   PetscErrorCode ierr;
1541 
1542   PetscFunctionBegin;
1543   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1544   pcis = (PC_IS*)mat_ctx->pc->data;
1545   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1546 
1547   /* apply B_delta^T */
1548   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1549   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1550   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
1551   /* compute rhs for BDDC application */
1552   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1553   if (pcbddc->switch_static) {
1554     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1555   }
1556   /* apply BDDC */
1557   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
1558   /* put values into standard global vector */
1559   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1560   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1561   if (!pcbddc->switch_static) {
1562     /* compute values into the interior if solved for the partially subassembled Schur complement */
1563     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
1564     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
1565     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1566   }
1567   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1568   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1569   /* final change of basis if needed
1570      Is also sums the dirichlet part removed during RHS assembling */
1571   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
1572   PetscFunctionReturn(0);
1573 }
1574 
1575 #undef __FUNCT__
1576 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
1577 /*@
1578  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
1579 
1580    Collective
1581 
1582    Input Parameters:
1583 +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
1584 .  fetidp_flux_sol - the solution of the FETIDP linear system
1585 
1586    Output Parameters:
1587 -  standard_sol      - the solution defined on the physical domain
1588 
1589    Level: developer
1590 
1591    Notes:
1592 
1593 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1594 @*/
1595 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1596 {
1597   FETIDPMat_ctx  mat_ctx;
1598   PetscErrorCode ierr;
1599 
1600   PetscFunctionBegin;
1601   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1602   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
1603   PetscFunctionReturn(0);
1604 }
1605 /* -------------------------------------------------------------------------- */
1606 
1607 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1608 extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1609 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1610 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1611 extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1612 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1613 
1614 #undef __FUNCT__
1615 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
1616 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1617 {
1618 
1619   FETIDPMat_ctx  fetidpmat_ctx;
1620   Mat            newmat;
1621   FETIDPPC_ctx   fetidppc_ctx;
1622   PC             newpc;
1623   MPI_Comm       comm;
1624   PetscErrorCode ierr;
1625 
1626   PetscFunctionBegin;
1627   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1628   /* FETIDP linear matrix */
1629   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
1630   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
1631   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
1632   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1633   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
1634   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
1635   ierr = MatSetUp(newmat);CHKERRQ(ierr);
1636   /* FETIDP preconditioner */
1637   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
1638   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
1639   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
1640   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
1641   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
1642   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1643   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
1644   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
1645   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
1646   ierr = PCSetUp(newpc);CHKERRQ(ierr);
1647   /* return pointers for objects created */
1648   *fetidp_mat=newmat;
1649   *fetidp_pc=newpc;
1650   PetscFunctionReturn(0);
1651 }
1652 
1653 #undef __FUNCT__
1654 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
1655 /*@
1656  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
1657 
1658    Collective
1659 
1660    Input Parameters:
1661 +  pc - the BDDC preconditioning context already setup
1662 
1663    Output Parameters:
1664 .  fetidp_mat - shell FETIDP matrix object
1665 .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
1666 
1667    Options Database Keys:
1668 -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
1669 
1670    Level: developer
1671 
1672    Notes:
1673      Currently the only operation provided for FETIDP matrix is MatMult
1674 
1675 .seealso: PCBDDC
1676 @*/
1677 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1678 {
1679   PetscErrorCode ierr;
1680 
1681   PetscFunctionBegin;
1682   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1683   if (pc->setupcalled) {
1684     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1685   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1686   PetscFunctionReturn(0);
1687 }
1688 /* -------------------------------------------------------------------------- */
1689 /*MC
1690    PCBDDC - Balancing Domain Decomposition by Constraints.
1691 
1692    An implementation of the BDDC preconditioner based on
1693 
1694 .vb
1695    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
1696    [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
1697    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
1698 .ve
1699 
1700    The matrix to be preconditioned (Pmat) must be of type MATIS.
1701 
1702    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
1703 
1704    It also works with unsymmetric and indefinite problems.
1705 
1706    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.
1707 
1708    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
1709 
1710    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()
1711 
1712    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace().
1713 
1714    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.
1715 
1716    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
1717 
1718    Options Database Keys:
1719 
1720 .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
1721 .    -pc_bddc_use_edges <1> - use or not edges in primal space
1722 .    -pc_bddc_use_faces <0> - use or not faces in primal space
1723 .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
1724 .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
1725 .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
1726 .    -pc_bddc_levels <0> - maximum number of levels for multilevel
1727 .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
1728 -    -pc_bddc_check_level <0> - set verbosity level of debugging output
1729 
1730    Options for Dirichlet, Neumann or coarse solver can be set with
1731 .vb
1732       -pc_bddc_dirichlet_
1733       -pc_bddc_neumann_
1734       -pc_bddc_coarse_
1735 .ve
1736    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
1737 
1738    When using a multilevel approach, solvers' options at the N-th level can be specified as
1739 .vb
1740       -pc_bddc_dirichlet_lN_
1741       -pc_bddc_neumann_lN_
1742       -pc_bddc_coarse_lN_
1743 .ve
1744    Note that level number ranges from the finest 0 to the coarsest N.
1745 
1746    Level: intermediate
1747 
1748    Developer notes:
1749 
1750    New deluxe scaling operator will be available soon.
1751 
1752    Contributed by Stefano Zampini
1753 
1754 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1755 M*/
1756 
1757 #undef __FUNCT__
1758 #define __FUNCT__ "PCCreate_BDDC"
1759 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1760 {
1761   PetscErrorCode      ierr;
1762   PC_BDDC             *pcbddc;
1763 
1764   PetscFunctionBegin;
1765   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1766   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
1767   pc->data  = (void*)pcbddc;
1768 
1769   /* create PCIS data structure */
1770   ierr = PCISCreate(pc);CHKERRQ(ierr);
1771 
1772   /* BDDC customization */
1773   pcbddc->use_local_adj       = PETSC_TRUE;
1774   pcbddc->use_vertices        = PETSC_TRUE;
1775   pcbddc->use_edges           = PETSC_TRUE;
1776   pcbddc->use_faces           = PETSC_FALSE;
1777   pcbddc->use_change_of_basis = PETSC_FALSE;
1778   pcbddc->use_change_on_faces = PETSC_FALSE;
1779   pcbddc->switch_static       = PETSC_FALSE;
1780   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
1781   pcbddc->dbg_flag            = 0;
1782   /* private */
1783   pcbddc->issym                      = PETSC_FALSE;
1784   pcbddc->BtoNmap                    = 0;
1785   pcbddc->local_primal_size          = 0;
1786   pcbddc->n_vertices                 = 0;
1787   pcbddc->n_actual_vertices          = 0;
1788   pcbddc->n_constraints              = 0;
1789   pcbddc->primal_indices_local_idxs  = 0;
1790   pcbddc->recompute_topography       = PETSC_FALSE;
1791   pcbddc->coarse_size                = -1;
1792   pcbddc->new_primal_space           = PETSC_FALSE;
1793   pcbddc->new_primal_space_local     = PETSC_FALSE;
1794   pcbddc->global_primal_indices      = 0;
1795   pcbddc->onearnullspace             = 0;
1796   pcbddc->onearnullvecs_state        = 0;
1797   pcbddc->user_primal_vertices       = 0;
1798   pcbddc->NullSpace                  = 0;
1799   pcbddc->temp_solution              = 0;
1800   pcbddc->original_rhs               = 0;
1801   pcbddc->local_mat                  = 0;
1802   pcbddc->ChangeOfBasisMatrix        = 0;
1803   pcbddc->user_ChangeOfBasisMatrix   = 0;
1804   pcbddc->new_global_mat             = 0;
1805   pcbddc->coarse_vec                 = 0;
1806   pcbddc->coarse_rhs                 = 0;
1807   pcbddc->coarse_ksp                 = 0;
1808   pcbddc->coarse_phi_B               = 0;
1809   pcbddc->coarse_phi_D               = 0;
1810   pcbddc->coarse_psi_B               = 0;
1811   pcbddc->coarse_psi_D               = 0;
1812   pcbddc->vec1_P                     = 0;
1813   pcbddc->vec1_R                     = 0;
1814   pcbddc->vec2_R                     = 0;
1815   pcbddc->local_auxmat1              = 0;
1816   pcbddc->local_auxmat2              = 0;
1817   pcbddc->R_to_B                     = 0;
1818   pcbddc->R_to_D                     = 0;
1819   pcbddc->ksp_D                      = 0;
1820   pcbddc->ksp_R                      = 0;
1821   pcbddc->NeumannBoundaries          = 0;
1822   pcbddc->NeumannBoundariesLocal     = 0;
1823   pcbddc->DirichletBoundaries        = 0;
1824   pcbddc->DirichletBoundariesLocal   = 0;
1825   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
1826   pcbddc->n_ISForDofs                = 0;
1827   pcbddc->n_ISForDofsLocal           = 0;
1828   pcbddc->ISForDofs                  = 0;
1829   pcbddc->ISForDofsLocal             = 0;
1830   pcbddc->ConstraintMatrix           = 0;
1831   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
1832   pcbddc->coarse_loc_to_glob         = 0;
1833   pcbddc->coarsening_ratio           = 8;
1834   pcbddc->current_level              = 0;
1835   pcbddc->max_levels                 = 0;
1836   pcbddc->use_coarse_estimates       = PETSC_FALSE;
1837   pcbddc->coarse_subassembling       = 0;
1838   pcbddc->coarse_subassembling_init  = 0;
1839 
1840   /* create local graph structure */
1841   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1842 
1843   /* scaling */
1844   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1845   pcbddc->deluxe_ctx                 = 0;
1846   pcbddc->work_scaling               = 0;
1847 
1848   /* function pointers */
1849   pc->ops->apply               = PCApply_BDDC;
1850   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
1851   pc->ops->setup               = PCSetUp_BDDC;
1852   pc->ops->destroy             = PCDestroy_BDDC;
1853   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1854   pc->ops->view                = 0;
1855   pc->ops->applyrichardson     = 0;
1856   pc->ops->applysymmetricleft  = 0;
1857   pc->ops->applysymmetricright = 0;
1858   pc->ops->presolve            = PCPreSolve_BDDC;
1859   pc->ops->postsolve           = PCPostSolve_BDDC;
1860 
1861   /* composing function */
1862   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
1863   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1864   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
1865   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1866   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
1867   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1868   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1869   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1870   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1871   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1872   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1873   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1874   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1875   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1876   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1877   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1878   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
1879   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1880   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1881   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1882   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1883   PetscFunctionReturn(0);
1884 }
1885 
1886