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