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