xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 06a4e24a4ea34074cd1dc021a3bf2cd6abf1fe12)
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   /* store partially computed solution and set initial guess */
1039   if (x && pcbddc->use_exact_dirichlet_trick) {
1040     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1041     ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1042     ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1043     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1044     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1045     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1046     if (ksp) {
1047       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1048     }
1049   }
1050 
1051   if (pcbddc->ChangeOfBasisMatrix) {
1052     PCBDDCChange_ctx change_ctx;
1053 
1054     /* get change ctx */
1055     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1056 
1057     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1058     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1059     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1060     change_ctx->original_mat = pc->mat;
1061 
1062     /* change iteration matrix */
1063     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1064     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1065     pc->mat = pcbddc->new_global_mat;
1066 
1067     /* store the original rhs */
1068     if (copy_rhs) {
1069       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1070       copy_rhs = PETSC_FALSE;
1071     }
1072 
1073     /* change rhs */
1074     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1075 
1076     /* change rhs on pressures (iteration matrix is of type MATIS) */
1077     if (pcbddc->benign_saddle_point) {
1078       Mat_IS *matis = (Mat_IS*)(change_ctx->original_mat->data);
1079 
1080       ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1081       ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1082       if (pcbddc->benign_change) {
1083         ierr = MatMultTranspose(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
1084         ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1085         ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1086         /* change local iteration matrix */
1087         ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_REUSE_MATRIX,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
1088         ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr);
1089         ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1090         pcbddc->benign_original_mat = matis->A;
1091         ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
1092         ierr = PetscObjectReference((PetscObject)pcbddc->local_mat);CHKERRQ(ierr);
1093         matis->A = pcbddc->local_mat;
1094         ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1095       } else {
1096         ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1097         ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1098       }
1099     }
1100 
1101     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
1102     pcbddc->rhs_change = PETSC_TRUE;
1103   }
1104 
1105   /* remove nullspace if present */
1106   if (ksp && x && pcbddc->NullSpace) {
1107     ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr);
1108     /* store the original rhs */
1109     if (copy_rhs) {
1110       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1111       copy_rhs = PETSC_FALSE;
1112     }
1113     pcbddc->rhs_change = PETSC_TRUE;
1114     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1115   }
1116   PetscFunctionReturn(0);
1117 }
1118 
1119 /* -------------------------------------------------------------------------- */
1120 #undef __FUNCT__
1121 #define __FUNCT__ "PCPostSolve_BDDC"
1122 /* -------------------------------------------------------------------------- */
1123 /*
1124    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1125                      approach has been selected. Also, restores rhs to its original state.
1126 
1127    Input Parameter:
1128 +  pc - the preconditioner contex
1129 
1130    Application Interface Routine: PCPostSolve()
1131 
1132    Notes:
1133      The interface routine PCPostSolve() is not usually called directly by
1134      the user, but instead is called by KSPSolve().
1135 */
1136 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1137 {
1138   PetscErrorCode ierr;
1139   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1140 
1141   PetscFunctionBegin;
1142   if (pcbddc->ChangeOfBasisMatrix) {
1143     PCBDDCChange_ctx change_ctx;
1144 
1145     /* get change ctx */
1146     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1147 
1148     /* restore iteration matrix */
1149     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1150     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1151     pc->mat = change_ctx->original_mat;
1152 
1153     /* get solution in original basis */
1154     if (x) {
1155       PC_IS *pcis = (PC_IS*)(pc->data);
1156       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1157       /* restore solution on pressures */
1158       if (pcbddc->benign_saddle_point) {
1159         Mat_IS *matis = (Mat_IS*)pc->mat->data;
1160 
1161         ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1162         ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1163         if (pcbddc->benign_change) {
1164 
1165           ierr = MatMult(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
1166           ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1167           ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1168           /* restore local iteration matrix */
1169           ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
1170           ierr = PetscObjectReference((PetscObject)pcbddc->benign_original_mat);CHKERRQ(ierr);
1171           matis->A = pcbddc->benign_original_mat;
1172           ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr);
1173         } else {
1174           ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1175           ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1176         }
1177       }
1178       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
1179     }
1180   }
1181 
1182   /* add solution removed in presolve */
1183   if (x && pcbddc->rhs_change) {
1184     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
1185   }
1186 
1187   /* restore rhs to its original state */
1188   if (rhs && pcbddc->rhs_change) {
1189     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1190   }
1191   pcbddc->rhs_change = PETSC_FALSE;
1192 
1193   /* restore ksp guess state */
1194   if (ksp) {
1195     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
1196   }
1197   PetscFunctionReturn(0);
1198 }
1199 /* -------------------------------------------------------------------------- */
1200 #undef __FUNCT__
1201 #define __FUNCT__ "PCSetUp_BDDC"
1202 /* -------------------------------------------------------------------------- */
1203 /*
1204    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1205                   by setting data structures and options.
1206 
1207    Input Parameter:
1208 +  pc - the preconditioner context
1209 
1210    Application Interface Routine: PCSetUp()
1211 
1212    Notes:
1213      The interface routine PCSetUp() is not usually called directly by
1214      the user, but instead is called by PCApply() if necessary.
1215 */
1216 PetscErrorCode PCSetUp_BDDC(PC pc)
1217 {
1218   PetscErrorCode ierr;
1219   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
1220   Mat_IS*        matis;
1221   MatNullSpace   nearnullspace;
1222   PetscInt       nrows,ncols;
1223   PetscBool      computetopography,computesolvers,computesubschurs;
1224   PetscBool      computeconstraintsmatrix;
1225   PetscBool      new_nearnullspace_provided,ismatis;
1226 
1227   PetscFunctionBegin;
1228   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr);
1229   if (!ismatis) {
1230     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
1231   }
1232   ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr);
1233   if (nrows != ncols) {
1234     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
1235   }
1236   matis = (Mat_IS*)pc->pmat->data;
1237   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1238   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1239      Also, BDDC directly build the Dirichlet problem */
1240   /* split work */
1241   if (pc->setupcalled) {
1242     if (pc->flag == SAME_NONZERO_PATTERN) {
1243       computetopography = PETSC_FALSE;
1244       computesolvers = PETSC_TRUE;
1245     } else { /* DIFFERENT_NONZERO_PATTERN */
1246       computetopography = PETSC_TRUE;
1247       computesolvers = PETSC_TRUE;
1248     }
1249   } else {
1250     computetopography = PETSC_TRUE;
1251     computesolvers = PETSC_TRUE;
1252   }
1253   if (pcbddc->recompute_topography) {
1254     computetopography = PETSC_TRUE;
1255   }
1256   computeconstraintsmatrix = PETSC_FALSE;
1257   if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) {
1258     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling");
1259   }
1260   pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling);
1261   if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;
1262 
1263   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
1264   if (pcbddc->faster_deluxe && pcbddc->adaptive_selection && pcbddc->use_change_of_basis) {
1265     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");
1266   }
1267 
1268   /* check if the iteration matrix is of type MATIS in case the benign trick has been requested */
1269   ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&ismatis);CHKERRQ(ierr);
1270   if (pcbddc->benign_saddle_point && !ismatis) {
1271     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner with benign subspace trick requires the iteration matrix to be of type MATIS");
1272   }
1273 
1274   /* Get stdout for dbg */
1275   if (pcbddc->dbg_flag) {
1276     if (!pcbddc->dbg_viewer) {
1277       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1278       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1279     }
1280     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1281   }
1282 
1283   if (pcbddc->user_ChangeOfBasisMatrix) {
1284     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1285     pcbddc->use_change_of_basis = PETSC_FALSE;
1286     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
1287   } else {
1288     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1289     if (!pcbddc->benign_saddle_point) {
1290       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1291       pcbddc->local_mat = matis->A;
1292     } else { /* TODO: handle user change of basis */
1293       PetscInt  nz;
1294       PetscBool sorted;
1295 
1296       ierr = ISDestroy(&pcbddc->zerodiag);CHKERRQ(ierr);
1297       /* should I also include nonzero pressures? */
1298       ierr = MatFindZeroDiagonals(matis->A,&pcbddc->zerodiag);CHKERRQ(ierr);
1299       ierr = ISSorted(pcbddc->zerodiag,&sorted);CHKERRQ(ierr);
1300       if (!sorted) {
1301         ierr = ISSort(pcbddc->zerodiag);CHKERRQ(ierr);
1302       }
1303       //if (!PetscGlobalRank) printf("ZERODIAG\n");
1304       //if (!PetscGlobalRank) ISView(pcbddc->zerodiag,NULL);
1305       ierr = ISGetLocalSize(pcbddc->zerodiag,&nz);CHKERRQ(ierr);
1306       if (nz) {
1307         IS                zerodiagc;
1308         PetscScalar       *array;
1309         const PetscInt    *idxs,*idxsc;
1310         PetscInt          i,n,*nnz;
1311 
1312         /* TODO: add check for shared dofs */
1313         pcbddc->use_change_of_basis = PETSC_TRUE;
1314         ierr = MatGetLocalSize(matis->A,&n,NULL);CHKERRQ(ierr);
1315         ierr = ISComplement(pcbddc->zerodiag,0,n,&zerodiagc);CHKERRQ(ierr);
1316         ierr = ISGetIndices(pcbddc->zerodiag,&idxs);CHKERRQ(ierr);
1317         ierr = ISGetIndices(zerodiagc,&idxsc);CHKERRQ(ierr);
1318         /* local change of basis for pressures */
1319         ierr = MatDestroy(&pcbddc->benign_change);CHKERRQ(ierr);
1320         ierr = MatCreate(PetscObjectComm((PetscObject)matis->A),&pcbddc->benign_change);CHKERRQ(ierr);
1321         ierr = MatSetType(pcbddc->benign_change,MATAIJ);CHKERRQ(ierr);
1322         ierr = MatSetSizes(pcbddc->benign_change,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
1323         ierr = PetscMalloc1(n,&nnz);CHKERRQ(ierr);
1324         for (i=0;i<n-nz;i++) nnz[idxsc[i]] = 1; /* identity on velocities */
1325         for (i=0;i<nz-1;i++) nnz[idxs[i]] = 2; /* change on pressures */
1326         nnz[idxs[nz-1]] = nz; /* last local pressure dof: _0 set */
1327         ierr = MatSeqAIJSetPreallocation(pcbddc->benign_change,0,nnz);CHKERRQ(ierr);
1328         ierr = PetscFree(nnz);CHKERRQ(ierr);
1329         /* set identity on velocities */
1330         for (i=0;i<n-nz;i++) {
1331           ierr = MatSetValue(pcbddc->benign_change,idxsc[i],idxsc[i],1.,INSERT_VALUES);CHKERRQ(ierr);
1332         }
1333         /* set change on pressures */
1334         for (i=0;i<nz-1;i++) {
1335           PetscScalar vals[2];
1336           PetscInt    cols[2];
1337 
1338           /* TODO: add quadrature */
1339           cols[0] = idxs[i];
1340           cols[1] = idxs[nz-1];
1341           vals[0] = 1.;
1342           vals[1] = 1./nz;
1343           ierr = MatSetValues(pcbddc->benign_change,1,idxs+i,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
1344         }
1345         ierr = PetscMalloc1(nz,&array);CHKERRQ(ierr);
1346         for (i=0;i<nz-1;i++) array[i] = -1.;
1347         array[nz-1] = 1./nz;
1348         ierr = MatSetValues(pcbddc->benign_change,1,idxs+nz-1,nz,idxs,array,INSERT_VALUES);CHKERRQ(ierr);
1349         ierr = PetscFree(array);CHKERRQ(ierr);
1350         ierr = MatAssemblyBegin(pcbddc->benign_change,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1351         ierr = MatAssemblyEnd(pcbddc->benign_change,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1352         /* TODO: need optimization? */
1353         ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_INITIAL_MATRIX,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
1354         //if (!PetscGlobalRank) printf("Local MAT\n");
1355         //if (!PetscGlobalRank) MatView(matis->A,NULL);
1356         //if (!PetscGlobalRank) printf("Local pressure change\n");
1357         //if (!PetscGlobalRank) MatView(pcbddc->benign_change,NULL);
1358         //if (!PetscGlobalRank) printf("Local AFTER CHANGE\n");
1359         //if (!PetscGlobalRank) MatView(pcbddc->local_mat,NULL);
1360         ierr = ISRestoreIndices(pcbddc->zerodiag,&idxs);CHKERRQ(ierr);
1361         ierr = ISRestoreIndices(zerodiagc,&idxsc);CHKERRQ(ierr);
1362         ierr = ISDestroy(&zerodiagc);CHKERRQ(ierr);
1363         /* pop B0 mat from pcbddc->local_mat */
1364         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
1365         //if (!PetscGlobalRank) printf("REMOVED AND B0\n");
1366         //if (!PetscGlobalRank) MatView(pcbddc->local_mat,NULL);
1367         //if (!PetscGlobalRank) {
1368         //  for (i=0;i<pcbddc->B0_ncol;i++) printf("%d %f\n",pcbddc->B0_cols[i],pcbddc->B0_vals[i]);
1369         //}
1370       } else { /* this is unlikely to happen */
1371         ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1372         pcbddc->local_mat = matis->A;
1373       }
1374     }
1375   }
1376 
1377   /* workaround for reals */
1378 #if !defined(PETSC_USE_COMPLEX)
1379   if (matis->A->symmetric_set) {
1380     ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr);
1381   }
1382 #endif
1383   if (matis->A->symmetric_set) {
1384     ierr = MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);CHKERRQ(ierr);
1385   }
1386   if (matis->A->spd_set) {
1387     ierr = MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);CHKERRQ(ierr);
1388   }
1389 
1390   /* Set up all the "iterative substructuring" common block without computing solvers */
1391   {
1392     Mat temp_mat;
1393 
1394     temp_mat = matis->A;
1395     matis->A = pcbddc->local_mat;
1396     ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr);
1397     pcbddc->local_mat = matis->A;
1398     matis->A = temp_mat;
1399       //if (!PetscGlobalRank) printf("IS_I_LOCAL\n");
1400       //if (!PetscGlobalRank) { PC_IS* pcis=(PC_IS*)(pc->data); ISView(pcis->is_I_local,NULL); }
1401   }
1402 
1403   /* Analyze interface */
1404   if (computetopography) {
1405     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1406     computeconstraintsmatrix = PETSC_TRUE;
1407   }
1408 
1409   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1410   if (computesolvers) {
1411     PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs;
1412 
1413     if (computesubschurs && computetopography) {
1414       ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr);
1415     }
1416     if (sub_schurs->use_mumps) {
1417       if (computesubschurs) {
1418         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1419       }
1420       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1421     } else {
1422       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1423       if (computesubschurs) {
1424         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1425       }
1426     }
1427     if (pcbddc->adaptive_selection) {
1428       ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr);
1429       computeconstraintsmatrix = PETSC_TRUE;
1430     }
1431   }
1432 
1433   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1434   new_nearnullspace_provided = PETSC_FALSE;
1435   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1436   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1437     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1438       new_nearnullspace_provided = PETSC_TRUE;
1439     } else {
1440       /* determine if the two nullspaces are different (should be lightweight) */
1441       if (nearnullspace != pcbddc->onearnullspace) {
1442         new_nearnullspace_provided = PETSC_TRUE;
1443       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1444         PetscInt         i;
1445         const Vec        *nearnullvecs;
1446         PetscObjectState state;
1447         PetscInt         nnsp_size;
1448         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1449         for (i=0;i<nnsp_size;i++) {
1450           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1451           if (pcbddc->onearnullvecs_state[i] != state) {
1452             new_nearnullspace_provided = PETSC_TRUE;
1453             break;
1454           }
1455         }
1456       }
1457     }
1458   } else {
1459     if (!nearnullspace) { /* both nearnullspaces are null */
1460       new_nearnullspace_provided = PETSC_FALSE;
1461     } else { /* nearnullspace attached later */
1462       new_nearnullspace_provided = PETSC_TRUE;
1463     }
1464   }
1465 
1466   /* Setup constraints and related work vectors */
1467   /* reset primal space flags */
1468   pcbddc->new_primal_space = PETSC_FALSE;
1469   pcbddc->new_primal_space_local = PETSC_FALSE;
1470   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1471     /* It also sets the primal space flags */
1472     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1473     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1474     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
1475   }
1476 
1477   if (computesolvers || pcbddc->new_primal_space) {
1478     if (pcbddc->use_change_of_basis) {
1479       PC_IS *pcis = (PC_IS*)(pc->data);
1480       Mat   temp_mat = NULL;
1481 
1482       if (pcbddc->zerodiag) {
1483         /* insert B0 in pcbddc->local_mat */
1484         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_FALSE);CHKERRQ(ierr);
1485         /* hack: swap pointers */
1486         temp_mat = matis->A;
1487         matis->A = pcbddc->local_mat;
1488         pcbddc->local_mat = NULL;
1489       }
1490       ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
1491       if (pcbddc->zerodiag) {
1492         /* restore original matrix */
1493         ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
1494         matis->A = temp_mat;
1495         /* pop B0 from pcbddc->local_mat */
1496         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
1497       }
1498       /* get submatrices */
1499       ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
1500       ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
1501       ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
1502       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
1503       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
1504       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
1505       /* set flag in pcis to not reuse submatrices during PCISCreate */
1506       pcis->reusesubmatrices = PETSC_FALSE;
1507     } else if (!pcbddc->user_ChangeOfBasisMatrix) {
1508       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1509       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1510       pcbddc->local_mat = matis->A;
1511     }
1512     /* SetUp coarse and local Neumann solvers */
1513     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1514     /* SetUp Scaling operator */
1515     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
1516   }
1517 
1518   if (pcbddc->dbg_flag) {
1519     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1520   }
1521   PetscFunctionReturn(0);
1522 }
1523 
1524 /* -------------------------------------------------------------------------- */
1525 /*
1526    PCApply_BDDC - Applies the BDDC operator to a vector.
1527 
1528    Input Parameters:
1529 +  pc - the preconditioner context
1530 -  r - input vector (global)
1531 
1532    Output Parameter:
1533 .  z - output vector (global)
1534 
1535    Application Interface Routine: PCApply()
1536  */
1537 #undef __FUNCT__
1538 #define __FUNCT__ "PCApply_BDDC"
1539 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1540 {
1541   PC_IS             *pcis = (PC_IS*)(pc->data);
1542   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1543   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
1544   PetscErrorCode    ierr;
1545   const PetscScalar one = 1.0;
1546   const PetscScalar m_one = -1.0;
1547   const PetscScalar zero = 0.0;
1548 
1549 /* This code is similar to that provided in nn.c for PCNN
1550    NN interface preconditioner changed to BDDC
1551    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
1552 
1553   PetscFunctionBegin;
1554   if (!pcbddc->use_exact_dirichlet_trick) {
1555     ierr = VecCopy(r,z);CHKERRQ(ierr);
1556     /* First Dirichlet solve */
1557     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1558     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1559     /*
1560       Assembling right hand side for BDDC operator
1561       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1562       - pcis->vec1_B the interface part of the global vector z
1563     */
1564     if (n_D) {
1565       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1566       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1567       if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1568       ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1569     } else {
1570       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1571     }
1572     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1573     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1574     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1575   } else {
1576     if (pcbddc->switch_static) {
1577       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1578     }
1579     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1580   }
1581 
1582   /* Apply interface preconditioner
1583      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1584   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
1585 
1586   /* Apply transpose of partition of unity operator */
1587   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1588 
1589   /* Second Dirichlet solve and assembling of output */
1590   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1591   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1592   if (n_B) {
1593     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1594     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1595   } else if (pcbddc->switch_static) {
1596     ierr = MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1597   }
1598   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1599   if (!pcbddc->use_exact_dirichlet_trick) {
1600     if (pcbddc->switch_static) {
1601       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1602     } else {
1603       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1604     }
1605     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1606     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1607   } else {
1608     if (pcbddc->switch_static) {
1609       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1610     } else {
1611       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1612     }
1613     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1614     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1615   }
1616   PetscFunctionReturn(0);
1617 }
1618 
1619 /* -------------------------------------------------------------------------- */
1620 /*
1621    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
1622 
1623    Input Parameters:
1624 +  pc - the preconditioner context
1625 -  r - input vector (global)
1626 
1627    Output Parameter:
1628 .  z - output vector (global)
1629 
1630    Application Interface Routine: PCApplyTranspose()
1631  */
1632 #undef __FUNCT__
1633 #define __FUNCT__ "PCApplyTranspose_BDDC"
1634 PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
1635 {
1636   PC_IS             *pcis = (PC_IS*)(pc->data);
1637   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1638   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
1639   PetscErrorCode    ierr;
1640   const PetscScalar one = 1.0;
1641   const PetscScalar m_one = -1.0;
1642   const PetscScalar zero = 0.0;
1643 
1644   PetscFunctionBegin;
1645   if (!pcbddc->use_exact_dirichlet_trick) {
1646     ierr = VecCopy(r,z);CHKERRQ(ierr);
1647     /* First Dirichlet solve */
1648     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1649     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1650     /*
1651       Assembling right hand side for BDDC operator
1652       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1653       - pcis->vec1_B the interface part of the global vector z
1654     */
1655     if (n_D) {
1656       ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1657       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
1658       if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1659       ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1660     } else {
1661       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1662     }
1663     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1664     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1665     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1666   } else {
1667     if (pcbddc->switch_static) {
1668       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1669     }
1670     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1671   }
1672 
1673   /* Apply interface preconditioner
1674      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1675   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
1676 
1677   /* Apply transpose of partition of unity operator */
1678   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
1679 
1680   /* Second Dirichlet solve and assembling of output */
1681   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1682   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1683   if (n_B) {
1684     ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
1685     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1686   } else if (pcbddc->switch_static) {
1687     ierr = MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1688   }
1689   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1690   if (!pcbddc->use_exact_dirichlet_trick) {
1691     if (pcbddc->switch_static) {
1692       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1693     } else {
1694       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1695     }
1696     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1697     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1698   } else {
1699     if (pcbddc->switch_static) {
1700       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1701     } else {
1702       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1703     }
1704     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1705     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1706   }
1707   PetscFunctionReturn(0);
1708 }
1709 /* -------------------------------------------------------------------------- */
1710 
1711 #undef __FUNCT__
1712 #define __FUNCT__ "PCDestroy_BDDC"
1713 PetscErrorCode PCDestroy_BDDC(PC pc)
1714 {
1715   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1716   PetscErrorCode ierr;
1717 
1718   PetscFunctionBegin;
1719   /* free data created by PCIS */
1720   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1721   /* free BDDC custom data  */
1722   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1723   /* destroy objects related to topography */
1724   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1725   /* free allocated graph structure */
1726   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1727   /* free allocated sub schurs structure */
1728   ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr);
1729   /* destroy objects for scaling operator */
1730   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1731   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
1732   /* free solvers stuff */
1733   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1734   /* free global vectors needed in presolve */
1735   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
1736   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1737   /* free stuff for change of basis hooks */
1738   if (pcbddc->new_global_mat) {
1739     PCBDDCChange_ctx change_ctx;
1740     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1741     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1742     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1743     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1744     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1745   }
1746   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
1747   /* remove functions */
1748   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1749   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1750   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
1751   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1752   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
1753   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1754   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1755   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1756   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1757   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1758   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1759   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1760   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1761   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1762   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1763   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1764   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1765   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1766   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1767   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1768   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1769   /* Free the private data structure */
1770   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1771   PetscFunctionReturn(0);
1772 }
1773 /* -------------------------------------------------------------------------- */
1774 
1775 #undef __FUNCT__
1776 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
1777 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1778 {
1779   FETIDPMat_ctx  mat_ctx;
1780   Vec            copy_standard_rhs;
1781   PC_IS*         pcis;
1782   PC_BDDC*       pcbddc;
1783   PetscErrorCode ierr;
1784 
1785   PetscFunctionBegin;
1786   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1787   pcis = (PC_IS*)mat_ctx->pc->data;
1788   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1789 
1790   /*
1791      change of basis for physical rhs if needed
1792      It also changes the rhs in case of dirichlet boundaries
1793      TODO: better management when FETIDP will have its own class
1794   */
1795   ierr = VecDuplicate(standard_rhs,&copy_standard_rhs);CHKERRQ(ierr);
1796   ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr);
1797   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr);
1798   /* store vectors for computation of fetidp final solution */
1799   ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1800   ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1801   /* scale rhs since it should be unassembled */
1802   /* TODO use counter scaling? (also below) */
1803   ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1804   ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1805   /* Apply partition of unity */
1806   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1807   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1808   if (!pcbddc->switch_static) {
1809     /* compute partially subassembled Schur complement right-hand side */
1810     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1811     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
1812     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1813     ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr);
1814     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1815     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1816     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1817     ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1818     ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1819     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1820   }
1821   ierr = VecDestroy(&copy_standard_rhs);CHKERRQ(ierr);
1822   /* BDDC rhs */
1823   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
1824   if (pcbddc->switch_static) {
1825     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1826   }
1827   /* apply BDDC */
1828   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
1829   /* Application of B_delta and assembling of rhs for fetidp fluxes */
1830   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
1831   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
1832   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1833   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1834   PetscFunctionReturn(0);
1835 }
1836 
1837 #undef __FUNCT__
1838 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
1839 /*@
1840  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side
1841 
1842    Collective
1843 
1844    Input Parameters:
1845 +  fetidp_mat      - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
1846 -  standard_rhs    - the right-hand side of the original linear system
1847 
1848    Output Parameters:
1849 .  fetidp_flux_rhs - the right-hand side for the FETI-DP linear system
1850 
1851    Level: developer
1852 
1853    Notes:
1854 
1855 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
1856 @*/
1857 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1858 {
1859   FETIDPMat_ctx  mat_ctx;
1860   PetscErrorCode ierr;
1861 
1862   PetscFunctionBegin;
1863   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1864   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
1865   PetscFunctionReturn(0);
1866 }
1867 /* -------------------------------------------------------------------------- */
1868 
1869 #undef __FUNCT__
1870 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
1871 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1872 {
1873   FETIDPMat_ctx  mat_ctx;
1874   PC_IS*         pcis;
1875   PC_BDDC*       pcbddc;
1876   PetscErrorCode ierr;
1877 
1878   PetscFunctionBegin;
1879   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1880   pcis = (PC_IS*)mat_ctx->pc->data;
1881   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1882 
1883   /* apply B_delta^T */
1884   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1885   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1886   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
1887   /* compute rhs for BDDC application */
1888   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1889   if (pcbddc->switch_static) {
1890     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1891   }
1892   /* apply BDDC */
1893   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
1894   /* put values into standard global vector */
1895   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1896   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1897   if (!pcbddc->switch_static) {
1898     /* compute values into the interior if solved for the partially subassembled Schur complement */
1899     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
1900     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
1901     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
1902   }
1903   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1904   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1905   /* final change of basis if needed
1906      Is also sums the dirichlet part removed during RHS assembling */
1907   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
1908   PetscFunctionReturn(0);
1909 }
1910 
1911 #undef __FUNCT__
1912 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
1913 /*@
1914  PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system
1915 
1916    Collective
1917 
1918    Input Parameters:
1919 +  fetidp_mat      - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
1920 -  fetidp_flux_sol - the solution of the FETI-DP linear system
1921 
1922    Output Parameters:
1923 .  standard_sol    - the solution defined on the physical domain
1924 
1925    Level: developer
1926 
1927    Notes:
1928 
1929 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
1930 @*/
1931 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1932 {
1933   FETIDPMat_ctx  mat_ctx;
1934   PetscErrorCode ierr;
1935 
1936   PetscFunctionBegin;
1937   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1938   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
1939   PetscFunctionReturn(0);
1940 }
1941 /* -------------------------------------------------------------------------- */
1942 
1943 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1944 extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1945 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1946 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1947 extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1948 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1949 
1950 #undef __FUNCT__
1951 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
1952 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1953 {
1954 
1955   FETIDPMat_ctx  fetidpmat_ctx;
1956   Mat            newmat;
1957   FETIDPPC_ctx   fetidppc_ctx;
1958   PC             newpc;
1959   MPI_Comm       comm;
1960   PetscErrorCode ierr;
1961 
1962   PetscFunctionBegin;
1963   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1964   /* FETIDP linear matrix */
1965   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
1966   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
1967   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
1968   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1969   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
1970   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
1971   ierr = MatSetUp(newmat);CHKERRQ(ierr);
1972   /* FETIDP preconditioner */
1973   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
1974   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
1975   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
1976   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
1977   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
1978   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1979   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
1980   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
1981   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
1982   ierr = PCSetUp(newpc);CHKERRQ(ierr);
1983   /* return pointers for objects created */
1984   *fetidp_mat=newmat;
1985   *fetidp_pc=newpc;
1986   PetscFunctionReturn(0);
1987 }
1988 
1989 #undef __FUNCT__
1990 #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
1991 /*@
1992  PCBDDCCreateFETIDPOperators - Create FETI-DP operators
1993 
1994    Collective
1995 
1996    Input Parameters:
1997 .  pc - the BDDC preconditioning context (setup should have been called before)
1998 
1999    Output Parameters:
2000 +  fetidp_mat - shell FETI-DP matrix object
2001 -  fetidp_pc  - shell Dirichlet preconditioner for FETI-DP matrix
2002 
2003    Options Database Keys:
2004 .    -fetidp_fullyredundant <false> - use or not a fully redundant set of Lagrange multipliers
2005 
2006    Level: developer
2007 
2008    Notes:
2009      Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose
2010 
2011 .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
2012 @*/
2013 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
2014 {
2015   PetscErrorCode ierr;
2016 
2017   PetscFunctionBegin;
2018   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2019   if (pc->setupcalled) {
2020     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
2021   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
2022   PetscFunctionReturn(0);
2023 }
2024 /* -------------------------------------------------------------------------- */
2025 /*MC
2026    PCBDDC - Balancing Domain Decomposition by Constraints.
2027 
2028    An implementation of the BDDC preconditioner based on
2029 
2030 .vb
2031    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
2032    [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
2033    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
2034    [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
2035 .ve
2036 
2037    The matrix to be preconditioned (Pmat) must be of type MATIS.
2038 
2039    Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
2040 
2041    It also works with unsymmetric and indefinite problems.
2042 
2043    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.
2044 
2045    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace()
2046 
2047    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()
2048    Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesLocalIS()
2049 
2050    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.
2051 
2052    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.
2053    User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()
2054 
2055    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.
2056 
2057    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.
2058 
2059    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.
2060    Deluxe scaling is not supported yet for FETI-DP.
2061 
2062    Options Database Keys (some of them, run with -h for a complete list):
2063 
2064 .    -pc_bddc_use_vertices <true> - use or not vertices in primal space
2065 .    -pc_bddc_use_edges <true> - use or not edges in primal space
2066 .    -pc_bddc_use_faces <false> - use or not faces in primal space
2067 .    -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
2068 .    -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
2069 .    -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
2070 .    -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
2071 .    -pc_bddc_levels <0> - maximum number of levels for multilevel
2072 .    -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)
2073 .    -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)
2074 .    -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
2075 .    -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)
2076 .    -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)
2077 -    -pc_bddc_check_level <0> - set verbosity level of debugging output
2078 
2079    Options for Dirichlet, Neumann or coarse solver can be set with
2080 .vb
2081       -pc_bddc_dirichlet_
2082       -pc_bddc_neumann_
2083       -pc_bddc_coarse_
2084 .ve
2085    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.
2086 
2087    When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
2088 .vb
2089       -pc_bddc_dirichlet_lN_
2090       -pc_bddc_neumann_lN_
2091       -pc_bddc_coarse_lN_
2092 .ve
2093    Note that level number ranges from the finest (0) to the coarsest (N).
2094    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.
2095 .vb
2096      -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
2097 .ve
2098    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
2099 
2100    Level: intermediate
2101 
2102    Developer notes:
2103 
2104    Contributed by Stefano Zampini
2105 
2106 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
2107 M*/
2108 
2109 #undef __FUNCT__
2110 #define __FUNCT__ "PCCreate_BDDC"
2111 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
2112 {
2113   PetscErrorCode      ierr;
2114   PC_BDDC             *pcbddc;
2115 
2116   PetscFunctionBegin;
2117   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
2118   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
2119   pc->data  = (void*)pcbddc;
2120 
2121   /* create PCIS data structure */
2122   ierr = PCISCreate(pc);CHKERRQ(ierr);
2123 
2124   /* BDDC customization */
2125   pcbddc->use_local_adj       = PETSC_TRUE;
2126   pcbddc->use_vertices        = PETSC_TRUE;
2127   pcbddc->use_edges           = PETSC_TRUE;
2128   pcbddc->use_faces           = PETSC_FALSE;
2129   pcbddc->use_change_of_basis = PETSC_FALSE;
2130   pcbddc->use_change_on_faces = PETSC_FALSE;
2131   pcbddc->switch_static       = PETSC_FALSE;
2132   pcbddc->use_nnsp_true       = PETSC_FALSE;
2133   pcbddc->use_qr_single       = PETSC_FALSE;
2134   pcbddc->symmetric_primal    = PETSC_TRUE;
2135   pcbddc->benign_saddle_point = PETSC_FALSE;
2136   pcbddc->dbg_flag            = 0;
2137   /* private */
2138   pcbddc->local_primal_size          = 0;
2139   pcbddc->local_primal_size_cc       = 0;
2140   pcbddc->local_primal_ref_node      = 0;
2141   pcbddc->local_primal_ref_mult      = 0;
2142   pcbddc->n_vertices                 = 0;
2143   pcbddc->primal_indices_local_idxs  = 0;
2144   pcbddc->recompute_topography       = PETSC_FALSE;
2145   pcbddc->coarse_size                = -1;
2146   pcbddc->new_primal_space           = PETSC_FALSE;
2147   pcbddc->new_primal_space_local     = PETSC_FALSE;
2148   pcbddc->global_primal_indices      = 0;
2149   pcbddc->onearnullspace             = 0;
2150   pcbddc->onearnullvecs_state        = 0;
2151   pcbddc->user_primal_vertices       = 0;
2152   pcbddc->NullSpace                  = 0;
2153   pcbddc->temp_solution              = 0;
2154   pcbddc->original_rhs               = 0;
2155   pcbddc->local_mat                  = 0;
2156   pcbddc->ChangeOfBasisMatrix        = 0;
2157   pcbddc->user_ChangeOfBasisMatrix   = 0;
2158   pcbddc->new_global_mat             = 0;
2159   pcbddc->coarse_vec                 = 0;
2160   pcbddc->coarse_ksp                 = 0;
2161   pcbddc->coarse_phi_B               = 0;
2162   pcbddc->coarse_phi_D               = 0;
2163   pcbddc->coarse_psi_B               = 0;
2164   pcbddc->coarse_psi_D               = 0;
2165   pcbddc->vec1_P                     = 0;
2166   pcbddc->vec1_R                     = 0;
2167   pcbddc->vec2_R                     = 0;
2168   pcbddc->local_auxmat1              = 0;
2169   pcbddc->local_auxmat2              = 0;
2170   pcbddc->R_to_B                     = 0;
2171   pcbddc->R_to_D                     = 0;
2172   pcbddc->ksp_D                      = 0;
2173   pcbddc->ksp_R                      = 0;
2174   pcbddc->NeumannBoundaries          = 0;
2175   pcbddc->NeumannBoundariesLocal     = 0;
2176   pcbddc->DirichletBoundaries        = 0;
2177   pcbddc->DirichletBoundariesLocal   = 0;
2178   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
2179   pcbddc->n_ISForDofs                = 0;
2180   pcbddc->n_ISForDofsLocal           = 0;
2181   pcbddc->ISForDofs                  = 0;
2182   pcbddc->ISForDofsLocal             = 0;
2183   pcbddc->ConstraintMatrix           = 0;
2184   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
2185   pcbddc->coarse_loc_to_glob         = 0;
2186   pcbddc->coarsening_ratio           = 8;
2187   pcbddc->coarse_adj_red             = 0;
2188   pcbddc->current_level              = 0;
2189   pcbddc->max_levels                 = 0;
2190   pcbddc->use_coarse_estimates       = PETSC_FALSE;
2191   pcbddc->redistribute_coarse        = 0;
2192   pcbddc->coarse_subassembling       = 0;
2193   pcbddc->coarse_subassembling_init  = 0;
2194 
2195   /* benign subspace trick */
2196   pcbddc->zerodiag                   = 0;
2197   pcbddc->B0_ncol                    = 0;
2198   pcbddc->B0_cols                    = NULL;
2199   pcbddc->B0_vals                    = NULL;
2200   pcbddc->benign_change              = 0;
2201 
2202   /* create local graph structure */
2203   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
2204 
2205   /* scaling */
2206   pcbddc->work_scaling          = 0;
2207   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
2208   pcbddc->faster_deluxe         = PETSC_FALSE;
2209 
2210   /* create sub schurs structure */
2211   ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr);
2212   pcbddc->sub_schurs_rebuild     = PETSC_FALSE;
2213   pcbddc->sub_schurs_layers      = -1;
2214   pcbddc->sub_schurs_use_useradj = PETSC_FALSE;
2215 
2216   pcbddc->computed_rowadj = PETSC_FALSE;
2217 
2218   /* adaptivity */
2219   pcbddc->adaptive_threshold      = 0.0;
2220   pcbddc->adaptive_nmax           = 0;
2221   pcbddc->adaptive_nmin           = 0;
2222 
2223   /* function pointers */
2224   pc->ops->apply               = PCApply_BDDC;
2225   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2226   pc->ops->setup               = PCSetUp_BDDC;
2227   pc->ops->destroy             = PCDestroy_BDDC;
2228   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
2229   pc->ops->view                = 0;
2230   pc->ops->applyrichardson     = 0;
2231   pc->ops->applysymmetricleft  = 0;
2232   pc->ops->applysymmetricright = 0;
2233   pc->ops->presolve            = PCPreSolve_BDDC;
2234   pc->ops->postsolve           = PCPostSolve_BDDC;
2235 
2236   /* composing function */
2237   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
2238   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
2239   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
2240   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
2241   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
2242   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
2243   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
2244   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
2245   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2246   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
2247   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2248   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
2249   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2250   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
2251   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2252   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
2253   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
2254   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
2255   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
2256   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
2257   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
2258   PetscFunctionReturn(0);
2259 }
2260 
2261