xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision c9ed8603f65e57fce3d31de1caa5774cbb25447b)
1 /* TODOLIST
2 
3    Solvers
4    - Add support for cholesky for coarse solver (similar to local solvers)
5    - Propagate ksp prefixes for solvers to mat objects?
6 
7    User interface
8    - ** DM attached to pc?
9 
10    Debugging output
11    - * Better management of verbosity levels of debugging output
12 
13    Extra
14    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
15    - BDDC with MG framework?
16 
17    FETIDP
18    - Move FETIDP code to its own classes
19 
20    MATIS related operations contained in BDDC code
21    - Provide general case for subassembling
22 
23 */
24 
25 #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
26 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
27 #include <petscblaslapack.h>
28 
29 /* -------------------------------------------------------------------------- */
30 #undef __FUNCT__
31 #define __FUNCT__ "PCSetFromOptions_BDDC"
32 PetscErrorCode PCSetFromOptions_BDDC(PetscOptions *PetscOptionsObject,PC pc)
33 {
34   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
35   PetscErrorCode ierr;
36 
37   PetscFunctionBegin;
38   ierr = PetscOptionsHead(PetscOptionsObject,"BDDC options");CHKERRQ(ierr);
39   /* Verbose debugging */
40   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
41   /* Primal space cumstomization */
42   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);
43   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
44   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
45   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
46   ierr = PetscOptionsBool("-pc_bddc_use_true_nnsp","Use near null space attached to the matrix without modifications","none",pcbddc->use_nnsp_true,&pcbddc->use_nnsp_true,NULL);CHKERRQ(ierr);
47   ierr = PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is used when multiple constraints are present)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);CHKERRQ(ierr);
48   /* Change of basis */
49   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);
50   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);
51   if (!pcbddc->use_change_of_basis) {
52     pcbddc->use_change_on_faces = PETSC_FALSE;
53   }
54   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
55   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);
56   ierr = PetscOptionsInt("-pc_bddc_coarse_redistribute","Number of procs where to redistribute coarse problem","none",pcbddc->redistribute_coarse,&pcbddc->redistribute_coarse,NULL);CHKERRQ(ierr);
57   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
58   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
59   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);
60   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
61   ierr = PetscOptionsBool("-pc_bddc_schur_rebuild","Whether or not the interface graph for Schur principal minors has to be rebuilt (i.e. define the interface without any adjacency)","none",pcbddc->sub_schurs_rebuild,&pcbddc->sub_schurs_rebuild,NULL);CHKERRQ(ierr);
62   ierr = PetscOptionsInt("-pc_bddc_schur_layers","Number of dofs' layers for the computation of principal minors (i.e. -1 uses all dofs)","none",pcbddc->sub_schurs_layers,&pcbddc->sub_schurs_layers,NULL);CHKERRQ(ierr);
63   ierr = PetscOptionsBool("-pc_bddc_schur_use_useradj","Whether or not the CSR graph specified by the user should be used for computing successive layers (default is to use adj of local mat)","none",pcbddc->sub_schurs_use_useradj,&pcbddc->sub_schurs_use_useradj,NULL);CHKERRQ(ierr);
64   ierr = PetscOptionsBool("-pc_bddc_deluxe_faster","Faster application of deluxe scaling (requires extra work during setup)","none",pcbddc->faster_deluxe,&pcbddc->faster_deluxe,NULL);CHKERRQ(ierr);
65   ierr = PetscOptionsReal("-pc_bddc_adaptive_threshold","Threshold to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&pcbddc->adaptive_threshold,NULL);CHKERRQ(ierr);
66   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);CHKERRQ(ierr);
67   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);CHKERRQ(ierr);
68   ierr = PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);CHKERRQ(ierr);
69   ierr = PetscOptionsInt("-pc_bddc_coarse_adj","Number of processors where to map the coarse adjacency list","none",pcbddc->coarse_adj_red,&pcbddc->coarse_adj_red,NULL);CHKERRQ(ierr);
70   ierr = PetscOptionsBool("-pc_bddc_benign_trick","Apply the benign subspace trick to a class of saddle point problems","none",pcbddc->benign_saddle_point,&pcbddc->benign_saddle_point,NULL);CHKERRQ(ierr);
71   ierr = PetscOptionsTail();CHKERRQ(ierr);
72   PetscFunctionReturn(0);
73 }
74 /* -------------------------------------------------------------------------- */
75 #undef __FUNCT__
76 #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
77 static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
78 {
79   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
80   PetscErrorCode ierr;
81 
82   PetscFunctionBegin;
83   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
84   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
85   pcbddc->user_ChangeOfBasisMatrix = change;
86   PetscFunctionReturn(0);
87 }
88 #undef __FUNCT__
89 #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
90 /*@
91  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
92 
93    Collective on PC
94 
95    Input Parameters:
96 +  pc - the preconditioning context
97 -  change - the change of basis matrix
98 
99    Level: intermediate
100 
101    Notes:
102 
103 .seealso: PCBDDC
104 @*/
105 PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
106 {
107   PetscErrorCode ierr;
108 
109   PetscFunctionBegin;
110   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
111   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
112   PetscCheckSameComm(pc,1,change,2);
113   if (pc->mat) {
114     PetscInt rows_c,cols_c,rows,cols;
115     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
116     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
117     if (rows_c != rows) {
118       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
119     }
120     if (cols_c != cols) {
121       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
122     }
123     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
124     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
125     if (rows_c != rows) {
126       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
127     }
128     if (cols_c != cols) {
129       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
130     }
131   }
132   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
133   PetscFunctionReturn(0);
134 }
135 /* -------------------------------------------------------------------------- */
136 #undef __FUNCT__
137 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
138 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
139 {
140   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
141   PetscErrorCode ierr;
142 
143   PetscFunctionBegin;
144   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
145   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
146   pcbddc->user_primal_vertices = PrimalVertices;
147   PetscFunctionReturn(0);
148 }
149 #undef __FUNCT__
150 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
151 /*@
152  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
153 
154    Collective
155 
156    Input Parameters:
157 +  pc - the preconditioning context
158 -  PrimalVertices - index set of primal vertices in local numbering (can be empty)
159 
160    Level: intermediate
161 
162    Notes:
163 
164 .seealso: PCBDDC
165 @*/
166 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
167 {
168   PetscErrorCode ierr;
169 
170   PetscFunctionBegin;
171   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
172   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
173   PetscCheckSameComm(pc,1,PrimalVertices,2);
174   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
175   PetscFunctionReturn(0);
176 }
177 /* -------------------------------------------------------------------------- */
178 #undef __FUNCT__
179 #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
180 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
181 {
182   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
183 
184   PetscFunctionBegin;
185   pcbddc->coarsening_ratio = k;
186   PetscFunctionReturn(0);
187 }
188 
189 #undef __FUNCT__
190 #define __FUNCT__ "PCBDDCSetCoarseningRatio"
191 /*@
192  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
193 
194    Logically collective on PC
195 
196    Input Parameters:
197 +  pc - the preconditioning context
198 -  k - coarsening ratio (H/h at the coarser level)
199 
200    Options Database Keys:
201 .    -pc_bddc_coarsening_ratio
202 
203    Level: intermediate
204 
205    Notes:
206      Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
207 
208 .seealso: PCBDDC, PCBDDCSetLevels()
209 @*/
210 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
211 {
212   PetscErrorCode ierr;
213 
214   PetscFunctionBegin;
215   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
216   PetscValidLogicalCollectiveInt(pc,k,2);
217   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
218   PetscFunctionReturn(0);
219 }
220 
221 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
222 #undef __FUNCT__
223 #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
224 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
225 {
226   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
227 
228   PetscFunctionBegin;
229   pcbddc->use_exact_dirichlet_trick = flg;
230   PetscFunctionReturn(0);
231 }
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
235 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
236 {
237   PetscErrorCode ierr;
238 
239   PetscFunctionBegin;
240   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
241   PetscValidLogicalCollectiveBool(pc,flg,2);
242   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "PCBDDCSetLevel_BDDC"
248 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
249 {
250   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
251 
252   PetscFunctionBegin;
253   pcbddc->current_level = level;
254   PetscFunctionReturn(0);
255 }
256 
257 #undef __FUNCT__
258 #define __FUNCT__ "PCBDDCSetLevel"
259 PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
260 {
261   PetscErrorCode ierr;
262 
263   PetscFunctionBegin;
264   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
265   PetscValidLogicalCollectiveInt(pc,level,2);
266   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
267   PetscFunctionReturn(0);
268 }
269 
270 #undef __FUNCT__
271 #define __FUNCT__ "PCBDDCSetLevels_BDDC"
272 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
273 {
274   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
275 
276   PetscFunctionBegin;
277   pcbddc->max_levels = levels;
278   PetscFunctionReturn(0);
279 }
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "PCBDDCSetLevels"
283 /*@
284  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
285 
286    Logically collective on PC
287 
288    Input Parameters:
289 +  pc - the preconditioning context
290 -  levels - the maximum number of levels (max 9)
291 
292    Options Database Keys:
293 .    -pc_bddc_levels
294 
295    Level: intermediate
296 
297    Notes:
298      Default value is 0, i.e. traditional one-level BDDC
299 
300 .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
301 @*/
302 PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
303 {
304   PetscErrorCode ierr;
305 
306   PetscFunctionBegin;
307   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
308   PetscValidLogicalCollectiveInt(pc,levels,2);
309   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
310   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
311   PetscFunctionReturn(0);
312 }
313 /* -------------------------------------------------------------------------- */
314 
315 #undef __FUNCT__
316 #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
317 static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
318 {
319   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
320   PetscErrorCode ierr;
321 
322   PetscFunctionBegin;
323   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
324   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
325   pcbddc->NullSpace = NullSpace;
326   PetscFunctionReturn(0);
327 }
328 
329 #undef __FUNCT__
330 #define __FUNCT__ "PCBDDCSetNullSpace"
331 /*@
332  PCBDDCSetNullSpace - Set nullspace for BDDC operator
333 
334    Logically collective on PC and MatNullSpace
335 
336    Input Parameters:
337 +  pc - the preconditioning context
338 -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
339 
340    Level: intermediate
341 
342    Notes:
343 
344 .seealso: PCBDDC
345 @*/
346 PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
347 {
348   PetscErrorCode ierr;
349 
350   PetscFunctionBegin;
351   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
352   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
353   PetscCheckSameComm(pc,1,NullSpace,2);
354   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
355   PetscFunctionReturn(0);
356 }
357 /* -------------------------------------------------------------------------- */
358 
359 #undef __FUNCT__
360 #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
361 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
362 {
363   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
364   PetscErrorCode ierr;
365 
366   PetscFunctionBegin;
367   /* last user setting takes precendence -> destroy any other customization */
368   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
369   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
370   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
371   pcbddc->DirichletBoundaries = DirichletBoundaries;
372   pcbddc->recompute_topography = PETSC_TRUE;
373   PetscFunctionReturn(0);
374 }
375 
376 #undef __FUNCT__
377 #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
378 /*@
379  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
380 
381    Collective
382 
383    Input Parameters:
384 +  pc - the preconditioning context
385 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
386 
387    Level: intermediate
388 
389    Notes:
390      Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node
391 
392 .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
393 @*/
394 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
395 {
396   PetscErrorCode ierr;
397 
398   PetscFunctionBegin;
399   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
400   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
401   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
402   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
403   PetscFunctionReturn(0);
404 }
405 /* -------------------------------------------------------------------------- */
406 
407 #undef __FUNCT__
408 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
409 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
410 {
411   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
412   PetscErrorCode ierr;
413 
414   PetscFunctionBegin;
415   /* last user setting takes precendence -> destroy any other customization */
416   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
417   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
418   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
419   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
420   pcbddc->recompute_topography = PETSC_TRUE;
421   PetscFunctionReturn(0);
422 }
423 
424 #undef __FUNCT__
425 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
426 /*@
427  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
428 
429    Collective
430 
431    Input Parameters:
432 +  pc - the preconditioning context
433 -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
434 
435    Level: intermediate
436 
437    Notes:
438 
439 .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
440 @*/
441 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
442 {
443   PetscErrorCode ierr;
444 
445   PetscFunctionBegin;
446   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
447   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
448   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
449   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
450   PetscFunctionReturn(0);
451 }
452 /* -------------------------------------------------------------------------- */
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
456 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
457 {
458   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
459   PetscErrorCode ierr;
460 
461   PetscFunctionBegin;
462   /* last user setting takes precendence -> destroy any other customization */
463   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
464   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
465   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
466   pcbddc->NeumannBoundaries = NeumannBoundaries;
467   pcbddc->recompute_topography = PETSC_TRUE;
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
473 /*@
474  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
475 
476    Collective
477 
478    Input Parameters:
479 +  pc - the preconditioning context
480 -  NeumannBoundaries - parallel IS defining the Neumann boundaries
481 
482    Level: intermediate
483 
484    Notes:
485      Any process can list any global node
486 
487 .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
488 @*/
489 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
490 {
491   PetscErrorCode ierr;
492 
493   PetscFunctionBegin;
494   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
495   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
496   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
497   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
498   PetscFunctionReturn(0);
499 }
500 /* -------------------------------------------------------------------------- */
501 
502 #undef __FUNCT__
503 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
504 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
505 {
506   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
507   PetscErrorCode ierr;
508 
509   PetscFunctionBegin;
510   /* last user setting takes precendence -> destroy any other customization */
511   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
512   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
513   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
514   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
515   pcbddc->recompute_topography = PETSC_TRUE;
516   PetscFunctionReturn(0);
517 }
518 
519 #undef __FUNCT__
520 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
521 /*@
522  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
523 
524    Collective
525 
526    Input Parameters:
527 +  pc - the preconditioning context
528 -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
529 
530    Level: intermediate
531 
532    Notes:
533 
534 .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
535 @*/
536 PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
537 {
538   PetscErrorCode ierr;
539 
540   PetscFunctionBegin;
541   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
542   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
543   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
544   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
545   PetscFunctionReturn(0);
546 }
547 /* -------------------------------------------------------------------------- */
548 
549 #undef __FUNCT__
550 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
551 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
552 {
553   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
554 
555   PetscFunctionBegin;
556   *DirichletBoundaries = pcbddc->DirichletBoundaries;
557   PetscFunctionReturn(0);
558 }
559 
560 #undef __FUNCT__
561 #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
562 /*@
563  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
564 
565    Collective
566 
567    Input Parameters:
568 .  pc - the preconditioning context
569 
570    Output Parameters:
571 .  DirichletBoundaries - index set defining the Dirichlet boundaries
572 
573    Level: intermediate
574 
575    Notes:
576      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
577 
578 .seealso: PCBDDC
579 @*/
580 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
581 {
582   PetscErrorCode ierr;
583 
584   PetscFunctionBegin;
585   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
586   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
587   PetscFunctionReturn(0);
588 }
589 /* -------------------------------------------------------------------------- */
590 
591 #undef __FUNCT__
592 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
593 static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
594 {
595   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
596 
597   PetscFunctionBegin;
598   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
599   PetscFunctionReturn(0);
600 }
601 
602 #undef __FUNCT__
603 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
604 /*@
605  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
606 
607    Collective
608 
609    Input Parameters:
610 .  pc - the preconditioning context
611 
612    Output Parameters:
613 .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
614 
615    Level: intermediate
616 
617    Notes:
618      The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetDirichletBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetDirichletBoundaries).
619           In the latter case, the IS will be available after PCSetUp.
620 
621 .seealso: PCBDDC
622 @*/
623 PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
624 {
625   PetscErrorCode ierr;
626 
627   PetscFunctionBegin;
628   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
629   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
630   PetscFunctionReturn(0);
631 }
632 /* -------------------------------------------------------------------------- */
633 
634 #undef __FUNCT__
635 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
636 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
637 {
638   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
639 
640   PetscFunctionBegin;
641   *NeumannBoundaries = pcbddc->NeumannBoundaries;
642   PetscFunctionReturn(0);
643 }
644 
645 #undef __FUNCT__
646 #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
647 /*@
648  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
649 
650    Collective
651 
652    Input Parameters:
653 .  pc - the preconditioning context
654 
655    Output Parameters:
656 .  NeumannBoundaries - index set defining the Neumann boundaries
657 
658    Level: intermediate
659 
660    Notes:
661      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
662 
663 .seealso: PCBDDC
664 @*/
665 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
666 {
667   PetscErrorCode ierr;
668 
669   PetscFunctionBegin;
670   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
671   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
672   PetscFunctionReturn(0);
673 }
674 /* -------------------------------------------------------------------------- */
675 
676 #undef __FUNCT__
677 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
678 static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
679 {
680   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
681 
682   PetscFunctionBegin;
683   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
684   PetscFunctionReturn(0);
685 }
686 
687 #undef __FUNCT__
688 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
689 /*@
690  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
691 
692    Collective
693 
694    Input Parameters:
695 .  pc - the preconditioning context
696 
697    Output Parameters:
698 .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
699 
700    Level: intermediate
701 
702    Notes:
703      The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetNeumannBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetNeumannBoundaries).
704           In the latter case, the IS will be available after PCSetUp.
705 
706 .seealso: PCBDDC
707 @*/
708 PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
709 {
710   PetscErrorCode ierr;
711 
712   PetscFunctionBegin;
713   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
714   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
715   PetscFunctionReturn(0);
716 }
717 /* -------------------------------------------------------------------------- */
718 
719 #undef __FUNCT__
720 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
721 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
722 {
723   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
724   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
725   PetscErrorCode ierr;
726 
727   PetscFunctionBegin;
728   /* free old CSR */
729   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
730   /* TODO: PCBDDCGraphSetAdjacency */
731   /* get CSR into graph structure */
732   if (copymode == PETSC_COPY_VALUES) {
733     ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr);
734     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
735     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
736     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
737   } else if (copymode == PETSC_OWN_POINTER) {
738     mat_graph->xadj = (PetscInt*)xadj;
739     mat_graph->adjncy = (PetscInt*)adjncy;
740   }
741   mat_graph->nvtxs_csr = nvtxs;
742   PetscFunctionReturn(0);
743 }
744 
745 #undef __FUNCT__
746 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
747 /*@
748  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local matrix
749 
750    Not collective
751 
752    Input Parameters:
753 +  pc - the preconditioning context
754 .  nvtxs - number of local vertices of the graph (i.e., the size of the local problem)
755 .  xadj, adjncy - the CSR graph
756 -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
757 
758    Level: intermediate
759 
760    Notes:
761 
762 .seealso: PCBDDC,PetscCopyMode
763 @*/
764 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
765 {
766   void (*f)(void) = 0;
767   PetscErrorCode ierr;
768 
769   PetscFunctionBegin;
770   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
771   PetscValidIntPointer(xadj,3);
772   PetscValidIntPointer(adjncy,4);
773   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
774     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
775   }
776   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
777   /* free arrays if PCBDDC is not the PC type */
778   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
779   if (!f && copymode == PETSC_OWN_POINTER) {
780     ierr = PetscFree(xadj);CHKERRQ(ierr);
781     ierr = PetscFree(adjncy);CHKERRQ(ierr);
782   }
783   PetscFunctionReturn(0);
784 }
785 /* -------------------------------------------------------------------------- */
786 
787 #undef __FUNCT__
788 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
789 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
790 {
791   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
792   PetscInt i;
793   PetscErrorCode ierr;
794 
795   PetscFunctionBegin;
796   /* Destroy ISes if they were already set */
797   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
798     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
799   }
800   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
801   /* last user setting takes precendence -> destroy any other customization */
802   for (i=0;i<pcbddc->n_ISForDofs;i++) {
803     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
804   }
805   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
806   pcbddc->n_ISForDofs = 0;
807   /* allocate space then set */
808   if (n_is) {
809     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
810   }
811   for (i=0;i<n_is;i++) {
812     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
813     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
814   }
815   pcbddc->n_ISForDofsLocal=n_is;
816   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
817   pcbddc->recompute_topography = PETSC_TRUE;
818   PetscFunctionReturn(0);
819 }
820 
821 #undef __FUNCT__
822 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
823 /*@
824  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
825 
826    Collective
827 
828    Input Parameters:
829 +  pc - the preconditioning context
830 .  n_is - number of index sets defining the fields
831 -  ISForDofs - array of IS describing the fields in local ordering
832 
833    Level: intermediate
834 
835    Notes:
836      n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
837 
838 .seealso: PCBDDC
839 @*/
840 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
841 {
842   PetscInt       i;
843   PetscErrorCode ierr;
844 
845   PetscFunctionBegin;
846   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
847   PetscValidLogicalCollectiveInt(pc,n_is,2);
848   for (i=0;i<n_is;i++) {
849     PetscCheckSameComm(pc,1,ISForDofs[i],3);
850     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
851   }
852   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
853   PetscFunctionReturn(0);
854 }
855 /* -------------------------------------------------------------------------- */
856 
857 #undef __FUNCT__
858 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
859 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
860 {
861   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
862   PetscInt i;
863   PetscErrorCode ierr;
864 
865   PetscFunctionBegin;
866   /* Destroy ISes if they were already set */
867   for (i=0;i<pcbddc->n_ISForDofs;i++) {
868     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
869   }
870   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
871   /* last user setting takes precendence -> destroy any other customization */
872   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
873     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
874   }
875   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
876   pcbddc->n_ISForDofsLocal = 0;
877   /* allocate space then set */
878   if (n_is) {
879     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
880   }
881   for (i=0;i<n_is;i++) {
882     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
883     pcbddc->ISForDofs[i]=ISForDofs[i];
884   }
885   pcbddc->n_ISForDofs=n_is;
886   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
887   pcbddc->recompute_topography = PETSC_TRUE;
888   PetscFunctionReturn(0);
889 }
890 
891 #undef __FUNCT__
892 #define __FUNCT__ "PCBDDCSetDofsSplitting"
893 /*@
894  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
895 
896    Collective
897 
898    Input Parameters:
899 +  pc - the preconditioning context
900 .  n_is - number of index sets defining the fields
901 -  ISForDofs - array of IS describing the fields in global ordering
902 
903    Level: intermediate
904 
905    Notes:
906      Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
907 
908 .seealso: PCBDDC
909 @*/
910 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
911 {
912   PetscInt       i;
913   PetscErrorCode ierr;
914 
915   PetscFunctionBegin;
916   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
917   PetscValidLogicalCollectiveInt(pc,n_is,2);
918   for (i=0;i<n_is;i++) {
919     PetscCheckSameComm(pc,1,ISForDofs[i],3);
920     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
921   }
922   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
923   PetscFunctionReturn(0);
924 }
925 
926 /* -------------------------------------------------------------------------- */
927 #undef __FUNCT__
928 #define __FUNCT__ "PCPreSolve_BDDC"
929 /* -------------------------------------------------------------------------- */
930 /*
931    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
932                      guess if a transformation of basis approach has been selected.
933 
934    Input Parameter:
935 +  pc - the preconditioner contex
936 
937    Application Interface Routine: PCPreSolve()
938 
939    Notes:
940      The interface routine PCPreSolve() is not usually called directly by
941    the user, but instead is called by KSPSolve().
942 */
943 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
944 {
945   PetscErrorCode ierr;
946   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
947   PC_IS          *pcis = (PC_IS*)(pc->data);
948   Vec            used_vec;
949   PetscBool      copy_rhs = PETSC_TRUE;
950 
951   PetscFunctionBegin;
952   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
953   if (ksp) {
954     PetscBool iscg;
955     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
956     if (!iscg) {
957       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
958     }
959   }
960   /* Creates parallel work vectors used in presolve */
961   if (!pcbddc->original_rhs) {
962     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
963   }
964   if (!pcbddc->temp_solution) {
965     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
966   }
967 
968   if (x) {
969     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
970     used_vec = x;
971   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
972     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
973     used_vec = pcbddc->temp_solution;
974     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
975   }
976 
977   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
978   if (ksp) {
979     /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
980     ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
981     if (!pcbddc->ksp_guess_nonzero) {
982       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
983     }
984   }
985 
986   pcbddc->rhs_change = PETSC_FALSE;
987   /* Take into account zeroed rows -> change rhs and store solution removed */
988   if (rhs) {
989     IS dirIS = NULL;
990 
991     /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
992     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
993     if (dirIS) {
994       Mat_IS            *matis = (Mat_IS*)pc->pmat->data;
995       PetscInt          dirsize,i,*is_indices;
996       PetscScalar       *array_x;
997       const PetscScalar *array_diagonal;
998 
999       ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
1000       ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
1001       ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1002       ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1003       ierr = VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1004       ierr = VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1005       ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
1006       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
1007       ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
1008       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
1009       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
1010       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
1011       ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
1012       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
1013       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1014       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1015       pcbddc->rhs_change = PETSC_TRUE;
1016       ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
1017     }
1018   }
1019 
1020   /* remove the computed solution or the initial guess from the rhs */
1021   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
1022     /* store the original rhs */
1023     if (copy_rhs) {
1024       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1025       copy_rhs = PETSC_FALSE;
1026     }
1027     pcbddc->rhs_change = PETSC_TRUE;
1028     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
1029     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
1030     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
1031     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
1032     if (ksp) {
1033       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr);
1034     }
1035   }
1036   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1037 
1038   /* If using the benign trick, change rhs on pressures (iteration matrix is surely of type MATIS) */
1039   if (pcbddc->benign_saddle_point) {
1040     Mat_IS *matis = (Mat_IS*)(pc->mat->data);
1041 
1042     /* store the original rhs */
1043     if (copy_rhs) {
1044       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1045       copy_rhs = PETSC_FALSE;
1046     }
1047 
1048     /* now change (locally) the basis */
1049     ierr = VecScatterBegin(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1050     ierr = VecScatterEnd(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1051     if (pcbddc->benign_change) {
1052       ierr = MatMultTranspose(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
1053       ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1054       ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1055       /* change local iteration matrix */
1056       ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_REUSE_MATRIX,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
1057       ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr);
1058       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1059       pcbddc->benign_original_mat = matis->A;
1060       ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
1061       ierr = PetscObjectReference((PetscObject)pcbddc->local_mat);CHKERRQ(ierr);
1062       matis->A = pcbddc->local_mat;
1063       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1064     } else {
1065       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1066       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1067     }
1068     pcbddc->rhs_change = PETSC_TRUE;
1069   }
1070 
1071   /* store partially computed solution and set initial guess */
1072   if (x && pcbddc->use_exact_dirichlet_trick) {
1073     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1074     ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1075     ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1076     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1077     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1078     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1079     if (ksp) {
1080       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1081     }
1082   }
1083 
1084   if (pcbddc->ChangeOfBasisMatrix) {
1085     PCBDDCChange_ctx change_ctx;
1086 
1087     /* get change ctx */
1088     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1089 
1090     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1091     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1092     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1093     change_ctx->original_mat = pc->mat;
1094 
1095     /* change iteration matrix */
1096     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1097     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1098     pc->mat = pcbddc->new_global_mat;
1099 
1100     /* store the original rhs */
1101     if (copy_rhs) {
1102       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1103       copy_rhs = PETSC_FALSE;
1104     }
1105 
1106     /* change rhs */
1107     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1108     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
1109     pcbddc->rhs_change = PETSC_TRUE;
1110   }
1111 
1112   /* remove nullspace if present */
1113   if (ksp && x && pcbddc->NullSpace) {
1114     ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr);
1115     /* store the original rhs */
1116     if (copy_rhs) {
1117       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1118       copy_rhs = PETSC_FALSE;
1119     }
1120     pcbddc->rhs_change = PETSC_TRUE;
1121     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1122   }
1123   PetscFunctionReturn(0);
1124 }
1125 
1126 /* -------------------------------------------------------------------------- */
1127 #undef __FUNCT__
1128 #define __FUNCT__ "PCPostSolve_BDDC"
1129 /* -------------------------------------------------------------------------- */
1130 /*
1131    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1132                      approach has been selected. Also, restores rhs to its original state.
1133 
1134    Input Parameter:
1135 +  pc - the preconditioner contex
1136 
1137    Application Interface Routine: PCPostSolve()
1138 
1139    Notes:
1140      The interface routine PCPostSolve() is not usually called directly by
1141      the user, but instead is called by KSPSolve().
1142 */
1143 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1144 {
1145   PetscErrorCode ierr;
1146   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1147 
1148   PetscFunctionBegin;
1149   if (pcbddc->ChangeOfBasisMatrix) {
1150     PCBDDCChange_ctx change_ctx;
1151 
1152     /* get change ctx */
1153     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1154 
1155     /* restore iteration matrix */
1156     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1157     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1158     pc->mat = change_ctx->original_mat;
1159 
1160     /* get solution in original basis */
1161     if (x) {
1162       PC_IS *pcis = (PC_IS*)(pc->data);
1163 
1164       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1165       /* restore solution on pressures */
1166       if (pcbddc->benign_saddle_point) {
1167         Mat_IS *matis = (Mat_IS*)pc->mat->data;
1168 
1169         ierr = VecScatterBegin(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1170         ierr = VecScatterEnd(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1171         if (pcbddc->benign_change) {
1172 
1173           ierr = MatMult(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
1174           ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1175           ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1176         } else {
1177           ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1178           ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1179         }
1180       }
1181       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
1182     }
1183   }
1184 
1185   /* add solution removed in presolve */
1186   if (x && pcbddc->rhs_change) {
1187     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
1188   }
1189 
1190   /* restore rhs to its original state */
1191   if (rhs && pcbddc->rhs_change) {
1192     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1193   }
1194   pcbddc->rhs_change = PETSC_FALSE;
1195 
1196   /* restore ksp guess state */
1197   if (ksp) {
1198     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
1199   }
1200   PetscFunctionReturn(0);
1201 }
1202 /* -------------------------------------------------------------------------- */
1203 #undef __FUNCT__
1204 #define __FUNCT__ "PCSetUp_BDDC"
1205 /* -------------------------------------------------------------------------- */
1206 /*
1207    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1208                   by setting data structures and options.
1209 
1210    Input Parameter:
1211 +  pc - the preconditioner context
1212 
1213    Application Interface Routine: PCSetUp()
1214 
1215    Notes:
1216      The interface routine PCSetUp() is not usually called directly by
1217      the user, but instead is called by PCApply() if necessary.
1218 */
1219 PetscErrorCode PCSetUp_BDDC(PC pc)
1220 {
1221   PetscErrorCode ierr;
1222   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
1223   Mat_IS*        matis;
1224   MatNullSpace   nearnullspace;
1225   IS             zerodiag = NULL;
1226   PetscInt       nrows,ncols;
1227   PetscBool      computetopography,computesolvers,computesubschurs;
1228   PetscBool      computeconstraintsmatrix;
1229   PetscBool      new_nearnullspace_provided,ismatis;
1230 
1231   PetscFunctionBegin;
1232   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr);
1233   if (!ismatis) {
1234     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
1235   }
1236   ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr);
1237   if (nrows != ncols) {
1238     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
1239   }
1240   matis = (Mat_IS*)pc->pmat->data;
1241   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1242   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1243      Also, BDDC directly build the Dirichlet problem */
1244   /* split work */
1245   if (pc->setupcalled) {
1246     if (pc->flag == SAME_NONZERO_PATTERN) {
1247       computetopography = PETSC_FALSE;
1248       computesolvers = PETSC_TRUE;
1249     } else { /* DIFFERENT_NONZERO_PATTERN */
1250       computetopography = PETSC_TRUE;
1251       computesolvers = PETSC_TRUE;
1252     }
1253   } else {
1254     computetopography = PETSC_TRUE;
1255     computesolvers = PETSC_TRUE;
1256   }
1257   if (pcbddc->recompute_topography) {
1258     computetopography = PETSC_TRUE;
1259   }
1260   computeconstraintsmatrix = PETSC_FALSE;
1261   if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) {
1262     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling");
1263   }
1264   pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling);
1265   if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;
1266 
1267   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
1268   if (pcbddc->faster_deluxe && pcbddc->adaptive_selection && pcbddc->use_change_of_basis) {
1269     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute faster deluxe if adaptivity and change of basis are both requested. Rerun with -pc_bddc_deluxe_faster false");
1270   }
1271 
1272   /* check if the iteration matrix is of type MATIS in case the benign trick has been requested */
1273   ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&ismatis);CHKERRQ(ierr);
1274   if (pcbddc->benign_saddle_point && !ismatis) {
1275     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner with benign subspace trick requires the iteration matrix to be of type MATIS");
1276   }
1277 
1278   /* Get stdout for dbg */
1279   if (pcbddc->dbg_flag) {
1280     if (!pcbddc->dbg_viewer) {
1281       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1282       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1283     }
1284     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1285   }
1286 
1287   if (pcbddc->user_ChangeOfBasisMatrix) {
1288     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1289     pcbddc->use_change_of_basis = PETSC_FALSE;
1290     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
1291   } else {
1292     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1293     if (!pcbddc->benign_saddle_point) {
1294       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1295       pcbddc->local_mat = matis->A;
1296     } else { /* TODO: handle user change of basis */
1297       PetscInt  nz;
1298       PetscBool sorted;
1299 
1300       /* should I also include nonzero pressures? */
1301       ierr = MatFindZeroDiagonals(matis->A,&zerodiag);CHKERRQ(ierr);
1302       ierr = ISSorted(zerodiag,&sorted);CHKERRQ(ierr);
1303       if (!sorted) {
1304         ierr = ISSort(zerodiag);CHKERRQ(ierr);
1305       }
1306       ierr = ISGetLocalSize(zerodiag,&nz);CHKERRQ(ierr);
1307       if (nz) {
1308         IS                zerodiagc;
1309         PetscScalar       *array;
1310         const PetscInt    *idxs,*idxsc;
1311         PetscInt          i,n,*nnz;
1312 
1313         /* TODO: add check for shared dofs */
1314         pcbddc->use_change_of_basis = PETSC_TRUE;
1315         pcbddc->use_change_on_faces = PETSC_TRUE;
1316         ierr = MatGetLocalSize(matis->A,&n,NULL);CHKERRQ(ierr);
1317         ierr = ISComplement(zerodiag,0,n,&zerodiagc);CHKERRQ(ierr);
1318         ierr = ISGetIndices(zerodiag,&idxs);CHKERRQ(ierr);
1319         ierr = ISGetIndices(zerodiagc,&idxsc);CHKERRQ(ierr);
1320         /* local change of basis for pressures */
1321         ierr = MatDestroy(&pcbddc->benign_change);CHKERRQ(ierr);
1322         ierr = MatCreate(PetscObjectComm((PetscObject)matis->A),&pcbddc->benign_change);CHKERRQ(ierr);
1323         ierr = MatSetType(pcbddc->benign_change,MATAIJ);CHKERRQ(ierr);
1324         ierr = MatSetSizes(pcbddc->benign_change,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
1325         ierr = PetscMalloc1(n,&nnz);CHKERRQ(ierr);
1326         for (i=0;i<n-nz;i++) nnz[idxsc[i]] = 1; /* identity on velocities */
1327         for (i=0;i<nz-1;i++) nnz[idxs[i]] = 2; /* change on pressures */
1328         nnz[idxs[nz-1]] = nz; /* last local pressure dof: _0 set */
1329         ierr = MatSeqAIJSetPreallocation(pcbddc->benign_change,0,nnz);CHKERRQ(ierr);
1330         ierr = PetscFree(nnz);CHKERRQ(ierr);
1331         /* set identity on velocities */
1332         for (i=0;i<n-nz;i++) {
1333           ierr = MatSetValue(pcbddc->benign_change,idxsc[i],idxsc[i],1.,INSERT_VALUES);CHKERRQ(ierr);
1334         }
1335         /* set change on pressures */
1336         for (i=0;i<nz-1;i++) {
1337           PetscScalar vals[2];
1338           PetscInt    cols[2];
1339 
1340           /* TODO: add quadrature */
1341           cols[0] = idxs[i];
1342           cols[1] = idxs[nz-1];
1343           vals[0] = 1.;
1344           vals[1] = 1./nz;
1345           ierr = MatSetValues(pcbddc->benign_change,1,idxs+i,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
1346         }
1347         ierr = PetscMalloc1(nz,&array);CHKERRQ(ierr);
1348         for (i=0;i<nz-1;i++) array[i] = -1.;
1349         array[nz-1] = 1./nz;
1350         ierr = MatSetValues(pcbddc->benign_change,1,idxs+nz-1,nz,idxs,array,INSERT_VALUES);CHKERRQ(ierr);
1351         ierr = PetscFree(array);CHKERRQ(ierr);
1352         ierr = MatAssemblyBegin(pcbddc->benign_change,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1353         ierr = MatAssemblyEnd(pcbddc->benign_change,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1354         /* TODO: need optimization? */
1355         ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_INITIAL_MATRIX,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
1356         /* store local and global idxs for p0 */
1357         pcbddc->benign_p0_lidx = idxs[nz-1];
1358         ierr = ISLocalToGlobalMappingApply(pc->pmat->rmap->mapping,1,&idxs[nz-1],&pcbddc->benign_p0_gidx);CHKERRQ(ierr);
1359         ierr = ISRestoreIndices(zerodiag,&idxs);CHKERRQ(ierr);
1360         ierr = ISRestoreIndices(zerodiagc,&idxsc);CHKERRQ(ierr);
1361         ierr = ISDestroy(&zerodiagc);CHKERRQ(ierr);
1362         /* pop B0 mat from pcbddc->local_mat */
1363         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
1364       } else { /* this is unlikely to happen but, just in case, destroy the empty IS */
1365         pcbddc->benign_p0_lidx = -1;
1366         pcbddc->benign_p0_gidx = -1;
1367         ierr = ISDestroy(&zerodiag);CHKERRQ(ierr);
1368         ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1369         pcbddc->local_mat = matis->A;
1370       }
1371     }
1372   }
1373 
1374   /* workaround for reals */
1375 #if !defined(PETSC_USE_COMPLEX)
1376   if (matis->A->symmetric_set) {
1377     ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr);
1378   }
1379 #endif
1380   if (matis->A->symmetric_set) {
1381     ierr = MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);CHKERRQ(ierr);
1382   }
1383   if (matis->A->spd_set) {
1384     ierr = MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);CHKERRQ(ierr);
1385   }
1386 
1387   /* Set up all the "iterative substructuring" common block without computing solvers */
1388   {
1389     Mat temp_mat;
1390 
1391     temp_mat = matis->A;
1392     matis->A = pcbddc->local_mat;
1393     ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr);
1394     pcbddc->local_mat = matis->A;
1395     matis->A = temp_mat;
1396   }
1397 
1398   /* Analyze interface */
1399   if (computetopography) {
1400     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1401     computeconstraintsmatrix = PETSC_TRUE;
1402   }
1403 
1404   /* check b(v_I,p_0) = 0 for all v_I */
1405   if (pcbddc->dbg_flag && zerodiag) {
1406     PC_IS          *pcis = (PC_IS*)(pc->data);
1407     IS             dirIS = NULL;
1408     PetscScalar    *vals;
1409     const PetscInt *idxs;
1410     PetscInt       i,nz;
1411 
1412     /* p0 */
1413     ierr = VecSet(pcis->vec1_N,0.);CHKERRQ(ierr);
1414     ierr = PetscMalloc1(pcis->n,&vals);CHKERRQ(ierr);
1415     ierr = ISGetLocalSize(zerodiag,&nz);CHKERRQ(ierr);
1416     ierr = ISGetIndices(zerodiag,&idxs);CHKERRQ(ierr);
1417     for (i=0;i<nz;i++) vals[i] = 1.;
1418     ierr = VecSetValues(pcis->vec1_N,nz,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1419     ierr = VecAssemblyBegin(pcis->vec1_N);CHKERRQ(ierr);
1420     ierr = VecAssemblyEnd(pcis->vec1_N);CHKERRQ(ierr);
1421     /* v_I */
1422     ierr = VecSetRandom(pcis->vec2_N,NULL);CHKERRQ(ierr);
1423     for (i=0;i<nz;i++) vals[i] = 0.;
1424     ierr = VecSetValues(pcis->vec2_N,nz,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1425     ierr = ISRestoreIndices(zerodiag,&idxs);CHKERRQ(ierr);
1426     ierr = ISGetIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
1427     for (i=0;i<pcis->n_B;i++) vals[i] = 0.;
1428     ierr = VecSetValues(pcis->vec2_N,pcis->n_B,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1429     ierr = ISRestoreIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
1430     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
1431     if (dirIS) {
1432       PetscInt n;
1433 
1434       ierr = ISGetLocalSize(dirIS,&n);CHKERRQ(ierr);
1435       ierr = ISGetIndices(dirIS,&idxs);CHKERRQ(ierr);
1436       for (i=0;i<n;i++) vals[i] = 0.;
1437       ierr = VecSetValues(pcis->vec2_N,n,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1438       ierr = ISRestoreIndices(dirIS,&idxs);CHKERRQ(ierr);
1439     }
1440     ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
1441     ierr = VecAssemblyBegin(pcis->vec2_N);CHKERRQ(ierr);
1442     ierr = VecAssemblyEnd(pcis->vec2_N);CHKERRQ(ierr);
1443     ierr = VecSet(matis->x,0.);CHKERRQ(ierr);
1444     ierr = MatMult(matis->A,pcis->vec1_N,matis->x);CHKERRQ(ierr);
1445     ierr = VecDot(matis->x,pcis->vec2_N,&vals[0]);CHKERRQ(ierr);
1446     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] check original mat: %1.4e\n",PetscGlobalRank,vals[0]);CHKERRQ(ierr);
1447     ierr = PetscFree(vals);CHKERRQ(ierr);
1448   }
1449 
1450   /* check PCBDDCBenignPopOrPush */
1451   if (pcbddc->dbg_flag && pcbddc->benign_saddle_point) {
1452     PC_IS  *pcis = (PC_IS*)(pc->data);
1453 
1454     ierr = VecSetRandom(pcis->vec1_global,NULL);CHKERRQ(ierr);
1455     pcbddc->benign_p0 = -PetscGlobalRank;
1456     ierr = PCBDDCBenignPopOrPushP0(pc,pcis->vec1_global,PETSC_FALSE);CHKERRQ(ierr);
1457     pcbddc->benign_p0 = 1;
1458     ierr = PCBDDCBenignPopOrPushP0(pc,pcis->vec1_global,PETSC_TRUE);CHKERRQ(ierr);
1459     if (pcbddc->benign_p0_gidx >=0 && pcbddc->benign_p0 != -PetscGlobalRank) {
1460       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error testing PCBDDCBenignPopOrPushP0! Found %1.4e instead of %1.4e\n",pcbddc->benign_p0,-PetscGlobalRank);CHKERRQ(ierr);
1461     }
1462   }
1463 
1464   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1465   if (computesolvers) {
1466     PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs;
1467 
1468     if (computesubschurs && computetopography) {
1469       ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr);
1470     }
1471     if (sub_schurs->use_mumps) {
1472       if (computesubschurs) {
1473         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1474       }
1475       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1476     } else {
1477       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1478       if (computesubschurs) {
1479         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1480       }
1481     }
1482     if (pcbddc->adaptive_selection) {
1483       ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr);
1484       computeconstraintsmatrix = PETSC_TRUE;
1485     }
1486   }
1487 
1488   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1489   new_nearnullspace_provided = PETSC_FALSE;
1490   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1491   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1492     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1493       new_nearnullspace_provided = PETSC_TRUE;
1494     } else {
1495       /* determine if the two nullspaces are different (should be lightweight) */
1496       if (nearnullspace != pcbddc->onearnullspace) {
1497         new_nearnullspace_provided = PETSC_TRUE;
1498       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1499         PetscInt         i;
1500         const Vec        *nearnullvecs;
1501         PetscObjectState state;
1502         PetscInt         nnsp_size;
1503         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1504         for (i=0;i<nnsp_size;i++) {
1505           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1506           if (pcbddc->onearnullvecs_state[i] != state) {
1507             new_nearnullspace_provided = PETSC_TRUE;
1508             break;
1509           }
1510         }
1511       }
1512     }
1513   } else {
1514     if (!nearnullspace) { /* both nearnullspaces are null */
1515       new_nearnullspace_provided = PETSC_FALSE;
1516     } else { /* nearnullspace attached later */
1517       new_nearnullspace_provided = PETSC_TRUE;
1518     }
1519   }
1520 
1521   /* Setup constraints and related work vectors */
1522   /* reset primal space flags */
1523   pcbddc->new_primal_space = PETSC_FALSE;
1524   pcbddc->new_primal_space_local = PETSC_FALSE;
1525   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1526     /* It also sets the primal space flags */
1527     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1528     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1529     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
1530   }
1531 
1532   if (computesolvers || pcbddc->new_primal_space) {
1533     if (pcbddc->use_change_of_basis) {
1534       PC_IS *pcis = (PC_IS*)(pc->data);
1535       Mat   temp_mat = NULL;
1536 
1537       if (zerodiag) {
1538         /* insert B0 in pcbddc->local_mat */
1539         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_FALSE);CHKERRQ(ierr);
1540         if (pcbddc->dbg_flag) {
1541           PC_IS          *pcis = (PC_IS*)(pc->data);
1542           IS             dirIS = NULL;
1543           PetscScalar    *vals;
1544           const PetscInt *idxs;
1545           PetscInt       i,nz;
1546 
1547           /* p0 */
1548           ierr = VecSet(pcis->vec1_N,0.);CHKERRQ(ierr);
1549           ierr = PetscMalloc1(pcis->n,&vals);CHKERRQ(ierr);
1550           ierr = ISGetLocalSize(zerodiag,&nz);CHKERRQ(ierr);
1551           ierr = ISGetIndices(zerodiag,&idxs);CHKERRQ(ierr);
1552           vals[0] = 1.;
1553           ierr = VecSetValues(pcis->vec1_N,1,&idxs[nz-1],vals,INSERT_VALUES);CHKERRQ(ierr);
1554           ierr = VecAssemblyBegin(pcis->vec1_N);CHKERRQ(ierr);
1555           ierr = VecAssemblyEnd(pcis->vec1_N);CHKERRQ(ierr);
1556           /* v_I */
1557           ierr = VecSetRandom(pcis->vec2_N,NULL);CHKERRQ(ierr);
1558           for (i=0;i<nz;i++) vals[i] = 0.;
1559           ierr = VecSetValues(pcis->vec2_N,nz,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1560           ierr = ISRestoreIndices(zerodiag,&idxs);CHKERRQ(ierr);
1561           ierr = ISGetIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
1562           for (i=0;i<pcis->n_B;i++) vals[i] = 0.;
1563           ierr = VecSetValues(pcis->vec2_N,pcis->n_B,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1564           ierr = ISRestoreIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
1565           ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
1566           if (dirIS) {
1567             PetscInt n;
1568 
1569             ierr = ISGetLocalSize(dirIS,&n);CHKERRQ(ierr);
1570             ierr = ISGetIndices(dirIS,&idxs);CHKERRQ(ierr);
1571             for (i=0;i<n;i++) vals[i] = 0.;
1572             ierr = VecSetValues(pcis->vec2_N,n,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
1573             ierr = ISRestoreIndices(dirIS,&idxs);CHKERRQ(ierr);
1574           }
1575           ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
1576           ierr = VecAssemblyBegin(pcis->vec2_N);CHKERRQ(ierr);
1577           ierr = VecAssemblyEnd(pcis->vec2_N);CHKERRQ(ierr);
1578           ierr = VecSet(matis->x,0.);CHKERRQ(ierr);
1579           ierr = MatMult(pcbddc->local_mat,pcis->vec1_N,matis->x);CHKERRQ(ierr);
1580           ierr = VecDot(matis->x,pcis->vec2_N,&vals[0]);CHKERRQ(ierr);
1581           ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] check new mat: %1.4e\n",PetscGlobalRank,vals[0]);CHKERRQ(ierr);
1582           ierr = PetscFree(vals);CHKERRQ(ierr);
1583         }
1584         /* hack: swap pointers */
1585         temp_mat = matis->A;
1586         matis->A = pcbddc->local_mat;
1587         pcbddc->local_mat = NULL;
1588       }
1589       ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
1590       if (zerodiag) {
1591         /* restore original matrix */
1592         ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
1593         matis->A = temp_mat;
1594         /* pop B0 from pcbddc->local_mat */
1595         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
1596       }
1597       /* get submatrices */
1598       ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
1599       ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
1600       ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
1601       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
1602       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
1603       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
1604       /* set flag in pcis to not reuse submatrices during PCISCreate */
1605       pcis->reusesubmatrices = PETSC_FALSE;
1606     } else if (!pcbddc->user_ChangeOfBasisMatrix) {
1607       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1608       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1609       pcbddc->local_mat = matis->A;
1610     }
1611     /* SetUp coarse and local Neumann solvers */
1612     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1613     /* SetUp Scaling operator */
1614     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
1615   }
1616 
1617   if (pcbddc->dbg_flag) {
1618     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1619   }
1620   PetscFunctionReturn(0);
1621 }
1622 
1623 /* -------------------------------------------------------------------------- */
1624 /*
1625    PCApply_BDDC - Applies the BDDC operator to a vector.
1626 
1627    Input Parameters:
1628 +  pc - the preconditioner context
1629 -  r - input vector (global)
1630 
1631    Output Parameter:
1632 .  z - output vector (global)
1633 
1634    Application Interface Routine: PCApply()
1635  */
1636 #undef __FUNCT__
1637 #define __FUNCT__ "PCApply_BDDC"
1638 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1639 {
1640   PC_IS             *pcis = (PC_IS*)(pc->data);
1641   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1642   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
1643   PetscErrorCode    ierr;
1644   const PetscScalar one = 1.0;
1645   const PetscScalar m_one = -1.0;
1646   const PetscScalar zero = 0.0;
1647 
1648 /* This code is similar to that provided in nn.c for PCNN
1649    NN interface preconditioner changed to BDDC
1650    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
1651 
1652   PetscFunctionBegin;
1653   if (pcbddc->benign_saddle_point) { /* extract p0 from r */
1654     ierr = PCBDDCBenignPopOrPushP0(pc,r,PETSC_TRUE);CHKERRQ(ierr);
1655   }
1656   if (!pcbddc->use_exact_dirichlet_trick) {
1657     ierr = VecCopy(r,z);CHKERRQ(ierr);
1658     /* First Dirichlet solve */
1659     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1660     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1661     /*
1662       Assembling right hand side for BDDC operator
1663       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1664       - pcis->vec1_B the interface part of the global vector z
1665     */
1666     if (n_D) {
1667       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1668       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1669       if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1670       ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1671     } else {
1672       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1673     }
1674     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1675     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1676     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1677   } else {
1678     if (pcbddc->switch_static) {
1679       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1680     }
1681     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1682   }
1683 
1684   /* Apply interface preconditioner
1685      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1686   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
1687 
1688   /* Apply transpose of partition of unity operator */
1689   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1690 
1691   /* Second Dirichlet solve and assembling of output */
1692   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1693   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1694   if (n_B) {
1695     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1696     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1697   } else if (pcbddc->switch_static) {
1698     ierr = MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1699   }
1700   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1701 
1702   if (!pcbddc->use_exact_dirichlet_trick) {
1703     if (pcbddc->switch_static) {
1704       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1705     } else {
1706       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1707     }
1708     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1709     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1710   } else {
1711     if (pcbddc->switch_static) {
1712       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1713     } else {
1714       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1715     }
1716     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1717     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1718   }
1719 
1720   if (pcbddc->benign_saddle_point) { /* push p0 (computed in PCBDDCApplyInterface) */
1721     ierr = PCBDDCBenignPopOrPushP0(pc,z,PETSC_FALSE);CHKERRQ(ierr);
1722   }
1723   PetscFunctionReturn(0);
1724 }
1725 
1726 /* -------------------------------------------------------------------------- */
1727 /*
1728    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
1729 
1730    Input Parameters:
1731 +  pc - the preconditioner context
1732 -  r - input vector (global)
1733 
1734    Output Parameter:
1735 .  z - output vector (global)
1736 
1737    Application Interface Routine: PCApplyTranspose()
1738  */
1739 #undef __FUNCT__
1740 #define __FUNCT__ "PCApplyTranspose_BDDC"
1741 PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
1742 {
1743   PC_IS             *pcis = (PC_IS*)(pc->data);
1744   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1745   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
1746   PetscErrorCode    ierr;
1747   const PetscScalar one = 1.0;
1748   const PetscScalar m_one = -1.0;
1749   const PetscScalar zero = 0.0;
1750 
1751   PetscFunctionBegin;
1752   if (!pcbddc->use_exact_dirichlet_trick) {
1753     ierr = VecCopy(r,z);CHKERRQ(ierr);
1754     /* First Dirichlet solve */
1755     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1756     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1757     /*
1758       Assembling right hand side for BDDC operator
1759       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1760       - pcis->vec1_B the interface part of the global vector z
1761     */
1762     if (n_D) {
1763       ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1764       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1765       if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1766       ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1767     } else {
1768       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1769     }
1770     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1771     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1772     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1773   } else {
1774     if (pcbddc->switch_static) {
1775       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1776     }
1777     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1778   }
1779 
1780   /* Apply interface preconditioner
1781      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1782   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
1783 
1784   /* Apply transpose of partition of unity operator */
1785   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1786 
1787   /* Second Dirichlet solve and assembling of output */
1788   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1789   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1790   if (n_B) {
1791     ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1792     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1793   } else if (pcbddc->switch_static) {
1794     ierr = MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1795   }
1796   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1797   if (!pcbddc->use_exact_dirichlet_trick) {
1798     if (pcbddc->switch_static) {
1799       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1800     } else {
1801       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1802     }
1803     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1804     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1805   } else {
1806     if (pcbddc->switch_static) {
1807       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1808     } else {
1809       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1810     }
1811     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1812     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1813   }
1814   PetscFunctionReturn(0);
1815 }
1816 /* -------------------------------------------------------------------------- */
1817 
1818 #undef __FUNCT__
1819 #define __FUNCT__ "PCDestroy_BDDC"
1820 PetscErrorCode PCDestroy_BDDC(PC pc)
1821 {
1822   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1823   PetscErrorCode ierr;
1824 
1825   PetscFunctionBegin;
1826   /* free data created by PCIS */
1827   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1828   /* free BDDC custom data  */
1829   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1830   /* destroy objects related to topography */
1831   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1832   /* free allocated graph structure */
1833   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1834   /* free allocated sub schurs structure */
1835   ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr);
1836   /* destroy objects for scaling operator */
1837   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1838   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
1839   /* free solvers stuff */
1840   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1841   /* free global vectors needed in presolve */
1842   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
1843   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1844   /* free stuff for change of basis hooks */
1845   if (pcbddc->new_global_mat) {
1846     PCBDDCChange_ctx change_ctx;
1847     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1848     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1849     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1850     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1851     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1852   }
1853   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
1854   /* remove functions */
1855   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1856   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1857   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
1858   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1859   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
1860   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1861   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1862   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1863   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1864   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1865   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1866   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1867   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1868   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1869   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1870   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1871   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1872   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1873   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1874   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1875   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1876   /* Free the private data structure */
1877   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1878   PetscFunctionReturn(0);
1879 }
1880 /* -------------------------------------------------------------------------- */
1881 
1882 #undef __FUNCT__
1883 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
1884 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1885 {
1886   FETIDPMat_ctx  mat_ctx;
1887   Vec            copy_standard_rhs;
1888   PC_IS*         pcis;
1889   PC_BDDC*       pcbddc;
1890   PetscErrorCode ierr;
1891 
1892   PetscFunctionBegin;
1893   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1894   pcis = (PC_IS*)mat_ctx->pc->data;
1895   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1896 
1897   /*
1898      change of basis for physical rhs if needed
1899      It also changes the rhs in case of dirichlet boundaries
1900      TODO: better management when FETIDP will have its own class
1901   */
1902   ierr = VecDuplicate(standard_rhs,&copy_standard_rhs);CHKERRQ(ierr);
1903   ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr);
1904   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr);
1905   /* store vectors for computation of fetidp final solution */
1906   ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1907   ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1908   /* scale rhs since it should be unassembled */
1909   /* TODO use counter scaling? (also below) */
1910   ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1911   ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1912   /* Apply partition of unity */
1913   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1914   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1915   if (!pcbddc->switch_static) {
1916     /* compute partially subassembled Schur complement right-hand side */
1917     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1918     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
1919     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1920     ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr);
1921     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1922     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1923     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1924     ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1925     ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1926     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1927   }
1928   ierr = VecDestroy(&copy_standard_rhs);CHKERRQ(ierr);
1929   /* BDDC rhs */
1930   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
1931   if (pcbddc->switch_static) {
1932     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1933   }
1934   /* apply BDDC */
1935   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
1936   /* Application of B_delta and assembling of rhs for fetidp fluxes */
1937   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
1938   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
1939   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1940   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1941   PetscFunctionReturn(0);
1942 }
1943 
1944 #undef __FUNCT__
1945 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
1946 /*@
1947  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side
1948 
1949    Collective
1950 
1951    Input Parameters:
1952 +  fetidp_mat      - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
1953 -  standard_rhs    - the right-hand side of the original linear system
1954 
1955    Output Parameters:
1956 .  fetidp_flux_rhs - the right-hand side for the FETI-DP linear system
1957 
1958    Level: developer
1959 
1960    Notes:
1961 
1962 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
1963 @*/
1964 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1965 {
1966   FETIDPMat_ctx  mat_ctx;
1967   PetscErrorCode ierr;
1968 
1969   PetscFunctionBegin;
1970   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1971   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
1972   PetscFunctionReturn(0);
1973 }
1974 /* -------------------------------------------------------------------------- */
1975 
1976 #undef __FUNCT__
1977 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
1978 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1979 {
1980   FETIDPMat_ctx  mat_ctx;
1981   PC_IS*         pcis;
1982   PC_BDDC*       pcbddc;
1983   PetscErrorCode ierr;
1984 
1985   PetscFunctionBegin;
1986   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1987   pcis = (PC_IS*)mat_ctx->pc->data;
1988   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1989 
1990   /* apply B_delta^T */
1991   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1992   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1993   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
1994   /* compute rhs for BDDC application */
1995   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1996   if (pcbddc->switch_static) {
1997     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1998   }
1999   /* apply BDDC */
2000   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
2001   /* put values into standard global vector */
2002   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2003   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2004   if (!pcbddc->switch_static) {
2005     /* compute values into the interior if solved for the partially subassembled Schur complement */
2006     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
2007     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
2008     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
2009   }
2010   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2011   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2012   /* final change of basis if needed
2013      Is also sums the dirichlet part removed during RHS assembling */
2014   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
2015   PetscFunctionReturn(0);
2016 }
2017 
2018 #undef __FUNCT__
2019 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
2020 /*@
2021  PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system
2022 
2023    Collective
2024 
2025    Input Parameters:
2026 +  fetidp_mat      - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
2027 -  fetidp_flux_sol - the solution of the FETI-DP linear system
2028 
2029    Output Parameters:
2030 .  standard_sol    - the solution defined on the physical domain
2031 
2032    Level: developer
2033 
2034    Notes:
2035 
2036 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
2037 @*/
2038 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
2039 {
2040   FETIDPMat_ctx  mat_ctx;
2041   PetscErrorCode ierr;
2042 
2043   PetscFunctionBegin;
2044   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
2045   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
2046   PetscFunctionReturn(0);
2047 }
2048 /* -------------------------------------------------------------------------- */
2049 
2050 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
2051 extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
2052 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
2053 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
2054 extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
2055 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
2056 
2057 #undef __FUNCT__
2058 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
2059 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
2060 {
2061 
2062   FETIDPMat_ctx  fetidpmat_ctx;
2063   Mat            newmat;
2064   FETIDPPC_ctx   fetidppc_ctx;
2065   PC             newpc;
2066   MPI_Comm       comm;
2067   PetscErrorCode ierr;
2068 
2069   PetscFunctionBegin;
2070   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
2071   /* FETIDP linear matrix */
2072   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
2073   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
2074   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
2075   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
2076   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
2077   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
2078   ierr = MatSetUp(newmat);CHKERRQ(ierr);
2079   /* FETIDP preconditioner */
2080   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
2081   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
2082   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
2083   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
2084   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
2085   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
2086   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
2087   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
2088   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
2089   ierr = PCSetUp(newpc);CHKERRQ(ierr);
2090   /* return pointers for objects created */
2091   *fetidp_mat=newmat;
2092   *fetidp_pc=newpc;
2093   PetscFunctionReturn(0);
2094 }
2095 
2096 #undef __FUNCT__
2097 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
2098 /*@
2099  PCBDDCCreateFETIDPOperators - Create FETI-DP operators
2100 
2101    Collective
2102 
2103    Input Parameters:
2104 .  pc - the BDDC preconditioning context (setup should have been called before)
2105 
2106    Output Parameters:
2107 +  fetidp_mat - shell FETI-DP matrix object
2108 -  fetidp_pc  - shell Dirichlet preconditioner for FETI-DP matrix
2109 
2110    Options Database Keys:
2111 .    -fetidp_fullyredundant <false> - use or not a fully redundant set of Lagrange multipliers
2112 
2113    Level: developer
2114 
2115    Notes:
2116      Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose
2117 
2118 .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
2119 @*/
2120 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
2121 {
2122   PetscErrorCode ierr;
2123 
2124   PetscFunctionBegin;
2125   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2126   if (pc->setupcalled) {
2127     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
2128   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
2129   PetscFunctionReturn(0);
2130 }
2131 /* -------------------------------------------------------------------------- */
2132 /*MC
2133    PCBDDC - Balancing Domain Decomposition by Constraints.
2134 
2135    An implementation of the BDDC preconditioner based on
2136 
2137 .vb
2138    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
2139    [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
2140    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
2141    [4] C. Pechstein and C. R. Dohrmann. "Modern domain decomposition methods BDDC, deluxe scaling, and an algebraic approach", Seminar talk, Linz, December 2013, http://people.ricam.oeaw.ac.at/c.pechstein/pechstein-bddc2013.pdf
2142 .ve
2143 
2144    The matrix to be preconditioned (Pmat) must be of type MATIS.
2145 
2146    Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
2147 
2148    It also works with unsymmetric and indefinite problems.
2149 
2150    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.
2151 
2152    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace()
2153 
2154    Boundary nodes are split in vertices, edges and faces classes 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()
2155    Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesLocalIS()
2156 
2157    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.
2158 
2159    Change of basis is performed similarly to [2] when requested. When more than 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.
2160    User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()
2161 
2162    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.
2163 
2164    Adaptive selection of primal constraints [4] is supported for SPD systems with high-contrast in the coefficients if MUMPS is present. Future versions of the code will also consider using MKL_PARDISO or PASTIX.
2165 
2166    An experimental interface to the FETI-DP method is available. FETI-DP operators could be created using PCBDDCCreateFETIDPOperators(). A stand-alone class for the FETI-DP method will be provided in the next releases.
2167    Deluxe scaling is not supported yet for FETI-DP.
2168 
2169    Options Database Keys (some of them, run with -h for a complete list):
2170 
2171 .    -pc_bddc_use_vertices <true> - use or not vertices in primal space
2172 .    -pc_bddc_use_edges <true> - use or not edges in primal space
2173 .    -pc_bddc_use_faces <false> - use or not faces in primal space
2174 .    -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
2175 .    -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
2176 .    -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
2177 .    -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
2178 .    -pc_bddc_levels <0> - maximum number of levels for multilevel
2179 .    -pc_bddc_coarsening_ratio <8> - number of subdomains which will be aggregated together at the coarser level (e.g. H/h ratio at the coarser level, significative only in the multilevel case)
2180 .    -pc_bddc_redistribute <0> - size of a subset of processors where the coarse problem will be remapped (the value is ignored if not at the coarsest level)
2181 .    -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
2182 .    -pc_bddc_schur_layers <-1> - select the economic version of deluxe scaling by specifying the number of layers (-1 corresponds to the original deluxe scaling)
2183 .    -pc_bddc_adaptive_threshold <0.0> - when a value greater than one is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS installed)
2184 -    -pc_bddc_check_level <0> - set verbosity level of debugging output
2185 
2186    Options for Dirichlet, Neumann or coarse solver can be set with
2187 .vb
2188       -pc_bddc_dirichlet_
2189       -pc_bddc_neumann_
2190       -pc_bddc_coarse_
2191 .ve
2192    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.
2193 
2194    When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
2195 .vb
2196       -pc_bddc_dirichlet_lN_
2197       -pc_bddc_neumann_lN_
2198       -pc_bddc_coarse_lN_
2199 .ve
2200    Note that level number ranges from the finest (0) to the coarsest (N).
2201    In order to specify options for the BDDC operators at the coarser levels (and not for the solvers), prepend -pc_bddc_coarse_ or -pc_bddc_coarse_l to the option, e.g.
2202 .vb
2203      -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
2204 .ve
2205    will use a threshold of 5 for constraints' selection at the first coarse level and will redistribute the coarse problem of the first coarse level on 3 processors
2206 
2207    Level: intermediate
2208 
2209    Developer notes:
2210 
2211    Contributed by Stefano Zampini
2212 
2213 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
2214 M*/
2215 
2216 #undef __FUNCT__
2217 #define __FUNCT__ "PCCreate_BDDC"
2218 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
2219 {
2220   PetscErrorCode      ierr;
2221   PC_BDDC             *pcbddc;
2222 
2223   PetscFunctionBegin;
2224   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
2225   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
2226   pc->data  = (void*)pcbddc;
2227 
2228   /* create PCIS data structure */
2229   ierr = PCISCreate(pc);CHKERRQ(ierr);
2230 
2231   /* BDDC customization */
2232   pcbddc->use_local_adj       = PETSC_TRUE;
2233   pcbddc->use_vertices        = PETSC_TRUE;
2234   pcbddc->use_edges           = PETSC_TRUE;
2235   pcbddc->use_faces           = PETSC_FALSE;
2236   pcbddc->use_change_of_basis = PETSC_FALSE;
2237   pcbddc->use_change_on_faces = PETSC_FALSE;
2238   pcbddc->switch_static       = PETSC_FALSE;
2239   pcbddc->use_nnsp_true       = PETSC_FALSE;
2240   pcbddc->use_qr_single       = PETSC_FALSE;
2241   pcbddc->symmetric_primal    = PETSC_TRUE;
2242   pcbddc->benign_saddle_point = PETSC_FALSE;
2243   pcbddc->dbg_flag            = 0;
2244   /* private */
2245   pcbddc->local_primal_size          = 0;
2246   pcbddc->local_primal_size_cc       = 0;
2247   pcbddc->local_primal_ref_node      = 0;
2248   pcbddc->local_primal_ref_mult      = 0;
2249   pcbddc->n_vertices                 = 0;
2250   pcbddc->primal_indices_local_idxs  = 0;
2251   pcbddc->recompute_topography       = PETSC_FALSE;
2252   pcbddc->coarse_size                = -1;
2253   pcbddc->new_primal_space           = PETSC_FALSE;
2254   pcbddc->new_primal_space_local     = PETSC_FALSE;
2255   pcbddc->global_primal_indices      = 0;
2256   pcbddc->onearnullspace             = 0;
2257   pcbddc->onearnullvecs_state        = 0;
2258   pcbddc->user_primal_vertices       = 0;
2259   pcbddc->NullSpace                  = 0;
2260   pcbddc->temp_solution              = 0;
2261   pcbddc->original_rhs               = 0;
2262   pcbddc->local_mat                  = 0;
2263   pcbddc->ChangeOfBasisMatrix        = 0;
2264   pcbddc->user_ChangeOfBasisMatrix   = 0;
2265   pcbddc->new_global_mat             = 0;
2266   pcbddc->coarse_vec                 = 0;
2267   pcbddc->coarse_ksp                 = 0;
2268   pcbddc->coarse_phi_B               = 0;
2269   pcbddc->coarse_phi_D               = 0;
2270   pcbddc->coarse_psi_B               = 0;
2271   pcbddc->coarse_psi_D               = 0;
2272   pcbddc->vec1_P                     = 0;
2273   pcbddc->vec1_R                     = 0;
2274   pcbddc->vec2_R                     = 0;
2275   pcbddc->local_auxmat1              = 0;
2276   pcbddc->local_auxmat2              = 0;
2277   pcbddc->R_to_B                     = 0;
2278   pcbddc->R_to_D                     = 0;
2279   pcbddc->ksp_D                      = 0;
2280   pcbddc->ksp_R                      = 0;
2281   pcbddc->NeumannBoundaries          = 0;
2282   pcbddc->NeumannBoundariesLocal     = 0;
2283   pcbddc->DirichletBoundaries        = 0;
2284   pcbddc->DirichletBoundariesLocal   = 0;
2285   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
2286   pcbddc->n_ISForDofs                = 0;
2287   pcbddc->n_ISForDofsLocal           = 0;
2288   pcbddc->ISForDofs                  = 0;
2289   pcbddc->ISForDofsLocal             = 0;
2290   pcbddc->ConstraintMatrix           = 0;
2291   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
2292   pcbddc->coarse_loc_to_glob         = 0;
2293   pcbddc->coarsening_ratio           = 8;
2294   pcbddc->coarse_adj_red             = 0;
2295   pcbddc->current_level              = 0;
2296   pcbddc->max_levels                 = 0;
2297   pcbddc->use_coarse_estimates       = PETSC_FALSE;
2298   pcbddc->redistribute_coarse        = 0;
2299   pcbddc->coarse_subassembling       = 0;
2300   pcbddc->coarse_subassembling_init  = 0;
2301 
2302   /* benign subspace trick */
2303   pcbddc->B0_ncol                    = 0;
2304   pcbddc->B0_cols                    = NULL;
2305   pcbddc->B0_vals                    = NULL;
2306   pcbddc->benign_change              = 0;
2307 
2308   /* create local graph structure */
2309   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
2310 
2311   /* scaling */
2312   pcbddc->work_scaling          = 0;
2313   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
2314   pcbddc->faster_deluxe         = PETSC_FALSE;
2315 
2316   /* create sub schurs structure */
2317   ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr);
2318   pcbddc->sub_schurs_rebuild     = PETSC_FALSE;
2319   pcbddc->sub_schurs_layers      = -1;
2320   pcbddc->sub_schurs_use_useradj = PETSC_FALSE;
2321 
2322   pcbddc->computed_rowadj = PETSC_FALSE;
2323 
2324   /* adaptivity */
2325   pcbddc->adaptive_threshold      = 0.0;
2326   pcbddc->adaptive_nmax           = 0;
2327   pcbddc->adaptive_nmin           = 0;
2328 
2329   /* function pointers */
2330   pc->ops->apply               = PCApply_BDDC;
2331   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2332   pc->ops->setup               = PCSetUp_BDDC;
2333   pc->ops->destroy             = PCDestroy_BDDC;
2334   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
2335   pc->ops->view                = 0;
2336   pc->ops->applyrichardson     = 0;
2337   pc->ops->applysymmetricleft  = 0;
2338   pc->ops->applysymmetricright = 0;
2339   pc->ops->presolve            = PCPreSolve_BDDC;
2340   pc->ops->postsolve           = PCPostSolve_BDDC;
2341 
2342   /* composing function */
2343   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
2344   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
2345   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
2346   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
2347   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
2348   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
2349   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
2350   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
2351   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2352   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
2353   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2354   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
2355   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2356   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
2357   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2358   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
2359   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
2360   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
2361   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
2362   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
2363   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
2364   PetscFunctionReturn(0);
2365 }
2366 
2367