xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 71783a16cd6b13dae51676bf855715a2a2b3171b)
153cdbc3dSStefano Zampini /* TODOLIST
2eb97c9d2SStefano Zampini 
3eb97c9d2SStefano Zampini    Solvers
4a0d3c3abSStefano Zampini    - Add support for cholesky for coarse solver (similar to local solvers)
5eb97c9d2SStefano Zampini    - Propagate ksp prefixes for solvers to mat objects?
6eb97c9d2SStefano Zampini 
7eb97c9d2SStefano Zampini    User interface
80f202f7eSStefano Zampini    - ** DM attached to pc?
9eb97c9d2SStefano Zampini 
10eb97c9d2SStefano Zampini    Debugging output
11b9b85e73SStefano Zampini    - * Better management of verbosity levels of debugging output
12eb97c9d2SStefano Zampini 
13eb97c9d2SStefano Zampini    Extra
14b9b85e73SStefano Zampini    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
15eb97c9d2SStefano Zampini    - BDDC with MG framework?
16eb97c9d2SStefano Zampini 
17eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
18eb97c9d2SStefano Zampini    - Provide general case for subassembling
19eb97c9d2SStefano Zampini 
2053cdbc3dSStefano Zampini */
210c7d97c5SJed Brown 
22ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
23ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
243b03a366Sstefano_zampini #include <petscblaslapack.h>
25674ae819SStefano Zampini 
26f3d41395Sstefano_zampini static PetscBool  cited = PETSC_FALSE;
27f3d41395Sstefano_zampini static const char citation[] =
28f3d41395Sstefano_zampini "@article{ZampiniPCBDDC,\n"
29f3d41395Sstefano_zampini "author = {Stefano Zampini},\n"
30f3d41395Sstefano_zampini "title = {{PCBDDC}: A Class of Robust Dual-Primal Methods in {PETS}c},\n"
31f3d41395Sstefano_zampini "journal = {SIAM Journal on Scientific Computing},\n"
32f3d41395Sstefano_zampini "volume = {38},\n"
33f3d41395Sstefano_zampini "number = {5},\n"
34f3d41395Sstefano_zampini "pages = {S282-S306},\n"
35f3d41395Sstefano_zampini "year = {2016},\n"
36f3d41395Sstefano_zampini "doi = {10.1137/15M1025785},\n"
37f3d41395Sstefano_zampini "URL = {http://dx.doi.org/10.1137/15M1025785},\n"
38f3d41395Sstefano_zampini "eprint = {http://dx.doi.org/10.1137/15M1025785}\n"
39f3d41395Sstefano_zampini "}\n";
40f3d41395Sstefano_zampini 
410369aaf7SStefano Zampini /* temporarily declare it */
420369aaf7SStefano Zampini PetscErrorCode PCApply_BDDC(PC,Vec,Vec);
430369aaf7SStefano Zampini 
444416b707SBarry Smith PetscErrorCode PCSetFromOptions_BDDC(PetscOptionItems *PetscOptionsObject,PC pc)
450c7d97c5SJed Brown {
460c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
47bd2a564bSStefano Zampini   PetscInt       nt;
480c7d97c5SJed Brown   PetscErrorCode ierr;
490c7d97c5SJed Brown 
500c7d97c5SJed Brown   PetscFunctionBegin;
51e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"BDDC options");CHKERRQ(ierr);
528eeda7d8SStefano Zampini   /* Verbose debugging */
53a13144ffSStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
54a13144ffSStefano Zampini   /* Approximate solvers */
55c7017625SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_dirichlet_approximate","Inform PCBDDC that we are using approximate Dirichlet solvers","none",pcbddc->NullSpace_corr[0],&pcbddc->NullSpace_corr[0],NULL);CHKERRQ(ierr);
56c7017625SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_dirichlet_approximate_scale","Inform PCBDDC that we need to scale the Dirichlet solve","none",pcbddc->NullSpace_corr[1],&pcbddc->NullSpace_corr[1],NULL);CHKERRQ(ierr);
57c7017625SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_neumann_approximate","Inform PCBDDC that we are using approximate Neumann solvers","none",pcbddc->NullSpace_corr[2],&pcbddc->NullSpace_corr[2],NULL);CHKERRQ(ierr);
58c7017625SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_neumann_approximate_scale","Inform PCBDDC that we need to scale the Neumann solve","none",pcbddc->NullSpace_corr[3],&pcbddc->NullSpace_corr[3],NULL);CHKERRQ(ierr);
596b78500eSPatrick Sanan   /* Primal space customization */
6008a5cf49SStefano Zampini   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);
61be12c134Sstefano_zampini   ierr = PetscOptionsInt("-pc_bddc_graph_maxcount","Maximum number of shared subdomains for a connected component","none",pcbddc->graphmaxcount,&pcbddc->graphmaxcount,NULL);CHKERRQ(ierr);
628eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
638eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
648eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
6514f95afaSStefano Zampini   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);
666661aa1dSStefano Zampini   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);
6714f95afaSStefano Zampini   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);
688eeda7d8SStefano Zampini   /* Change of basis */
69b9b85e73SStefano Zampini   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);
70b9b85e73SStefano Zampini   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);
71674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
72674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
73674ae819SStefano Zampini   }
748eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
758eeda7d8SStefano Zampini   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);
7657de7509SStefano Zampini   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);
770298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
782b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
79323d395dSStefano Zampini   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);
80674ae819SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
81b96c3477SStefano Zampini   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);
82b96c3477SStefano Zampini   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);
83b96c3477SStefano Zampini   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);
84683d3df6SStefano Zampini   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);
85bf3a8328SStefano Zampini   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);
86839e9adbSstefano_zampini   ierr = PetscOptionsBool("-pc_bddc_deluxe_singlemat","Collapse deluxe operators","none",pcbddc->deluxe_singlemat,&pcbddc->deluxe_singlemat,NULL);CHKERRQ(ierr);
87bf3a8328SStefano Zampini   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);
88bd2a564bSStefano Zampini   nt   = 2;
89bd2a564bSStefano Zampini   ierr = PetscOptionsRealArray("-pc_bddc_adaptive_threshold","Thresholds to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&nt,NULL);CHKERRQ(ierr);
90bd2a564bSStefano Zampini   if (nt == 1) pcbddc->adaptive_threshold[1] = pcbddc->adaptive_threshold[0];
9108122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);CHKERRQ(ierr);
9208122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);CHKERRQ(ierr);
933301b35fSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);CHKERRQ(ierr);
94b0c7d250SStefano Zampini   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);
95b3afcdbeSStefano Zampini   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);
96e9627c49SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_benign_change","Compute the pressure change of basis explicitly","none",pcbddc->benign_change_explicit,&pcbddc->benign_change_explicit,NULL);CHKERRQ(ierr);
9727b6a85dSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_benign_compute_correction","Compute the benign correction during PreSolve","none",pcbddc->benign_compute_correction,&pcbddc->benign_compute_correction,NULL);CHKERRQ(ierr);
98a198735bSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_nonetflux","Automatic computation of no-net-flux quadrature weights","none",pcbddc->compute_nonetflux,&pcbddc->compute_nonetflux,NULL);CHKERRQ(ierr);
994f1b2e48SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_detect_disconnected","Detects disconnected subdomains","none",pcbddc->detect_disconnected,&pcbddc->detect_disconnected,NULL);CHKERRQ(ierr);
10070c64980SStefano Zampini   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);
1010c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
1020c7d97c5SJed Brown   PetscFunctionReturn(0);
1030c7d97c5SJed Brown }
1046b78500eSPatrick Sanan 
1056b78500eSPatrick Sanan static PetscErrorCode PCView_BDDC(PC pc,PetscViewer viewer)
1066b78500eSPatrick Sanan {
1076b78500eSPatrick Sanan   PC_BDDC              *pcbddc = (PC_BDDC*)pc->data;
108e9627c49SStefano Zampini   PC_IS                *pcis = (PC_IS*)pc->data;
1096b78500eSPatrick Sanan   PetscErrorCode       ierr;
110*71783a16SStefano Zampini   PetscBool            isascii;
111e9627c49SStefano Zampini   PetscSubcomm         subcomm;
112e9627c49SStefano Zampini   PetscViewer          subviewer;
1136b78500eSPatrick Sanan 
1146b78500eSPatrick Sanan   PetscFunctionBegin;
1156b78500eSPatrick Sanan   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1166b78500eSPatrick Sanan   /* ASCII viewer */
1176b78500eSPatrick Sanan   if (isascii) {
1184b2aedd3SStefano Zampini     PetscMPIInt   color,rank,size;
119fbad9177SStefano Zampini     PetscInt64    loc[7],gsum[6],gmax[6],gmin[6],totbenign;
120e9627c49SStefano Zampini     PetscScalar   interface_size;
121e9627c49SStefano Zampini     PetscReal     ratio1=0.,ratio2=0.;
122e9627c49SStefano Zampini     Vec           counter;
1236b78500eSPatrick Sanan 
124b74ba07aSstefano_zampini     if (!pc->setupcalled) {
125b74ba07aSstefano_zampini       ierr = PetscViewerASCIIPrintf(viewer,"  Partial information available: preconditioner has not been setup yet\n");CHKERRQ(ierr);
126b74ba07aSstefano_zampini     }
127efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use verbose output: %D\n",pcbddc->dbg_flag);CHKERRQ(ierr);
1286f0c0a6aSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Use user-defined CSR: %d\n",!!pcbddc->mat_graph->nvtxs_csr);CHKERRQ(ierr);
1296f0c0a6aSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Use local mat graph: %d\n",pcbddc->use_local_adj && !pcbddc->mat_graph->nvtxs_csr);CHKERRQ(ierr);
130e9627c49SStefano Zampini     if (pcbddc->mat_graph->twodim) {
131efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Connectivity graph topological dimension: 2\n");CHKERRQ(ierr);
132e9627c49SStefano Zampini     } else {
133efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Connectivity graph topological dimension: 3\n");CHKERRQ(ierr);
134e9627c49SStefano Zampini     }
135aefa1729SStefano Zampini     if (pcbddc->graphmaxcount != PETSC_MAX_INT) {
136efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Graph max count: %D\n",pcbddc->graphmaxcount);CHKERRQ(ierr);
137aefa1729SStefano Zampini     }
13850e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Use vertices: %d (vertex size %D)\n",pcbddc->use_vertices,pcbddc->vertex_size);CHKERRQ(ierr);
139efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use edges: %d\n",pcbddc->use_edges);CHKERRQ(ierr);
140efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use faces: %d\n",pcbddc->use_faces);CHKERRQ(ierr);
141efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use true near null space: %d\n",pcbddc->use_nnsp_true);CHKERRQ(ierr);
142efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use QR for single constraints on cc: %d\n",pcbddc->use_qr_single);CHKERRQ(ierr);
143efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use change of basis on local edge nodes: %d\n",pcbddc->use_change_of_basis);CHKERRQ(ierr);
144efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use change of basis on local face nodes: %d\n",pcbddc->use_change_on_faces);CHKERRQ(ierr);
145efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  User defined change of basis matrix: %d\n",!!pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
146efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Has change of basis matrix: %d\n",!!pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
147efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Eliminate dirichlet boundary dofs: %d\n",pcbddc->eliminate_dirdofs);CHKERRQ(ierr);
148efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Switch on static condensation ops around the interface preconditioner: %d\n",pcbddc->switch_static);CHKERRQ(ierr);
149efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use exact dirichlet trick: %d\n",pcbddc->use_exact_dirichlet_trick);CHKERRQ(ierr);
15050e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Multilevel max levels: %D\n",pcbddc->max_levels);CHKERRQ(ierr);
15150e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Multilevel coarsening ratio: %D\n",pcbddc->coarsening_ratio);CHKERRQ(ierr);
152efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use estimated eigs for coarse problem: %d\n",pcbddc->use_coarse_estimates);CHKERRQ(ierr);
153efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use deluxe scaling: %d\n",pcbddc->use_deluxe_scaling);CHKERRQ(ierr);
154efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use deluxe zerorows: %d\n",pcbddc->deluxe_zerorows);CHKERRQ(ierr);
155efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use deluxe singlemat: %d\n",pcbddc->deluxe_singlemat);CHKERRQ(ierr);
156efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Rebuild interface graph for Schur principal minors: %d\n",pcbddc->sub_schurs_rebuild);CHKERRQ(ierr);
15750e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Number of dofs' layers for the computation of principal minors: %D\n",pcbddc->sub_schurs_layers);CHKERRQ(ierr);
158efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Use user CSR graph to compute successive layers: %d\n",pcbddc->sub_schurs_use_useradj);CHKERRQ(ierr);
159bd2a564bSStefano Zampini     if (pcbddc->adaptive_threshold[1] != pcbddc->adaptive_threshold[0]) {
160bd2a564bSStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"  Adaptive constraint selection thresholds (active %d, userdefined %d): %g,%g\n",pcbddc->adaptive_selection,pcbddc->adaptive_userdefined,pcbddc->adaptive_threshold[0],pcbddc->adaptive_threshold[1]);CHKERRQ(ierr);
161bd2a564bSStefano Zampini     } else {
162bd2a564bSStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"  Adaptive constraint selection threshold (active %d, userdefined %d): %g\n",pcbddc->adaptive_selection,pcbddc->adaptive_userdefined,pcbddc->adaptive_threshold[0]);CHKERRQ(ierr);
163bd2a564bSStefano Zampini     }
16450e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Min constraints / connected component: %D\n",pcbddc->adaptive_nmin);CHKERRQ(ierr);
16550e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Max constraints / connected component: %D\n",pcbddc->adaptive_nmax);CHKERRQ(ierr);
166efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Invert exact Schur complement for adaptive selection: %d\n",pcbddc->sub_schurs_exact_schur);CHKERRQ(ierr);
167efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Symmetric computation of primal basis functions: %d\n",pcbddc->symmetric_primal);CHKERRQ(ierr);
16850e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Num. Procs. to map coarse adjacency list: %D\n",pcbddc->coarse_adj_red);CHKERRQ(ierr);
16950e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Coarse eqs per proc (significant at the coarsest level): %D\n",pcbddc->coarse_eqs_per_proc);CHKERRQ(ierr);
170efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Detect disconnected: %d\n",pcbddc->detect_disconnected);CHKERRQ(ierr);
171efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Benign subspace trick: %d (change explicit %d)\n",pcbddc->benign_saddle_point,pcbddc->benign_change_explicit);CHKERRQ(ierr);
172efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Benign subspace trick is active: %d\n",pcbddc->benign_have_null);CHKERRQ(ierr);
173efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Algebraic computation of no-net-flux %d\n",pcbddc->compute_nonetflux);CHKERRQ(ierr);
174b74ba07aSstefano_zampini     if (!pc->setupcalled) PetscFunctionReturn(0);
1756b78500eSPatrick Sanan 
176fbad9177SStefano Zampini     /* compute interface size */
177e9627c49SStefano Zampini     ierr = VecSet(pcis->vec1_B,1.0);CHKERRQ(ierr);
178e9627c49SStefano Zampini     ierr = MatCreateVecs(pc->pmat,&counter,0);CHKERRQ(ierr);
179e9627c49SStefano Zampini     ierr = VecSet(counter,0.0);CHKERRQ(ierr);
180e9627c49SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
181e9627c49SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,counter,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
182e9627c49SStefano Zampini     ierr = VecSum(counter,&interface_size);CHKERRQ(ierr);
183e9627c49SStefano Zampini     ierr = VecDestroy(&counter);CHKERRQ(ierr);
184fbad9177SStefano Zampini 
185fbad9177SStefano Zampini     /* compute some statistics on the domain decomposition */
186e9627c49SStefano Zampini     gsum[0] = 1;
187fbad9177SStefano Zampini     gsum[1] = gsum[2] = gsum[3] = gsum[4] = gsum[5] = 0;
188e9627c49SStefano Zampini     loc[0]  = !!pcis->n;
189e9627c49SStefano Zampini     loc[1]  = pcis->n - pcis->n_B;
190e9627c49SStefano Zampini     loc[2]  = pcis->n_B;
191e9627c49SStefano Zampini     loc[3]  = pcbddc->local_primal_size;
192345ecf6cSStefano Zampini     loc[4]  = pcis->n;
193fbad9177SStefano Zampini     loc[5]  = pcbddc->n_local_subs > 0 ? pcbddc->n_local_subs : (pcis->n ? 1 : 0);
194fbad9177SStefano Zampini     loc[6]  = pcbddc->benign_n;
195fbad9177SStefano Zampini     ierr = MPI_Reduce(loc,gsum,6,MPIU_INT64,MPI_SUM,0,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
196fbad9177SStefano Zampini     if (!loc[0]) loc[1] = loc[2] = loc[3] = loc[4] = loc[5] = -1;
197fbad9177SStefano Zampini     ierr = MPI_Reduce(loc,gmax,6,MPIU_INT64,MPI_MAX,0,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
198fbad9177SStefano Zampini     if (!loc[0]) loc[1] = loc[2] = loc[3] = loc[4] = loc[5] = PETSC_MAX_INT;
199fbad9177SStefano Zampini     ierr = MPI_Reduce(loc,gmin,6,MPIU_INT64,MPI_MIN,0,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
200fbad9177SStefano Zampini     ierr = MPI_Reduce(&loc[6],&totbenign,1,MPIU_INT64,MPI_SUM,0,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
201e9627c49SStefano Zampini     if (pcbddc->coarse_size) {
202e9627c49SStefano Zampini       ratio1 = pc->pmat->rmap->N/(1.*pcbddc->coarse_size);
203e9627c49SStefano Zampini       ratio2 = PetscRealPart(interface_size)/pcbddc->coarse_size;
204e9627c49SStefano Zampini     }
205efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  ********************************** STATISTICS AT LEVEL %d **********************************\n",pcbddc->current_level);CHKERRQ(ierr);
20650e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Global dofs sizes: all %D interface %D coarse %D\n",pc->pmat->rmap->N,(PetscInt)PetscRealPart(interface_size),pcbddc->coarse_size);CHKERRQ(ierr);
20750e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Coarsening ratios: all/coarse %D interface/coarse %D\n",(PetscInt)ratio1,(PetscInt)ratio2);CHKERRQ(ierr);
20850e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Active processes : %D\n",(PetscInt)gsum[0]);CHKERRQ(ierr);
20950e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Total subdomains : %D\n",(PetscInt)gsum[5]);CHKERRQ(ierr);
210345ecf6cSStefano Zampini     if (pcbddc->benign_have_null) {
21150e0721cSStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"  Benign subs      : %D\n",(PetscInt)totbenign);CHKERRQ(ierr);
212345ecf6cSStefano Zampini     }
21350e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Dofs type        :\tMIN\tMAX\tMEAN\n");CHKERRQ(ierr);
21450e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Interior  dofs   :\t%D\t%D\t%D\n",(PetscInt)gmin[1],(PetscInt)gmax[1],(PetscInt)(gsum[1]/gsum[0]));CHKERRQ(ierr);
21550e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Interface dofs   :\t%D\t%D\t%D\n",(PetscInt)gmin[2],(PetscInt)gmax[2],(PetscInt)(gsum[2]/gsum[0]));CHKERRQ(ierr);
21650e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Primal    dofs   :\t%D\t%D\t%D\n",(PetscInt)gmin[3],(PetscInt)gmax[3],(PetscInt)(gsum[3]/gsum[0]));CHKERRQ(ierr);
21750e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Local     dofs   :\t%D\t%D\t%D\n",(PetscInt)gmin[4],(PetscInt)gmax[4],(PetscInt)(gsum[4]/gsum[0]));CHKERRQ(ierr);
21850e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  Local     subs   :\t%D\t%D\n"    ,(PetscInt)gmin[5],(PetscInt)gmax[5]);CHKERRQ(ierr);
21950e0721cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  ********************************** COARSE PROBLEM DETAILS *********************************\n");CHKERRQ(ierr);
22027b6a85dSStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
221e9627c49SStefano Zampini 
222fbad9177SStefano Zampini     /* the coarse problem can be handled by a different communicator */
223e9627c49SStefano Zampini     if (pcbddc->coarse_ksp) color = 1;
224e9627c49SStefano Zampini     else color = 0;
225e9627c49SStefano Zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
2264b2aedd3SStefano Zampini     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
227e9627c49SStefano Zampini     ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)pc),&subcomm);CHKERRQ(ierr);
2284b2aedd3SStefano Zampini     ierr = PetscSubcommSetNumber(subcomm,PetscMin(size,2));CHKERRQ(ierr);
229e9627c49SStefano Zampini     ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
230e9627c49SStefano Zampini     ierr = PetscViewerGetSubViewer(viewer,PetscSubcommChild(subcomm),&subviewer);CHKERRQ(ierr);
231e9627c49SStefano Zampini     if (color == 1) {
232e9627c49SStefano Zampini       ierr = KSPView(pcbddc->coarse_ksp,subviewer);CHKERRQ(ierr);
233e9627c49SStefano Zampini       ierr = PetscViewerFlush(subviewer);CHKERRQ(ierr);
234e9627c49SStefano Zampini     }
235e9627c49SStefano Zampini     ierr = PetscViewerRestoreSubViewer(viewer,PetscSubcommChild(subcomm),&subviewer);CHKERRQ(ierr);
236e9627c49SStefano Zampini     ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
237e9627c49SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
238e9627c49SStefano Zampini   }
2396b78500eSPatrick Sanan   PetscFunctionReturn(0);
2406b78500eSPatrick Sanan }
241a13144ffSStefano Zampini 
2421e0482f5SStefano Zampini static PetscErrorCode PCBDDCSetDiscreteGradient_BDDC(PC pc, Mat G, PetscInt order, PetscInt field, PetscBool global, PetscBool conforming)
243a13144ffSStefano Zampini {
244a13144ffSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
245a13144ffSStefano Zampini   PetscErrorCode ierr;
246a13144ffSStefano Zampini 
247a13144ffSStefano Zampini   PetscFunctionBegin;
248a13144ffSStefano Zampini   ierr = PetscObjectReference((PetscObject)G);CHKERRQ(ierr);
249a13144ffSStefano Zampini   ierr = MatDestroy(&pcbddc->discretegradient);CHKERRQ(ierr);
250a13144ffSStefano Zampini   pcbddc->discretegradient = G;
251a13144ffSStefano Zampini   pcbddc->nedorder         = order > 0 ? order : -order;
252495a2a07SStefano Zampini   pcbddc->nedfield         = field;
2531e0482f5SStefano Zampini   pcbddc->nedglobal        = global;
2541e0482f5SStefano Zampini   pcbddc->conforming       = conforming;
255a13144ffSStefano Zampini   PetscFunctionReturn(0);
256a13144ffSStefano Zampini }
257a13144ffSStefano Zampini 
258a13144ffSStefano Zampini /*@
259a13144ffSStefano Zampini  PCBDDCSetDiscreteGradient - Sets the discrete gradient
260a13144ffSStefano Zampini 
261a13144ffSStefano Zampini    Collective on PC
262a13144ffSStefano Zampini 
263a13144ffSStefano Zampini    Input Parameters:
264a13144ffSStefano Zampini +  pc         - the preconditioning context
265a13144ffSStefano Zampini .  G          - the discrete gradient matrix (should be in AIJ format)
266a13144ffSStefano Zampini .  order      - the order of the Nedelec space (1 for the lowest order)
267495a2a07SStefano Zampini .  field      - the field id of the Nedelec dofs (not used if the fields have not been specified)
2681e0482f5SStefano Zampini .  global     - the type of global ordering for the rows of G
269a13144ffSStefano Zampini -  conforming - whether the mesh is conforming or not
270a13144ffSStefano Zampini 
271a13144ffSStefano Zampini    Level: advanced
272a13144ffSStefano Zampini 
273a13144ffSStefano Zampini    Notes: The discrete gradient matrix G is used to analyze the subdomain edges, and it should not contain any zero entry.
274495a2a07SStefano Zampini           For variable order spaces, the order should be set to zero.
2751e0482f5SStefano Zampini           If global is true, the rows of G should be given in global ordering for the whole dofs;
2761e0482f5SStefano Zampini           if false, the ordering should be global for the Nedelec field.
2771e0482f5SStefano Zampini           In the latter case, it should hold gid[i] < gid[j] iff geid[i] < geid[j], with gid the global orderding for all the dofs
2781e0482f5SStefano Zampini           and geid the one for the Nedelec field.
279a13144ffSStefano Zampini 
280495a2a07SStefano Zampini .seealso: PCBDDC,PCBDDCSetDofsSplitting(),PCBDDCSetDofsSplittingLocal()
281a13144ffSStefano Zampini @*/
2821e0482f5SStefano Zampini PetscErrorCode PCBDDCSetDiscreteGradient(PC pc, Mat G, PetscInt order, PetscInt field, PetscBool global, PetscBool conforming)
283a13144ffSStefano Zampini {
284a13144ffSStefano Zampini   PetscErrorCode ierr;
285a13144ffSStefano Zampini 
286a13144ffSStefano Zampini   PetscFunctionBegin;
287a13144ffSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
288a13144ffSStefano Zampini   PetscValidHeaderSpecific(G,MAT_CLASSID,2);
289a13144ffSStefano Zampini   PetscValidLogicalCollectiveInt(pc,order,3);
2901e0482f5SStefano Zampini   PetscValidLogicalCollectiveInt(pc,field,4);
2911e0482f5SStefano Zampini   PetscValidLogicalCollectiveBool(pc,global,5);
2921e0482f5SStefano Zampini   PetscValidLogicalCollectiveBool(pc,conforming,6);
2931e0482f5SStefano Zampini   PetscCheckSameComm(pc,1,G,2);
2941e0482f5SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDiscreteGradient_C",(PC,Mat,PetscInt,PetscInt,PetscBool,PetscBool),(pc,G,order,field,global,conforming));CHKERRQ(ierr);
295a13144ffSStefano Zampini   PetscFunctionReturn(0);
296a13144ffSStefano Zampini }
297a13144ffSStefano Zampini 
2988ae0ca82SStefano Zampini static PetscErrorCode PCBDDCSetDivergenceMat_BDDC(PC pc, Mat divudotp, PetscBool trans, IS vl2l)
299a198735bSStefano Zampini {
300a198735bSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
301a198735bSStefano Zampini   PetscErrorCode ierr;
3026b78500eSPatrick Sanan 
303a198735bSStefano Zampini   PetscFunctionBegin;
304a198735bSStefano Zampini   ierr = PetscObjectReference((PetscObject)divudotp);CHKERRQ(ierr);
305a198735bSStefano Zampini   ierr = MatDestroy(&pcbddc->divudotp);CHKERRQ(ierr);
306a198735bSStefano Zampini   pcbddc->divudotp = divudotp;
3078ae0ca82SStefano Zampini   pcbddc->divudotp_trans = trans;
308a198735bSStefano Zampini   pcbddc->compute_nonetflux = PETSC_TRUE;
309a198735bSStefano Zampini   if (vl2l) {
310a198735bSStefano Zampini     ierr = PetscObjectReference((PetscObject)vl2l);CHKERRQ(ierr);
311fa23a32eSStefano Zampini     ierr = ISDestroy(&pcbddc->divudotp_vl2l);CHKERRQ(ierr);
312a198735bSStefano Zampini     pcbddc->divudotp_vl2l = vl2l;
313a198735bSStefano Zampini   }
314a198735bSStefano Zampini   PetscFunctionReturn(0);
315a198735bSStefano Zampini }
3163d996552SStefano Zampini 
317a198735bSStefano Zampini /*@
318a198735bSStefano Zampini  PCBDDCSetDivergenceMat - Sets the linear operator representing \int_\Omega \div {\bf u} \cdot p dx
319a198735bSStefano Zampini 
320a198735bSStefano Zampini    Collective on PC
321a198735bSStefano Zampini 
322a198735bSStefano Zampini    Input Parameters:
323a198735bSStefano Zampini +  pc - the preconditioning context
324a198735bSStefano Zampini .  divudotp - the matrix (must be of type MATIS)
3258ae0ca82SStefano Zampini .  trans - if trans if false (resp. true), then pressures are in the test (trial) space and velocities are in the trial (test) space.
32605a3bf82SStefano Zampini -  vl2l - optional index set describing the local (wrt the local matrix in divudotp) to local (wrt the local matrix in the preconditioning matrix) map for the velocities
327a198735bSStefano Zampini 
328a198735bSStefano Zampini    Level: advanced
329a198735bSStefano Zampini 
3308ae0ca82SStefano Zampini    Notes: This auxiliary matrix is used to compute quadrature weights representing the net-flux across subdomain boundaries
33105a3bf82SStefano Zampini           If vl2l is NULL, the local ordering for velocities in divudotp should match that of the preconditioning matrix
332a198735bSStefano Zampini 
333a198735bSStefano Zampini .seealso: PCBDDC
334a198735bSStefano Zampini @*/
3358ae0ca82SStefano Zampini PetscErrorCode PCBDDCSetDivergenceMat(PC pc, Mat divudotp, PetscBool trans, IS vl2l)
336a198735bSStefano Zampini {
337a198735bSStefano Zampini   PetscBool      ismatis;
338a198735bSStefano Zampini   PetscErrorCode ierr;
339a198735bSStefano Zampini 
340a198735bSStefano Zampini   PetscFunctionBegin;
341a198735bSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
342a198735bSStefano Zampini   PetscValidHeaderSpecific(divudotp,MAT_CLASSID,2);
343a198735bSStefano Zampini   PetscCheckSameComm(pc,1,divudotp,2);
3448ae0ca82SStefano Zampini   PetscValidLogicalCollectiveBool(pc,trans,3);
3458ae0ca82SStefano Zampini   if (vl2l) PetscValidHeaderSpecific(divudotp,IS_CLASSID,4);
346a198735bSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)divudotp,MATIS,&ismatis);CHKERRQ(ierr);
347a198735bSStefano Zampini   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Divergence matrix needs to be of type MATIS");
3488ae0ca82SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDivergenceMat_C",(PC,Mat,PetscBool,IS),(pc,divudotp,trans,vl2l));CHKERRQ(ierr);
349a198735bSStefano Zampini   PetscFunctionReturn(0);
350a198735bSStefano Zampini }
3512d505d7fSStefano Zampini 
3521dd7afcfSStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change, PetscBool interior)
353b9b85e73SStefano Zampini {
354b9b85e73SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
355b9b85e73SStefano Zampini   PetscErrorCode ierr;
356b9b85e73SStefano Zampini 
357b9b85e73SStefano Zampini   PetscFunctionBegin;
358b9b85e73SStefano Zampini   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
35956282151SStefano Zampini   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
360b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix = change;
3611dd7afcfSStefano Zampini   pcbddc->change_interior = interior;
362b9b85e73SStefano Zampini   PetscFunctionReturn(0);
363b9b85e73SStefano Zampini }
364b9b85e73SStefano Zampini /*@
365906d46d4SStefano Zampini  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
366b9b85e73SStefano Zampini 
367b9b85e73SStefano Zampini    Collective on PC
368b9b85e73SStefano Zampini 
369b9b85e73SStefano Zampini    Input Parameters:
370b9b85e73SStefano Zampini +  pc - the preconditioning context
3711dd7afcfSStefano Zampini .  change - the change of basis matrix
3721dd7afcfSStefano Zampini -  interior - whether or not the change of basis modifies interior dofs
373b9b85e73SStefano Zampini 
374b9b85e73SStefano Zampini    Level: intermediate
375b9b85e73SStefano Zampini 
376b9b85e73SStefano Zampini    Notes:
377b9b85e73SStefano Zampini 
378b9b85e73SStefano Zampini .seealso: PCBDDC
379b9b85e73SStefano Zampini @*/
3801dd7afcfSStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change, PetscBool interior)
381b9b85e73SStefano Zampini {
382b9b85e73SStefano Zampini   PetscErrorCode ierr;
383b9b85e73SStefano Zampini 
384b9b85e73SStefano Zampini   PetscFunctionBegin;
385b9b85e73SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
386b9b85e73SStefano Zampini   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
387906d46d4SStefano Zampini   PetscCheckSameComm(pc,1,change,2);
388906d46d4SStefano Zampini   if (pc->mat) {
389906d46d4SStefano Zampini     PetscInt rows_c,cols_c,rows,cols;
390906d46d4SStefano Zampini     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
391906d46d4SStefano Zampini     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
392e0fe2d75SToby Isaac     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);
393e0fe2d75SToby Isaac     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);
394906d46d4SStefano Zampini     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
395906d46d4SStefano Zampini     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
396e0fe2d75SToby Isaac     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);
397e0fe2d75SToby Isaac     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);
398906d46d4SStefano Zampini   }
3991dd7afcfSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat,PetscBool),(pc,change,interior));CHKERRQ(ierr);
400b9b85e73SStefano Zampini   PetscFunctionReturn(0);
401b9b85e73SStefano Zampini }
4022d505d7fSStefano Zampini 
40330368db7SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesIS_BDDC(PC pc, IS PrimalVertices)
40430368db7SStefano Zampini {
40530368db7SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
40656282151SStefano Zampini   PetscBool      isequal = PETSC_FALSE;
40730368db7SStefano Zampini   PetscErrorCode ierr;
40830368db7SStefano Zampini 
40930368db7SStefano Zampini   PetscFunctionBegin;
41056282151SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
41156282151SStefano Zampini   if (pcbddc->user_primal_vertices) {
41256282151SStefano Zampini     ierr = ISEqual(PrimalVertices,pcbddc->user_primal_vertices,&isequal);CHKERRQ(ierr);
41356282151SStefano Zampini   }
41430368db7SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
41530368db7SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices_local);CHKERRQ(ierr);
41630368db7SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
41756282151SStefano Zampini   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
41830368db7SStefano Zampini   PetscFunctionReturn(0);
41930368db7SStefano Zampini }
420ab8c8b98SStefano Zampini 
42130368db7SStefano Zampini /*@
42230368db7SStefano Zampini  PCBDDCSetPrimalVerticesIS - Set additional user defined primal vertices in PCBDDC
42330368db7SStefano Zampini 
42430368db7SStefano Zampini    Collective
42530368db7SStefano Zampini 
42630368db7SStefano Zampini    Input Parameters:
42730368db7SStefano Zampini +  pc - the preconditioning context
42830368db7SStefano Zampini -  PrimalVertices - index set of primal vertices in global numbering (can be empty)
42930368db7SStefano Zampini 
43030368db7SStefano Zampini    Level: intermediate
43130368db7SStefano Zampini 
43230368db7SStefano Zampini    Notes:
43330368db7SStefano Zampini      Any process can list any global node
43430368db7SStefano Zampini 
43530368db7SStefano Zampini .seealso: PCBDDC
43630368db7SStefano Zampini @*/
43730368db7SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesIS(PC pc, IS PrimalVertices)
43830368db7SStefano Zampini {
43930368db7SStefano Zampini   PetscErrorCode ierr;
44030368db7SStefano Zampini 
44130368db7SStefano Zampini   PetscFunctionBegin;
44230368db7SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
44330368db7SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
44430368db7SStefano Zampini   PetscCheckSameComm(pc,1,PrimalVertices,2);
44530368db7SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
44630368db7SStefano Zampini   PetscFunctionReturn(0);
44730368db7SStefano Zampini }
4482d505d7fSStefano Zampini 
449674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
450674ae819SStefano Zampini {
451674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
45256282151SStefano Zampini   PetscBool      isequal = PETSC_FALSE;
453674ae819SStefano Zampini   PetscErrorCode ierr;
4541e6b0712SBarry Smith 
455674ae819SStefano Zampini   PetscFunctionBegin;
45656282151SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
45756282151SStefano Zampini   if (pcbddc->user_primal_vertices_local) {
45856282151SStefano Zampini     ierr = ISEqual(PrimalVertices,pcbddc->user_primal_vertices_local,&isequal);CHKERRQ(ierr);
45956282151SStefano Zampini   }
460674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
46130368db7SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices_local);CHKERRQ(ierr);
46230368db7SStefano Zampini   pcbddc->user_primal_vertices_local = PrimalVertices;
46356282151SStefano Zampini   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
464674ae819SStefano Zampini   PetscFunctionReturn(0);
465674ae819SStefano Zampini }
466674ae819SStefano Zampini /*@
46728509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
468674ae819SStefano Zampini 
46917eb1463SStefano Zampini    Collective
470674ae819SStefano Zampini 
471674ae819SStefano Zampini    Input Parameters:
472674ae819SStefano Zampini +  pc - the preconditioning context
47317eb1463SStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering (can be empty)
474674ae819SStefano Zampini 
475674ae819SStefano Zampini    Level: intermediate
476674ae819SStefano Zampini 
477674ae819SStefano Zampini    Notes:
478674ae819SStefano Zampini 
479674ae819SStefano Zampini .seealso: PCBDDC
480674ae819SStefano Zampini @*/
481674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
482674ae819SStefano Zampini {
483674ae819SStefano Zampini   PetscErrorCode ierr;
484674ae819SStefano Zampini 
485674ae819SStefano Zampini   PetscFunctionBegin;
486674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
487674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
48817eb1463SStefano Zampini   PetscCheckSameComm(pc,1,PrimalVertices,2);
489674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
490674ae819SStefano Zampini   PetscFunctionReturn(0);
491674ae819SStefano Zampini }
4922d505d7fSStefano Zampini 
4934fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
4944fad6a16SStefano Zampini {
4954fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4964fad6a16SStefano Zampini 
4974fad6a16SStefano Zampini   PetscFunctionBegin;
4984fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
4994fad6a16SStefano Zampini   PetscFunctionReturn(0);
5004fad6a16SStefano Zampini }
5011e6b0712SBarry Smith 
5024fad6a16SStefano Zampini /*@
50328509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
5044fad6a16SStefano Zampini 
5054fad6a16SStefano Zampini    Logically collective on PC
5064fad6a16SStefano Zampini 
5074fad6a16SStefano Zampini    Input Parameters:
5084fad6a16SStefano Zampini +  pc - the preconditioning context
50928509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
5104fad6a16SStefano Zampini 
5110f202f7eSStefano Zampini    Options Database Keys:
5120f202f7eSStefano Zampini .    -pc_bddc_coarsening_ratio
5134fad6a16SStefano Zampini 
5144fad6a16SStefano Zampini    Level: intermediate
5154fad6a16SStefano Zampini 
5164fad6a16SStefano Zampini    Notes:
5170f202f7eSStefano Zampini      Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
5184fad6a16SStefano Zampini 
5190f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetLevels()
5204fad6a16SStefano Zampini @*/
5214fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
5224fad6a16SStefano Zampini {
5234fad6a16SStefano Zampini   PetscErrorCode ierr;
5244fad6a16SStefano Zampini 
5254fad6a16SStefano Zampini   PetscFunctionBegin;
5264fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5272b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
5284fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
5294fad6a16SStefano Zampini   PetscFunctionReturn(0);
5304fad6a16SStefano Zampini }
5312b510759SStefano Zampini 
532b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
533b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
534b8ffe317SStefano Zampini {
535b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
536b8ffe317SStefano Zampini 
537b8ffe317SStefano Zampini   PetscFunctionBegin;
53885c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
539b8ffe317SStefano Zampini   PetscFunctionReturn(0);
540b8ffe317SStefano Zampini }
541b8ffe317SStefano Zampini 
542b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
5432b510759SStefano Zampini {
5442b510759SStefano Zampini   PetscErrorCode ierr;
5452b510759SStefano Zampini 
5462b510759SStefano Zampini   PetscFunctionBegin;
5472b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
548b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
549b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
5502b510759SStefano Zampini   PetscFunctionReturn(0);
5512b510759SStefano Zampini }
5521e6b0712SBarry Smith 
5532b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
5544fad6a16SStefano Zampini {
5554fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5564fad6a16SStefano Zampini 
5574fad6a16SStefano Zampini   PetscFunctionBegin;
5582b510759SStefano Zampini   pcbddc->current_level = level;
5594fad6a16SStefano Zampini   PetscFunctionReturn(0);
5604fad6a16SStefano Zampini }
5611e6b0712SBarry Smith 
562b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
563b8ffe317SStefano Zampini {
564b8ffe317SStefano Zampini   PetscErrorCode ierr;
565b8ffe317SStefano Zampini 
566b8ffe317SStefano Zampini   PetscFunctionBegin;
567b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
568b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
569b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
570b8ffe317SStefano Zampini   PetscFunctionReturn(0);
571b8ffe317SStefano Zampini }
572b8ffe317SStefano Zampini 
5732b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
5742b510759SStefano Zampini {
5752b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5762b510759SStefano Zampini 
5772b510759SStefano Zampini   PetscFunctionBegin;
5782b510759SStefano Zampini   pcbddc->max_levels = levels;
5792b510759SStefano Zampini   PetscFunctionReturn(0);
5802b510759SStefano Zampini }
5812b510759SStefano Zampini 
5824fad6a16SStefano Zampini /*@
58328509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
5844fad6a16SStefano Zampini 
5854fad6a16SStefano Zampini    Logically collective on PC
5864fad6a16SStefano Zampini 
5874fad6a16SStefano Zampini    Input Parameters:
5884fad6a16SStefano Zampini +  pc - the preconditioning context
58928509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
5904fad6a16SStefano Zampini 
5910f202f7eSStefano Zampini    Options Database Keys:
5920f202f7eSStefano Zampini .    -pc_bddc_levels
5934fad6a16SStefano Zampini 
5944fad6a16SStefano Zampini    Level: intermediate
5954fad6a16SStefano Zampini 
5964fad6a16SStefano Zampini    Notes:
5970f202f7eSStefano Zampini      Default value is 0, i.e. traditional one-level BDDC
5984fad6a16SStefano Zampini 
5990f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
6004fad6a16SStefano Zampini @*/
6012b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
6024fad6a16SStefano Zampini {
6034fad6a16SStefano Zampini   PetscErrorCode ierr;
6044fad6a16SStefano Zampini 
6054fad6a16SStefano Zampini   PetscFunctionBegin;
6064fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
6072b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
608312be037SStefano Zampini   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
6092b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
6104fad6a16SStefano Zampini   PetscFunctionReturn(0);
6114fad6a16SStefano Zampini }
6121e6b0712SBarry Smith 
6133b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
6143b03a366Sstefano_zampini {
6153b03a366Sstefano_zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
61656282151SStefano Zampini   PetscBool      isequal = PETSC_FALSE;
6173b03a366Sstefano_zampini   PetscErrorCode ierr;
6183b03a366Sstefano_zampini 
6193b03a366Sstefano_zampini   PetscFunctionBegin;
62056282151SStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
62156282151SStefano Zampini   if (pcbddc->DirichletBoundaries) {
62256282151SStefano Zampini     ierr = ISEqual(DirichletBoundaries,pcbddc->DirichletBoundaries,&isequal);CHKERRQ(ierr);
62356282151SStefano Zampini   }
624785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
625785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
6263b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
62736e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
62856282151SStefano Zampini   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
6293b03a366Sstefano_zampini   PetscFunctionReturn(0);
6303b03a366Sstefano_zampini }
6311e6b0712SBarry Smith 
6323b03a366Sstefano_zampini /*@
63328509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
6343b03a366Sstefano_zampini 
635785d1243SStefano Zampini    Collective
6363b03a366Sstefano_zampini 
6373b03a366Sstefano_zampini    Input Parameters:
6383b03a366Sstefano_zampini +  pc - the preconditioning context
639785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
6403b03a366Sstefano_zampini 
6413b03a366Sstefano_zampini    Level: intermediate
6423b03a366Sstefano_zampini 
6430f202f7eSStefano Zampini    Notes:
6440f202f7eSStefano Zampini      Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node
6453b03a366Sstefano_zampini 
6460f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
6473b03a366Sstefano_zampini @*/
6483b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
6493b03a366Sstefano_zampini {
6503b03a366Sstefano_zampini   PetscErrorCode ierr;
6513b03a366Sstefano_zampini 
6523b03a366Sstefano_zampini   PetscFunctionBegin;
6533b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
654674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
655785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
6563b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
6573b03a366Sstefano_zampini   PetscFunctionReturn(0);
6583b03a366Sstefano_zampini }
6591e6b0712SBarry Smith 
66082ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
6613b03a366Sstefano_zampini {
6623b03a366Sstefano_zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
66356282151SStefano Zampini   PetscBool      isequal = PETSC_FALSE;
6643b03a366Sstefano_zampini   PetscErrorCode ierr;
6653b03a366Sstefano_zampini 
6663b03a366Sstefano_zampini   PetscFunctionBegin;
66756282151SStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
66856282151SStefano Zampini   if (pcbddc->DirichletBoundariesLocal) {
66956282151SStefano Zampini     ierr = ISEqual(DirichletBoundaries,pcbddc->DirichletBoundariesLocal,&isequal);CHKERRQ(ierr);
67056282151SStefano Zampini   }
671785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
672785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
6733b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
674785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
67556282151SStefano Zampini   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
6763b03a366Sstefano_zampini   PetscFunctionReturn(0);
6773b03a366Sstefano_zampini }
6783b03a366Sstefano_zampini 
6793b03a366Sstefano_zampini /*@
68082ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
6813b03a366Sstefano_zampini 
682785d1243SStefano Zampini    Collective
6833b03a366Sstefano_zampini 
6843b03a366Sstefano_zampini    Input Parameters:
6853b03a366Sstefano_zampini +  pc - the preconditioning context
68682ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
6873b03a366Sstefano_zampini 
6883b03a366Sstefano_zampini    Level: intermediate
6893b03a366Sstefano_zampini 
6903b03a366Sstefano_zampini    Notes:
6913b03a366Sstefano_zampini 
6920f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
6933b03a366Sstefano_zampini @*/
69482ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
6953b03a366Sstefano_zampini {
6963b03a366Sstefano_zampini   PetscErrorCode ierr;
6973b03a366Sstefano_zampini 
6983b03a366Sstefano_zampini   PetscFunctionBegin;
6993b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
7003b03a366Sstefano_zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
70182ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
70282ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
7033b03a366Sstefano_zampini   PetscFunctionReturn(0);
7043b03a366Sstefano_zampini }
7053b03a366Sstefano_zampini 
70653cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
7070c7d97c5SJed Brown {
7080c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
70956282151SStefano Zampini   PetscBool      isequal = PETSC_FALSE;
71053cdbc3dSStefano Zampini   PetscErrorCode ierr;
7110c7d97c5SJed Brown 
7120c7d97c5SJed Brown   PetscFunctionBegin;
71356282151SStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
71456282151SStefano Zampini   if (pcbddc->NeumannBoundaries) {
71556282151SStefano Zampini     ierr = ISEqual(NeumannBoundaries,pcbddc->NeumannBoundaries,&isequal);CHKERRQ(ierr);
71656282151SStefano Zampini   }
717785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
718785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
71953cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
72036e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
72156282151SStefano Zampini   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
7220c7d97c5SJed Brown   PetscFunctionReturn(0);
7230c7d97c5SJed Brown }
7241e6b0712SBarry Smith 
72557527edcSJed Brown /*@
72628509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
72757527edcSJed Brown 
728785d1243SStefano Zampini    Collective
72957527edcSJed Brown 
73057527edcSJed Brown    Input Parameters:
73157527edcSJed Brown +  pc - the preconditioning context
732785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
73357527edcSJed Brown 
73457527edcSJed Brown    Level: intermediate
73557527edcSJed Brown 
7360f202f7eSStefano Zampini    Notes:
7370f202f7eSStefano Zampini      Any process can list any global node
73857527edcSJed Brown 
7390f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
74057527edcSJed Brown @*/
74153cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
7420c7d97c5SJed Brown {
7430c7d97c5SJed Brown   PetscErrorCode ierr;
7440c7d97c5SJed Brown 
7450c7d97c5SJed Brown   PetscFunctionBegin;
7460c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
747674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
748785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
74953cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
75053cdbc3dSStefano Zampini   PetscFunctionReturn(0);
75153cdbc3dSStefano Zampini }
7521e6b0712SBarry Smith 
75382ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
7540c7d97c5SJed Brown {
7550c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
75656282151SStefano Zampini   PetscBool      isequal = PETSC_FALSE;
7570c7d97c5SJed Brown   PetscErrorCode ierr;
7580c7d97c5SJed Brown 
7590c7d97c5SJed Brown   PetscFunctionBegin;
76056282151SStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
76156282151SStefano Zampini   if (pcbddc->NeumannBoundariesLocal) {
76256282151SStefano Zampini     ierr = ISEqual(NeumannBoundaries,pcbddc->NeumannBoundariesLocal,&isequal);CHKERRQ(ierr);
76356282151SStefano Zampini   }
764785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
765785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
7660c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
767785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
76856282151SStefano Zampini   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
7690c7d97c5SJed Brown   PetscFunctionReturn(0);
7700c7d97c5SJed Brown }
7710c7d97c5SJed Brown 
7720c7d97c5SJed Brown /*@
77382ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
7740c7d97c5SJed Brown 
775785d1243SStefano Zampini    Collective
7760c7d97c5SJed Brown 
7770c7d97c5SJed Brown    Input Parameters:
7780c7d97c5SJed Brown +  pc - the preconditioning context
77982ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
7800c7d97c5SJed Brown 
7810c7d97c5SJed Brown    Level: intermediate
7820c7d97c5SJed Brown 
7830c7d97c5SJed Brown    Notes:
7840c7d97c5SJed Brown 
7850f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
7860c7d97c5SJed Brown @*/
78782ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
7880c7d97c5SJed Brown {
7890c7d97c5SJed Brown   PetscErrorCode ierr;
7900c7d97c5SJed Brown 
7910c7d97c5SJed Brown   PetscFunctionBegin;
7920c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
7930c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
79482ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
79582ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
79653cdbc3dSStefano Zampini   PetscFunctionReturn(0);
79753cdbc3dSStefano Zampini }
79853cdbc3dSStefano Zampini 
799da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
800da1bb401SStefano Zampini {
801da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
802da1bb401SStefano Zampini 
803da1bb401SStefano Zampini   PetscFunctionBegin;
804da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
805da1bb401SStefano Zampini   PetscFunctionReturn(0);
806da1bb401SStefano Zampini }
8071e6b0712SBarry Smith 
808da1bb401SStefano Zampini /*@
809785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
810da1bb401SStefano Zampini 
811785d1243SStefano Zampini    Collective
812785d1243SStefano Zampini 
813785d1243SStefano Zampini    Input Parameters:
814785d1243SStefano Zampini .  pc - the preconditioning context
815785d1243SStefano Zampini 
816785d1243SStefano Zampini    Output Parameters:
817785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
818785d1243SStefano Zampini 
819785d1243SStefano Zampini    Level: intermediate
820785d1243SStefano Zampini 
8210f202f7eSStefano Zampini    Notes:
8220f202f7eSStefano Zampini      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
823785d1243SStefano Zampini 
824785d1243SStefano Zampini .seealso: PCBDDC
825785d1243SStefano Zampini @*/
826785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
827785d1243SStefano Zampini {
828785d1243SStefano Zampini   PetscErrorCode ierr;
829785d1243SStefano Zampini 
830785d1243SStefano Zampini   PetscFunctionBegin;
831785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
832785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
833785d1243SStefano Zampini   PetscFunctionReturn(0);
834785d1243SStefano Zampini }
835785d1243SStefano Zampini 
836785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
837785d1243SStefano Zampini {
838785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
839785d1243SStefano Zampini 
840785d1243SStefano Zampini   PetscFunctionBegin;
841785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
842785d1243SStefano Zampini   PetscFunctionReturn(0);
843785d1243SStefano Zampini }
844785d1243SStefano Zampini 
845da1bb401SStefano Zampini /*@
84682ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
847da1bb401SStefano Zampini 
848785d1243SStefano Zampini    Collective
849da1bb401SStefano Zampini 
850da1bb401SStefano Zampini    Input Parameters:
85128509bceSStefano Zampini .  pc - the preconditioning context
852da1bb401SStefano Zampini 
853da1bb401SStefano Zampini    Output Parameters:
85428509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
855da1bb401SStefano Zampini 
856da1bb401SStefano Zampini    Level: intermediate
857da1bb401SStefano Zampini 
858da1bb401SStefano Zampini    Notes:
8590f202f7eSStefano Zampini      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).
8600f202f7eSStefano Zampini           In the latter case, the IS will be available after PCSetUp.
861da1bb401SStefano Zampini 
862da1bb401SStefano Zampini .seealso: PCBDDC
863da1bb401SStefano Zampini @*/
86482ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
865da1bb401SStefano Zampini {
866da1bb401SStefano Zampini   PetscErrorCode ierr;
867da1bb401SStefano Zampini 
868da1bb401SStefano Zampini   PetscFunctionBegin;
869da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
87082ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
871da1bb401SStefano Zampini   PetscFunctionReturn(0);
872da1bb401SStefano Zampini }
8731e6b0712SBarry Smith 
87453cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
87553cdbc3dSStefano Zampini {
87653cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
87753cdbc3dSStefano Zampini 
87853cdbc3dSStefano Zampini   PetscFunctionBegin;
87953cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
88053cdbc3dSStefano Zampini   PetscFunctionReturn(0);
88153cdbc3dSStefano Zampini }
8821e6b0712SBarry Smith 
88353cdbc3dSStefano Zampini /*@
884785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
88553cdbc3dSStefano Zampini 
886785d1243SStefano Zampini    Collective
887785d1243SStefano Zampini 
888785d1243SStefano Zampini    Input Parameters:
889785d1243SStefano Zampini .  pc - the preconditioning context
890785d1243SStefano Zampini 
891785d1243SStefano Zampini    Output Parameters:
892785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
893785d1243SStefano Zampini 
894785d1243SStefano Zampini    Level: intermediate
895785d1243SStefano Zampini 
8960f202f7eSStefano Zampini    Notes:
8970f202f7eSStefano Zampini      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
898785d1243SStefano Zampini 
899785d1243SStefano Zampini .seealso: PCBDDC
900785d1243SStefano Zampini @*/
901785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
902785d1243SStefano Zampini {
903785d1243SStefano Zampini   PetscErrorCode ierr;
904785d1243SStefano Zampini 
905785d1243SStefano Zampini   PetscFunctionBegin;
906785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
907785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
908785d1243SStefano Zampini   PetscFunctionReturn(0);
909785d1243SStefano Zampini }
910785d1243SStefano Zampini 
911785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
912785d1243SStefano Zampini {
913785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
914785d1243SStefano Zampini 
915785d1243SStefano Zampini   PetscFunctionBegin;
916785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
917785d1243SStefano Zampini   PetscFunctionReturn(0);
918785d1243SStefano Zampini }
919785d1243SStefano Zampini 
92053cdbc3dSStefano Zampini /*@
92182ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
92253cdbc3dSStefano Zampini 
923785d1243SStefano Zampini    Collective
92453cdbc3dSStefano Zampini 
92553cdbc3dSStefano Zampini    Input Parameters:
92628509bceSStefano Zampini .  pc - the preconditioning context
92753cdbc3dSStefano Zampini 
92853cdbc3dSStefano Zampini    Output Parameters:
92928509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
93053cdbc3dSStefano Zampini 
93153cdbc3dSStefano Zampini    Level: intermediate
93253cdbc3dSStefano Zampini 
93353cdbc3dSStefano Zampini    Notes:
9340f202f7eSStefano Zampini      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).
9350f202f7eSStefano Zampini           In the latter case, the IS will be available after PCSetUp.
93653cdbc3dSStefano Zampini 
93753cdbc3dSStefano Zampini .seealso: PCBDDC
93853cdbc3dSStefano Zampini @*/
93982ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
94053cdbc3dSStefano Zampini {
94153cdbc3dSStefano Zampini   PetscErrorCode ierr;
94253cdbc3dSStefano Zampini 
94353cdbc3dSStefano Zampini   PetscFunctionBegin;
94453cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
94582ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
9460c7d97c5SJed Brown   PetscFunctionReturn(0);
9470c7d97c5SJed Brown }
9481e6b0712SBarry Smith 
9491a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
95036e030ebSStefano Zampini {
95136e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
952da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
95356282151SStefano Zampini   PetscBool      same_data = PETSC_FALSE;
954da1bb401SStefano Zampini   PetscErrorCode ierr;
95536e030ebSStefano Zampini 
95636e030ebSStefano Zampini   PetscFunctionBegin;
9578687889aSStefano Zampini   if (!nvtxs) {
95804194a47SStefano Zampini     if (copymode == PETSC_OWN_POINTER) {
95904194a47SStefano Zampini       ierr = PetscFree(xadj);CHKERRQ(ierr);
96004194a47SStefano Zampini       ierr = PetscFree(adjncy);CHKERRQ(ierr);
96104194a47SStefano Zampini     }
9628687889aSStefano Zampini     ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
9638687889aSStefano Zampini     PetscFunctionReturn(0);
9648687889aSStefano Zampini   }
96566da6bd7Sstefano_zampini   if (mat_graph->nvtxs == nvtxs && mat_graph->freecsr) { /* we own the data */
96656282151SStefano Zampini     if (mat_graph->xadj == xadj && mat_graph->adjncy == adjncy) same_data = PETSC_TRUE;
96756282151SStefano Zampini     if (!same_data && mat_graph->xadj[nvtxs] == xadj[nvtxs]) {
9682d505d7fSStefano Zampini       ierr = PetscMemcmp(xadj,mat_graph->xadj,(nvtxs+1)*sizeof(PetscInt),&same_data);CHKERRQ(ierr);
9692d505d7fSStefano Zampini       if (same_data) {
9702d505d7fSStefano Zampini         ierr = PetscMemcmp(adjncy,mat_graph->adjncy,xadj[nvtxs]*sizeof(PetscInt),&same_data);CHKERRQ(ierr);
9712d505d7fSStefano Zampini       }
97256282151SStefano Zampini     }
97356282151SStefano Zampini   }
97456282151SStefano Zampini   if (!same_data) {
975674ae819SStefano Zampini     /* free old CSR */
976674ae819SStefano Zampini     ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
977674ae819SStefano Zampini     /* get CSR into graph structure */
978da1bb401SStefano Zampini     if (copymode == PETSC_COPY_VALUES) {
979854ce69bSBarry Smith       ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr);
980785e854fSJed Brown       ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
981674ae819SStefano Zampini       ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
982674ae819SStefano Zampini       ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
983a1dbd327SStefano Zampini       mat_graph->freecsr = PETSC_TRUE;
984da1bb401SStefano Zampini     } else if (copymode == PETSC_OWN_POINTER) {
9851a83f524SJed Brown       mat_graph->xadj    = (PetscInt*)xadj;
9861a83f524SJed Brown       mat_graph->adjncy  = (PetscInt*)adjncy;
987a1dbd327SStefano Zampini       mat_graph->freecsr = PETSC_TRUE;
988a1dbd327SStefano Zampini     } else if (copymode == PETSC_USE_POINTER) {
989a1dbd327SStefano Zampini       mat_graph->xadj    = (PetscInt*)xadj;
990a1dbd327SStefano Zampini       mat_graph->adjncy  = (PetscInt*)adjncy;
991a1dbd327SStefano Zampini       mat_graph->freecsr = PETSC_FALSE;
992e0fe2d75SToby Isaac     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %D",copymode);
993575ad6abSStefano Zampini     mat_graph->nvtxs_csr = nvtxs;
99456282151SStefano Zampini     pcbddc->recompute_topography = PETSC_TRUE;
99556282151SStefano Zampini   }
99636e030ebSStefano Zampini   PetscFunctionReturn(0);
99736e030ebSStefano Zampini }
9981e6b0712SBarry Smith 
99936e030ebSStefano Zampini /*@
100054fffbccSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local degrees of freedom.
100136e030ebSStefano Zampini 
100236e030ebSStefano Zampini    Not collective
100336e030ebSStefano Zampini 
100436e030ebSStefano Zampini    Input Parameters:
100554fffbccSStefano Zampini +  pc - the preconditioning context.
100654fffbccSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the number of local dofs).
100754fffbccSStefano Zampini .  xadj, adjncy - the connectivity of the dofs in CSR format.
100854fffbccSStefano Zampini -  copymode - supported modes are PETSC_COPY_VALUES, PETSC_USE_POINTER or PETSC_OWN_POINTER.
100936e030ebSStefano Zampini 
101036e030ebSStefano Zampini    Level: intermediate
101136e030ebSStefano Zampini 
101254fffbccSStefano Zampini    Notes: A dof is considered connected with all local dofs if xadj[dof+1]-xadj[dof] == 1 and adjncy[xadj[dof]] is negative.
101336e030ebSStefano Zampini 
101428509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
101536e030ebSStefano Zampini @*/
10161a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
101736e030ebSStefano Zampini {
1018575ad6abSStefano Zampini   void (*f)(void) = 0;
101936e030ebSStefano Zampini   PetscErrorCode ierr;
102036e030ebSStefano Zampini 
102136e030ebSStefano Zampini   PetscFunctionBegin;
102236e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10238687889aSStefano Zampini   if (nvtxs) {
1024674ae819SStefano Zampini     PetscValidIntPointer(xadj,3);
10251633d1f0SStefano Zampini     if (xadj[nvtxs]) PetscValidIntPointer(adjncy,4);
10268687889aSStefano Zampini   }
10271a83f524SJed Brown   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
1028575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
1029575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
1030575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
1031575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
1032575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
1033da1bb401SStefano Zampini   }
103436e030ebSStefano Zampini   PetscFunctionReturn(0);
103536e030ebSStefano Zampini }
10361e6b0712SBarry Smith 
103763602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
103863602bcaSStefano Zampini {
103963602bcaSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
104063602bcaSStefano Zampini   PetscInt       i;
104156282151SStefano Zampini   PetscBool      isequal = PETSC_FALSE;
104263602bcaSStefano Zampini   PetscErrorCode ierr;
104363602bcaSStefano Zampini 
104463602bcaSStefano Zampini   PetscFunctionBegin;
104556282151SStefano Zampini   if (pcbddc->n_ISForDofsLocal == n_is) {
104656282151SStefano Zampini     for (i=0;i<n_is;i++) {
104756282151SStefano Zampini       PetscBool isequalt;
104856282151SStefano Zampini       ierr = ISEqual(ISForDofs[i],pcbddc->ISForDofsLocal[i],&isequalt);CHKERRQ(ierr);
104956282151SStefano Zampini       if (!isequalt) break;
105056282151SStefano Zampini     }
105156282151SStefano Zampini     if (i == n_is) isequal = PETSC_TRUE;
105256282151SStefano Zampini   }
105356282151SStefano Zampini   for (i=0;i<n_is;i++) {
105456282151SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
105556282151SStefano Zampini   }
105663602bcaSStefano Zampini   /* Destroy ISes if they were already set */
105763602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
105863602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
105963602bcaSStefano Zampini   }
106063602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
106163602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
106263602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
106363602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
106463602bcaSStefano Zampini   }
106563602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
106663602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
106763602bcaSStefano Zampini   /* allocate space then set */
1068d02579f5SStefano Zampini   if (n_is) {
1069d02579f5SStefano Zampini     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
1070d02579f5SStefano Zampini   }
107163602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
107263602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i] = ISForDofs[i];
107363602bcaSStefano Zampini   }
107463602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = n_is;
107563602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
107656282151SStefano Zampini   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
107763602bcaSStefano Zampini   PetscFunctionReturn(0);
107863602bcaSStefano Zampini }
107963602bcaSStefano Zampini 
108063602bcaSStefano Zampini /*@
108163602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
108263602bcaSStefano Zampini 
108363602bcaSStefano Zampini    Collective
108463602bcaSStefano Zampini 
108563602bcaSStefano Zampini    Input Parameters:
108663602bcaSStefano Zampini +  pc - the preconditioning context
10870f202f7eSStefano Zampini .  n_is - number of index sets defining the fields
10880f202f7eSStefano Zampini -  ISForDofs - array of IS describing the fields in local ordering
108963602bcaSStefano Zampini 
109063602bcaSStefano Zampini    Level: intermediate
109163602bcaSStefano Zampini 
10920f202f7eSStefano Zampini    Notes:
10930f202f7eSStefano Zampini      n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
109463602bcaSStefano Zampini 
109563602bcaSStefano Zampini .seealso: PCBDDC
109663602bcaSStefano Zampini @*/
109763602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
109863602bcaSStefano Zampini {
109963602bcaSStefano Zampini   PetscInt       i;
110063602bcaSStefano Zampini   PetscErrorCode ierr;
110163602bcaSStefano Zampini 
110263602bcaSStefano Zampini   PetscFunctionBegin;
110363602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
110463602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
110563602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
110663602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
110763602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
110863602bcaSStefano Zampini   }
1109e71e7a71SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
111063602bcaSStefano Zampini   PetscFunctionReturn(0);
111163602bcaSStefano Zampini }
111263602bcaSStefano Zampini 
11139c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
11149c0446d6SStefano Zampini {
11159c0446d6SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
11169c0446d6SStefano Zampini   PetscInt       i;
111756282151SStefano Zampini   PetscBool      isequal = PETSC_FALSE;
11189c0446d6SStefano Zampini   PetscErrorCode ierr;
11199c0446d6SStefano Zampini 
11209c0446d6SStefano Zampini   PetscFunctionBegin;
112156282151SStefano Zampini   if (pcbddc->n_ISForDofs == n_is) {
112256282151SStefano Zampini     for (i=0;i<n_is;i++) {
112356282151SStefano Zampini       PetscBool isequalt;
112456282151SStefano Zampini       ierr = ISEqual(ISForDofs[i],pcbddc->ISForDofs[i],&isequalt);CHKERRQ(ierr);
112556282151SStefano Zampini       if (!isequalt) break;
112656282151SStefano Zampini     }
112756282151SStefano Zampini     if (i == n_is) isequal = PETSC_TRUE;
112856282151SStefano Zampini   }
112956282151SStefano Zampini   for (i=0;i<n_is;i++) {
113056282151SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
113156282151SStefano Zampini   }
1132da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
11339c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
11349c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
11359c0446d6SStefano Zampini   }
1136d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
113763602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
113863602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
113963602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
114063602bcaSStefano Zampini   }
114163602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
114263602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
1143da1bb401SStefano Zampini   /* allocate space then set */
1144d02579f5SStefano Zampini   if (n_is) {
1145785e854fSJed Brown     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
1146d02579f5SStefano Zampini   }
11479c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
1148da1bb401SStefano Zampini     pcbddc->ISForDofs[i] = ISForDofs[i];
11499c0446d6SStefano Zampini   }
11509c0446d6SStefano Zampini   pcbddc->n_ISForDofs = n_is;
115163602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
115256282151SStefano Zampini   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
11539c0446d6SStefano Zampini   PetscFunctionReturn(0);
11549c0446d6SStefano Zampini }
11551e6b0712SBarry Smith 
11569c0446d6SStefano Zampini /*@
115763602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
11589c0446d6SStefano Zampini 
115963602bcaSStefano Zampini    Collective
11609c0446d6SStefano Zampini 
11619c0446d6SStefano Zampini    Input Parameters:
11629c0446d6SStefano Zampini +  pc - the preconditioning context
11630f202f7eSStefano Zampini .  n_is - number of index sets defining the fields
11640f202f7eSStefano Zampini -  ISForDofs - array of IS describing the fields in global ordering
11659c0446d6SStefano Zampini 
11669c0446d6SStefano Zampini    Level: intermediate
11679c0446d6SStefano Zampini 
11680f202f7eSStefano Zampini    Notes:
11690f202f7eSStefano Zampini      Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
11709c0446d6SStefano Zampini 
11719c0446d6SStefano Zampini .seealso: PCBDDC
11729c0446d6SStefano Zampini @*/
11739c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
11749c0446d6SStefano Zampini {
11752b510759SStefano Zampini   PetscInt       i;
11769c0446d6SStefano Zampini   PetscErrorCode ierr;
11779c0446d6SStefano Zampini 
11789c0446d6SStefano Zampini   PetscFunctionBegin;
11799c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
118063602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
11812b510759SStefano Zampini   for (i=0;i<n_is;i++) {
118263602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
1183a011d5a7Sstefano_zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
11842b510759SStefano Zampini   }
11859c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
11869c0446d6SStefano Zampini   PetscFunctionReturn(0);
11879c0446d6SStefano Zampini }
1188906d46d4SStefano Zampini 
1189534831adSStefano Zampini /*
1190534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
1191534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
11929c0446d6SStefano Zampini 
1193534831adSStefano Zampini    Input Parameter:
1194534831adSStefano Zampini +  pc - the preconditioner contex
1195534831adSStefano Zampini 
1196534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
1197534831adSStefano Zampini 
1198534831adSStefano Zampini    Notes:
1199534831adSStefano Zampini      The interface routine PCPreSolve() is not usually called directly by
1200534831adSStefano Zampini    the user, but instead is called by KSPSolve().
1201534831adSStefano Zampini */
1202534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1203534831adSStefano Zampini {
1204534831adSStefano Zampini   PetscErrorCode ierr;
1205534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1206534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
12073972b0daSStefano Zampini   Vec            used_vec;
120827b6a85dSStefano Zampini   PetscBool      save_rhs = PETSC_TRUE, benign_correction_computed;
1209534831adSStefano Zampini 
1210534831adSStefano Zampini   PetscFunctionBegin;
12111f4df5f7SStefano Zampini   /* if we are working with CG, one dirichlet solve can be avoided during Krylov iterations */
121285c4d303SStefano Zampini   if (ksp) {
12131f4df5f7SStefano Zampini     PetscBool iscg, isgroppcg, ispipecg, ispipecgrr;
121485c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
121527b6a85dSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPGROPPCG,&isgroppcg);CHKERRQ(ierr);
121627b6a85dSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPPIPECG,&ispipecg);CHKERRQ(ierr);
1217f94e96cbSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPPIPECGRR,&ispipecgrr);CHKERRQ(ierr);
12181f4df5f7SStefano Zampini     if (pcbddc->benign_apply_coarse_only || pcbddc->switch_static || (!iscg && !isgroppcg && !ispipecg && !ispipecgrr)) {
121985c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
122085c4d303SStefano Zampini     }
122185c4d303SStefano Zampini   }
1222fc17d649SStefano Zampini   if (pcbddc->benign_apply_coarse_only || pcbddc->switch_static) {
1223fc17d649SStefano Zampini     ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
1224fc17d649SStefano Zampini   }
12251f4df5f7SStefano Zampini 
122685c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
122762a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
122862a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
122962a6ff1dSStefano Zampini   }
123062a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
123162a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
123262a6ff1dSStefano Zampini   }
12338d00608fSStefano Zampini 
123427b6a85dSStefano Zampini   pcbddc->temp_solution_used = PETSC_FALSE;
12353972b0daSStefano Zampini   if (x) {
12363972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
12373972b0daSStefano Zampini     used_vec = x;
12388d00608fSStefano Zampini   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
12393972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
12403972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
12413972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
124227b6a85dSStefano Zampini     pcbddc->temp_solution_used = PETSC_TRUE;
1243266e20e9SStefano Zampini     ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
1244266e20e9SStefano Zampini     save_rhs = PETSC_FALSE;
1245266e20e9SStefano Zampini     pcbddc->eliminate_dirdofs = PETSC_TRUE;
12463972b0daSStefano Zampini   }
12478efcfb23SStefano Zampini 
12488efcfb23SStefano Zampini   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
12493972b0daSStefano Zampini   if (ksp) {
1250a0cb1b98SStefano Zampini     /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
12518efcfb23SStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
12528efcfb23SStefano Zampini     if (!pcbddc->ksp_guess_nonzero) {
12533972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
12543972b0daSStefano Zampini     }
12553972b0daSStefano Zampini   }
12563308cffdSStefano Zampini 
12578d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
12583972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
125970c64980SStefano Zampini   if (rhs && pcbddc->eliminate_dirdofs) {
12603975b054SStefano Zampini     IS dirIS = NULL;
12613975b054SStefano Zampini 
1262a07ea27aSStefano Zampini     /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
12633975b054SStefano Zampini     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
12643975b054SStefano Zampini     if (dirIS) {
1265906d46d4SStefano Zampini       Mat_IS            *matis = (Mat_IS*)pc->pmat->data;
1266785d1243SStefano Zampini       PetscInt          dirsize,i,*is_indices;
12672b095fd8SStefano Zampini       PetscScalar       *array_x;
12682b095fd8SStefano Zampini       const PetscScalar *array_diagonal;
1269785d1243SStefano Zampini 
12703972b0daSStefano Zampini       ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
12713972b0daSStefano Zampini       ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
1272e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1273e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1274e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1275e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
127682ba6b80SStefano Zampini       ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
12773972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
12782b095fd8SStefano Zampini       ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
12793972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
12802fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
12813972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
12822b095fd8SStefano Zampini       ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
12833972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
1284e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1285e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12868d00608fSStefano Zampini       pcbddc->rhs_change = PETSC_TRUE;
12871b968477SStefano Zampini       ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
12888efcfb23SStefano Zampini     }
1289a07ea27aSStefano Zampini   }
1290b76ba322SStefano Zampini 
12918efcfb23SStefano Zampini   /* remove the computed solution or the initial guess from the rhs */
12928d00608fSStefano Zampini   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
129327b6a85dSStefano Zampini     /* save the original rhs */
129427b6a85dSStefano Zampini     if (save_rhs) {
129527b6a85dSStefano Zampini       ierr = VecSwap(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
129627b6a85dSStefano Zampini       save_rhs = PETSC_FALSE;
12978d00608fSStefano Zampini     }
12988d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
12993972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
130027b6a85dSStefano Zampini     ierr = MatMultAdd(pc->mat,used_vec,pcbddc->original_rhs,rhs);CHKERRQ(ierr);
13013972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
13028efcfb23SStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
130327b6a85dSStefano Zampini     pcbddc->temp_solution_used = PETSC_TRUE;
13047acc28cbSStefano Zampini     if (ksp) {
13057acc28cbSStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr);
13067acc28cbSStefano Zampini     }
13073308cffdSStefano Zampini   }
13088efcfb23SStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1309b76ba322SStefano Zampini 
1310fc17d649SStefano Zampini   /* compute initial vector in benign space if needed
131127b6a85dSStefano Zampini      and remove non-benign solution from the rhs */
131227b6a85dSStefano Zampini   benign_correction_computed = PETSC_FALSE;
131327b6a85dSStefano Zampini   if (rhs && pcbddc->benign_compute_correction && pcbddc->benign_have_null) {
13141f4df5f7SStefano Zampini     /* compute u^*_h using ideas similar to those in Xuemin Tu's PhD thesis (see Section 4.8.1)
13151f4df5f7SStefano Zampini        Recursively apply BDDC in the multilevel case */
13160369aaf7SStefano Zampini     if (!pcbddc->benign_vec) {
13170369aaf7SStefano Zampini       ierr = VecDuplicate(rhs,&pcbddc->benign_vec);CHKERRQ(ierr);
13180369aaf7SStefano Zampini     }
13194fee134fSStefano Zampini     pcbddc->benign_apply_coarse_only = PETSC_TRUE;
132027b6a85dSStefano Zampini     if (!pcbddc->benign_skip_correction) {
13211dd7afcfSStefano Zampini       ierr = PCApply_BDDC(pc,rhs,pcbddc->benign_vec);CHKERRQ(ierr);
13223bca92a6SStefano Zampini       benign_correction_computed = PETSC_TRUE;
13231f4df5f7SStefano Zampini       if (pcbddc->temp_solution_used) {
13241f4df5f7SStefano Zampini         ierr = VecAXPY(pcbddc->temp_solution,1.0,pcbddc->benign_vec);CHKERRQ(ierr);
13251f4df5f7SStefano Zampini       }
13261f4df5f7SStefano Zampini       ierr = VecScale(pcbddc->benign_vec,-1.0);CHKERRQ(ierr);
132727b6a85dSStefano Zampini       /* store the original rhs if not done earlier */
132827b6a85dSStefano Zampini       if (save_rhs) {
132927b6a85dSStefano Zampini         ierr = VecSwap(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
133092e3dcfbSStefano Zampini       }
133127b6a85dSStefano Zampini       if (pcbddc->rhs_change) {
13320369aaf7SStefano Zampini         ierr = MatMultAdd(pc->mat,pcbddc->benign_vec,rhs,rhs);CHKERRQ(ierr);
133327b6a85dSStefano Zampini       } else {
133427b6a85dSStefano Zampini         ierr = MatMultAdd(pc->mat,pcbddc->benign_vec,pcbddc->original_rhs,rhs);CHKERRQ(ierr);
133527b6a85dSStefano Zampini       }
13360369aaf7SStefano Zampini       pcbddc->rhs_change = PETSC_TRUE;
133727b6a85dSStefano Zampini     }
133827b6a85dSStefano Zampini     pcbddc->benign_apply_coarse_only = PETSC_FALSE;
13390369aaf7SStefano Zampini   }
13402d4c4fecSStefano Zampini 
13412d4c4fecSStefano Zampini   /* dbg output */
1342a198735bSStefano Zampini   if (pcbddc->dbg_flag && benign_correction_computed) {
13431f4df5f7SStefano Zampini     Vec v;
13441f4df5f7SStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&v);CHKERRQ(ierr);
13451f4df5f7SStefano Zampini     ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,rhs,v);CHKERRQ(ierr);
13461f4df5f7SStefano Zampini     ierr = PCBDDCBenignGetOrSetP0(pc,v,PETSC_TRUE);CHKERRQ(ierr);
1347e0fe2d75SToby Isaac     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"LEVEL %D: is the correction benign?\n",pcbddc->current_level);CHKERRQ(ierr);
13482d4c4fecSStefano Zampini     ierr = PetscScalarView(pcbddc->benign_n,pcbddc->benign_p0,PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)));CHKERRQ(ierr);
13491f4df5f7SStefano Zampini     ierr = VecDestroy(&v);CHKERRQ(ierr);
13501f4df5f7SStefano Zampini   }
13510369aaf7SStefano Zampini 
13520369aaf7SStefano Zampini   /* set initial guess if using PCG */
13538ae0ca82SStefano Zampini   pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
13540369aaf7SStefano Zampini   if (x && pcbddc->use_exact_dirichlet_trick) {
13550369aaf7SStefano Zampini     ierr = VecSet(x,0.0);CHKERRQ(ierr);
13561dd7afcfSStefano Zampini     if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) {
135727b6a85dSStefano Zampini       if (benign_correction_computed) { /* we have already saved the changed rhs */
13581dd7afcfSStefano Zampini         ierr = VecLockPop(pcis->vec1_global);CHKERRQ(ierr);
13591dd7afcfSStefano Zampini       } else {
13601dd7afcfSStefano Zampini         ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,rhs,pcis->vec1_global);CHKERRQ(ierr);
13611dd7afcfSStefano Zampini       }
13621dd7afcfSStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_global,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13631dd7afcfSStefano Zampini       ierr = VecScatterEnd(pcis->global_to_D,pcis->vec1_global,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13641dd7afcfSStefano Zampini     } else {
13650369aaf7SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13660369aaf7SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13671dd7afcfSStefano Zampini     }
13680369aaf7SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
13691dd7afcfSStefano Zampini     if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) {
13701dd7afcfSStefano Zampini       ierr = VecSet(pcis->vec1_global,0.);CHKERRQ(ierr);
13711dd7afcfSStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13721dd7afcfSStefano Zampini       ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13731dd7afcfSStefano Zampini       ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_global,x);CHKERRQ(ierr);
13741dd7afcfSStefano Zampini     } else {
13750369aaf7SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13760369aaf7SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13771dd7afcfSStefano Zampini     }
13780369aaf7SStefano Zampini     if (ksp) {
13790369aaf7SStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
13800369aaf7SStefano Zampini     }
13818ae0ca82SStefano Zampini     pcbddc->exact_dirichlet_trick_app = PETSC_TRUE;
1382266e20e9SStefano Zampini   } else if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior && benign_correction_computed && pcbddc->use_exact_dirichlet_trick) {
1383266e20e9SStefano Zampini     ierr = VecLockPop(pcis->vec1_global);CHKERRQ(ierr);
13840369aaf7SStefano Zampini   }
1385534831adSStefano Zampini   PetscFunctionReturn(0);
1386534831adSStefano Zampini }
1387906d46d4SStefano Zampini 
1388534831adSStefano Zampini /*
1389534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1390534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
1391534831adSStefano Zampini 
1392534831adSStefano Zampini    Input Parameter:
1393534831adSStefano Zampini +  pc - the preconditioner contex
1394534831adSStefano Zampini 
1395534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
1396534831adSStefano Zampini 
1397534831adSStefano Zampini    Notes:
1398534831adSStefano Zampini      The interface routine PCPostSolve() is not usually called directly by
1399534831adSStefano Zampini      the user, but instead is called by KSPSolve().
1400534831adSStefano Zampini */
1401534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1402534831adSStefano Zampini {
1403534831adSStefano Zampini   PetscErrorCode ierr;
1404534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1405534831adSStefano Zampini 
1406534831adSStefano Zampini   PetscFunctionBegin;
14073972b0daSStefano Zampini   /* add solution removed in presolve */
14086bcfc461SStefano Zampini   if (x && pcbddc->rhs_change) {
140927b6a85dSStefano Zampini     if (pcbddc->temp_solution_used) {
14103425bc38SStefano Zampini       ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
1411af140850Sstefano_zampini     } else if (pcbddc->benign_compute_correction && pcbddc->benign_vec) {
141227b6a85dSStefano Zampini       ierr = VecAXPY(x,-1.0,pcbddc->benign_vec);CHKERRQ(ierr);
14133425bc38SStefano Zampini     }
1414af140850Sstefano_zampini     /* restore to original state (not for FETI-DP) */
1415af140850Sstefano_zampini     if (ksp) pcbddc->temp_solution_used = PETSC_FALSE;
141627b6a85dSStefano Zampini   }
141727b6a85dSStefano Zampini 
1418266e20e9SStefano Zampini   /* restore rhs to its original state (not needed for FETI-DP) */
14198d00608fSStefano Zampini   if (rhs && pcbddc->rhs_change) {
142027b6a85dSStefano Zampini     ierr = VecSwap(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
14218d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_FALSE;
1422af140850Sstefano_zampini   }
14238efcfb23SStefano Zampini   /* restore ksp guess state */
14248efcfb23SStefano Zampini   if (ksp) {
14258efcfb23SStefano Zampini     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
14268ae0ca82SStefano Zampini     /* reset flag for exact dirichlet trick */
14278ae0ca82SStefano Zampini     pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
1428af140850Sstefano_zampini   }
1429534831adSStefano Zampini   PetscFunctionReturn(0);
1430534831adSStefano Zampini }
1431af140850Sstefano_zampini 
14320c7d97c5SJed Brown /*
14330c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
14340c7d97c5SJed Brown                   by setting data structures and options.
14350c7d97c5SJed Brown 
14360c7d97c5SJed Brown    Input Parameter:
143753cdbc3dSStefano Zampini +  pc - the preconditioner context
14380c7d97c5SJed Brown 
14390c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
14400c7d97c5SJed Brown 
14410c7d97c5SJed Brown    Notes:
14420c7d97c5SJed Brown      The interface routine PCSetUp() is not usually called directly by
14430c7d97c5SJed Brown      the user, but instead is called by PCApply() if necessary.
14440c7d97c5SJed Brown */
144553cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
14460c7d97c5SJed Brown {
14470c7d97c5SJed Brown   PC_BDDC*        pcbddc = (PC_BDDC*)pc->data;
1448c703fcc7SStefano Zampini   PCBDDCSubSchurs sub_schurs;
14495e8657edSStefano Zampini   Mat_IS*         matis;
145008122e43SStefano Zampini   MatNullSpace    nearnullspace;
145135509ce9Sstefano_zampini   Mat             lA;
145235509ce9Sstefano_zampini   IS              lP,zerodiag = NULL;
145391e8d312SStefano Zampini   PetscInt        nrows,ncols;
1454c703fcc7SStefano Zampini   PetscBool       computesubschurs;
14558de1fae6SStefano Zampini   PetscBool       computeconstraintsmatrix;
14565e8657edSStefano Zampini   PetscBool       new_nearnullspace_provided,ismatis;
1457c703fcc7SStefano Zampini   PetscErrorCode  ierr;
14580c7d97c5SJed Brown 
14590c7d97c5SJed Brown   PetscFunctionBegin;
14605e8657edSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr);
14616c4ed002SBarry Smith   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
146291e8d312SStefano Zampini   ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr);
14636c4ed002SBarry Smith   if (nrows != ncols) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
14645e8657edSStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
1465f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
14663b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
146771582508SStefano Zampini      Also, BDDC builds its own KSP for the Dirichlet problem */
1468a485753aSstefano_zampini   if (!pc->setupcalled || pc->flag == DIFFERENT_NONZERO_PATTERN) pcbddc->recompute_topography = PETSC_TRUE;
1469c83e1ba7SStefano Zampini   if (pcbddc->recompute_topography) {
1470c83e1ba7SStefano Zampini     pcbddc->graphanalyzed    = PETSC_FALSE;
1471c83e1ba7SStefano Zampini     computeconstraintsmatrix = PETSC_TRUE;
1472c83e1ba7SStefano Zampini   } else {
14738de1fae6SStefano Zampini     computeconstraintsmatrix = PETSC_FALSE;
1474c83e1ba7SStefano Zampini   }
1475b087196eSStefano Zampini 
1476b087196eSStefano Zampini   /* check parameters' compatibility */
1477b7ab4a40SStefano Zampini   if (!pcbddc->use_deluxe_scaling) pcbddc->deluxe_zerorows = PETSC_FALSE;
1478bd2a564bSStefano Zampini   pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold[0] != 0.0 || pcbddc->adaptive_threshold[1] != 0.0);
1479bf3a8328SStefano Zampini   pcbddc->adaptive_userdefined = (PetscBool)(pcbddc->adaptive_selection && pcbddc->adaptive_userdefined);
1480862806e4SStefano Zampini   if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;
1481862806e4SStefano Zampini 
14825a95e1ceSStefano Zampini   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
148316909a7fSStefano Zampini   if (pcbddc->switch_static) {
148416909a7fSStefano Zampini     PetscBool ismatis;
148516909a7fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&ismatis);CHKERRQ(ierr);
148616909a7fSStefano Zampini     if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"When the static switch is one, the iteration matrix should be of type MATIS");
148716909a7fSStefano Zampini   }
148816909a7fSStefano Zampini 
148971582508SStefano Zampini   /* activate all connected components if the netflux has been requested */
1490bb05f991SStefano Zampini   if (pcbddc->compute_nonetflux) {
1491bb05f991SStefano Zampini     pcbddc->use_vertices = PETSC_TRUE;
1492bb05f991SStefano Zampini     pcbddc->use_edges = PETSC_TRUE;
1493bb05f991SStefano Zampini     pcbddc->use_faces = PETSC_TRUE;
1494bb05f991SStefano Zampini   }
1495bb05f991SStefano Zampini 
1496f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
149770cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
149870cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
149958a03d70SStefano Zampini       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
15001575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(pcbddc->dbg_viewer);CHKERRQ(ierr);
1501f4ddd8eeSStefano Zampini     }
150258a03d70SStefano Zampini     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1503f4ddd8eeSStefano Zampini   }
1504f4ddd8eeSStefano Zampini 
1505c703fcc7SStefano Zampini   /* process topology information */
150671582508SStefano Zampini   if (pcbddc->recompute_topography) {
150771582508SStefano Zampini     ierr = PCBDDCComputeLocalTopologyInfo(pc);CHKERRQ(ierr);
1508c703fcc7SStefano Zampini     if (pcbddc->discretegradient) {
1509a13144ffSStefano Zampini       ierr = PCBDDCNedelecSupport(pc);CHKERRQ(ierr);
1510a13144ffSStefano Zampini     }
1511c703fcc7SStefano Zampini   }
1512a13144ffSStefano Zampini 
1513c703fcc7SStefano Zampini   /* change basis if requested by the user */
15145e8657edSStefano Zampini   if (pcbddc->user_ChangeOfBasisMatrix) {
15155e8657edSStefano Zampini     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
15165e8657edSStefano Zampini     pcbddc->use_change_of_basis = PETSC_FALSE;
15175e8657edSStefano Zampini     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
15185e8657edSStefano Zampini   } else {
1519b96c3477SStefano Zampini     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
15205e8657edSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
15215e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
1522d16cbb6bSStefano Zampini   }
1523d16cbb6bSStefano Zampini 
15244f1b2e48SStefano Zampini   /*
1525c703fcc7SStefano Zampini      Compute change of basis on local pressures (aka zerodiag dofs) with the benign trick
15264f1b2e48SStefano Zampini      This should come earlier then PCISSetUp for extracting the correct subdomain matrices
15274f1b2e48SStefano Zampini   */
15281dd7afcfSStefano Zampini   ierr = PCBDDCBenignShellMat(pc,PETSC_TRUE);CHKERRQ(ierr);
1529d16cbb6bSStefano Zampini   if (pcbddc->benign_saddle_point) {
15309f47a83aSStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
15319f47a83aSStefano Zampini 
153205b28244SStefano Zampini     if (pcbddc->user_ChangeOfBasisMatrix || pcbddc->use_change_of_basis || !computesubschurs) pcbddc->benign_change_explicit = PETSC_TRUE;
1533669cc0f4SStefano Zampini     /* detect local saddle point and change the basis in pcbddc->local_mat (TODO: reuse case) */
1534339f8db1SStefano Zampini     ierr = PCBDDCBenignDetectSaddlePoint(pc,&zerodiag);CHKERRQ(ierr);
1535a3df083aSStefano Zampini     /* pop B0 mat from local mat */
1536c263805aSStefano Zampini     ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
15371dd7afcfSStefano Zampini     /* give pcis a hint to not reuse submatrices during PCISCreate */
15381dd7afcfSStefano Zampini     if (pc->flag == SAME_NONZERO_PATTERN && pcis->reusesubmatrices == PETSC_TRUE) {
15391dd7afcfSStefano Zampini       if (pcbddc->benign_n && (pcbddc->benign_change_explicit || pcbddc->dbg_flag)) {
15401dd7afcfSStefano Zampini         pcis->reusesubmatrices = PETSC_FALSE;
15411dd7afcfSStefano Zampini       } else {
1542a3df083aSStefano Zampini         pcis->reusesubmatrices = PETSC_TRUE;
15431dd7afcfSStefano Zampini       }
1544a3df083aSStefano Zampini     } else {
15459f47a83aSStefano Zampini       pcis->reusesubmatrices = PETSC_FALSE;
1546674ae819SStefano Zampini     }
1547a3df083aSStefano Zampini   }
154827b6a85dSStefano Zampini 
1549c703fcc7SStefano Zampini   /* propagate relevant information -> TODO remove*/
15504f1b2e48SStefano Zampini #if !defined(PETSC_USE_COMPLEX) /* workaround for reals */
15513301b35fSStefano Zampini   if (matis->A->symmetric_set) {
15523301b35fSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr);
15533301b35fSStefano Zampini   }
1554e496cd5dSStefano Zampini #endif
155506a4e24aSStefano Zampini   if (matis->A->symmetric_set) {
155606a4e24aSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);CHKERRQ(ierr);
155706a4e24aSStefano Zampini   }
155806a4e24aSStefano Zampini   if (matis->A->spd_set) {
155906a4e24aSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);CHKERRQ(ierr);
156006a4e24aSStefano Zampini   }
1561e496cd5dSStefano Zampini 
15625e8657edSStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
15635e8657edSStefano Zampini   {
15645e8657edSStefano Zampini     Mat temp_mat;
15655e8657edSStefano Zampini 
15665e8657edSStefano Zampini     temp_mat = matis->A;
15675e8657edSStefano Zampini     matis->A = pcbddc->local_mat;
15685e8657edSStefano Zampini     ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr);
15695e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
15705e8657edSStefano Zampini     matis->A = temp_mat;
15715e8657edSStefano Zampini   }
1572684f6988SStefano Zampini 
157381d14e9dSStefano Zampini   /* Analyze interface */
157464ac59b8SStefano Zampini   if (!pcbddc->graphanalyzed) {
1575674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
15768de1fae6SStefano Zampini     computeconstraintsmatrix = PETSC_TRUE;
1577345ecf6cSStefano Zampini     if (pcbddc->adaptive_selection && !pcbddc->use_deluxe_scaling && !pcbddc->mat_graph->twodim) {
15784247aa23Sstefano_zampini       SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute the adaptive primal space for a problem with 3D edges without deluxe scaling");
1579345ecf6cSStefano Zampini     }
1580a198735bSStefano Zampini     if (pcbddc->compute_nonetflux) {
1581669cc0f4SStefano Zampini       MatNullSpace nnfnnsp;
1582669cc0f4SStefano Zampini 
158321ef3d20SStefano Zampini       if (!pcbddc->divudotp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Missing divudotp operator");
15848ae0ca82SStefano Zampini       ierr = PCBDDCComputeNoNetFlux(pc->pmat,pcbddc->divudotp,pcbddc->divudotp_trans,pcbddc->divudotp_vl2l,pcbddc->mat_graph,&nnfnnsp);CHKERRQ(ierr);
158571582508SStefano Zampini       /* TODO what if a nearnullspace is already attached? */
1586669cc0f4SStefano Zampini       ierr = MatSetNearNullSpace(pc->pmat,nnfnnsp);CHKERRQ(ierr);
1587669cc0f4SStefano Zampini       ierr = MatNullSpaceDestroy(&nnfnnsp);CHKERRQ(ierr);
1588669cc0f4SStefano Zampini     }
1589674ae819SStefano Zampini   }
1590fb8d54d4SStefano Zampini 
15915408967cSStefano Zampini   /* check existence of a divergence free extension, i.e.
15925408967cSStefano Zampini      b(v_I,p_0) = 0 for all v_I (raise error if not).
15935408967cSStefano Zampini      Also, check that PCBDDCBenignGetOrSetP0 works */
1594ff1f7e73Sstefano_zampini   if (pcbddc->benign_saddle_point && pcbddc->dbg_flag > 1) {
15955408967cSStefano Zampini     ierr = PCBDDCBenignCheck(pc,zerodiag);CHKERRQ(ierr);
159609f581a4SStefano Zampini   }
15974f1b2e48SStefano Zampini   ierr = ISDestroy(&zerodiag);CHKERRQ(ierr);
159806f24817SStefano Zampini 
1599b96c3477SStefano Zampini   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1600c703fcc7SStefano Zampini   if (computesubschurs && pcbddc->recompute_topography) {
160108122e43SStefano Zampini     ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr);
1602b1b3d7a2SStefano Zampini   }
16039d54b7f4SStefano Zampini   /* SetUp Scaling operator (scaling matrices could be needed in SubSchursSetUp)*/
16049d54b7f4SStefano Zampini   if (!pcbddc->use_deluxe_scaling) {
16059d54b7f4SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
16069d54b7f4SStefano Zampini   }
1607c703fcc7SStefano Zampini 
1608c703fcc7SStefano Zampini   /* finish setup solvers and do adaptive selection of constraints */
1609b334f244SStefano Zampini   sub_schurs = pcbddc->sub_schurs;
1610b334f244SStefano Zampini   if (sub_schurs && sub_schurs->schur_explicit) {
16112070dbb6SStefano Zampini     if (computesubschurs) {
161208122e43SStefano Zampini       ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
16132070dbb6SStefano Zampini     }
1614d5574798SStefano Zampini     ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1615d5574798SStefano Zampini   } else {
1616d5574798SStefano Zampini     ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
16172070dbb6SStefano Zampini     if (computesubschurs) {
1618d5574798SStefano Zampini       ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1619d5574798SStefano Zampini     }
16202070dbb6SStefano Zampini   }
162108122e43SStefano Zampini   if (pcbddc->adaptive_selection) {
162208122e43SStefano Zampini     ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr);
16238de1fae6SStefano Zampini     computeconstraintsmatrix = PETSC_TRUE;
1624b7eb3628SStefano Zampini   }
1625684f6988SStefano Zampini 
1626f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1627fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1628f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1629f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1630f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1631f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1632f4ddd8eeSStefano Zampini     } else {
1633f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1634f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1635f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1636165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1637f4ddd8eeSStefano Zampini         PetscInt         i;
1638165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1639165b64e2SStefano Zampini         PetscObjectState state;
1640165b64e2SStefano Zampini         PetscInt         nnsp_size;
1641165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1642f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1643f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1644165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1645f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1646f4ddd8eeSStefano Zampini             break;
1647f4ddd8eeSStefano Zampini           }
1648f4ddd8eeSStefano Zampini         }
1649f4ddd8eeSStefano Zampini       }
1650f4ddd8eeSStefano Zampini     }
1651f4ddd8eeSStefano Zampini   } else {
1652f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1653f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1654f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1655f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1656f4ddd8eeSStefano Zampini     }
1657f4ddd8eeSStefano Zampini   }
1658f4ddd8eeSStefano Zampini 
1659f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1660727cdba6SStefano Zampini   /* reset primal space flags */
1661f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1662727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
16638de1fae6SStefano Zampini   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1664727cdba6SStefano Zampini     /* It also sets the primal space flags */
1665674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
16669543d0ffSStefano Zampini   }
1667e7b262bdSStefano Zampini   /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1668f4ddd8eeSStefano Zampini   ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
16695e8657edSStefano Zampini 
16705e8657edSStefano Zampini   if (pcbddc->use_change_of_basis) {
16715e8657edSStefano Zampini     PC_IS *pcis = (PC_IS*)(pc->data);
16725e8657edSStefano Zampini 
16735e8657edSStefano Zampini     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
16744f1b2e48SStefano Zampini     if (pcbddc->benign_change) {
16751dd7afcfSStefano Zampini       ierr = MatDestroy(&pcbddc->benign_B0);CHKERRQ(ierr);
1676c263805aSStefano Zampini       /* pop B0 from pcbddc->local_mat */
1677c263805aSStefano Zampini       ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
1678c263805aSStefano Zampini     }
16795e8657edSStefano Zampini     /* get submatrices */
16805e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
16815e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
16825e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
16837dae84e0SHong Zhang     ierr = MatCreateSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
16847dae84e0SHong Zhang     ierr = MatCreateSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
16857dae84e0SHong Zhang     ierr = MatCreateSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
16863975b054SStefano Zampini     /* set flag in pcis to not reuse submatrices during PCISCreate */
16873975b054SStefano Zampini     pcis->reusesubmatrices = PETSC_FALSE;
16889c6a02ceSStefano Zampini   } else if (!pcbddc->user_ChangeOfBasisMatrix && !pcbddc->benign_change) {
1689b96c3477SStefano Zampini     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
16905e8657edSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
16915e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
16925e8657edSStefano Zampini   }
169335509ce9Sstefano_zampini 
169435509ce9Sstefano_zampini   /* interface pressure block row for B_C */
169535509ce9Sstefano_zampini   ierr = PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_lP" ,(PetscObject*)&lP);CHKERRQ(ierr);
169635509ce9Sstefano_zampini   ierr = PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_lA" ,(PetscObject*)&lA);CHKERRQ(ierr);
169735509ce9Sstefano_zampini   if (lA && lP) {
169835509ce9Sstefano_zampini     PC_IS*    pcis = (PC_IS*)pc->data;
169935509ce9Sstefano_zampini     Mat       B_BI,B_BB,Bt_BI,Bt_BB;
170035509ce9Sstefano_zampini     PetscBool issym;
170135509ce9Sstefano_zampini     ierr = MatIsSymmetric(lA,PETSC_SMALL,&issym);CHKERRQ(ierr);
17026cc1294bSstefano_zampini     if (issym) {
17037dae84e0SHong Zhang       ierr = MatCreateSubMatrix(lA,lP,pcis->is_I_local,MAT_INITIAL_MATRIX,&B_BI);CHKERRQ(ierr);
17047dae84e0SHong Zhang       ierr = MatCreateSubMatrix(lA,lP,pcis->is_B_local,MAT_INITIAL_MATRIX,&B_BB);CHKERRQ(ierr);
170535509ce9Sstefano_zampini       ierr = MatCreateTranspose(B_BI,&Bt_BI);CHKERRQ(ierr);
170635509ce9Sstefano_zampini       ierr = MatCreateTranspose(B_BB,&Bt_BB);CHKERRQ(ierr);
170735509ce9Sstefano_zampini     } else {
17087dae84e0SHong Zhang       ierr = MatCreateSubMatrix(lA,lP,pcis->is_I_local,MAT_INITIAL_MATRIX,&B_BI);CHKERRQ(ierr);
17097dae84e0SHong Zhang       ierr = MatCreateSubMatrix(lA,lP,pcis->is_B_local,MAT_INITIAL_MATRIX,&B_BB);CHKERRQ(ierr);
17107dae84e0SHong Zhang       ierr = MatCreateSubMatrix(lA,pcis->is_I_local,lP,MAT_INITIAL_MATRIX,&Bt_BI);CHKERRQ(ierr);
17117dae84e0SHong Zhang       ierr = MatCreateSubMatrix(lA,pcis->is_B_local,lP,MAT_INITIAL_MATRIX,&Bt_BB);CHKERRQ(ierr);
171235509ce9Sstefano_zampini     }
171335509ce9Sstefano_zampini     ierr = PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_B_BI",(PetscObject)B_BI);CHKERRQ(ierr);
171435509ce9Sstefano_zampini     ierr = PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_B_BB",(PetscObject)B_BB);CHKERRQ(ierr);
171535509ce9Sstefano_zampini     ierr = PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_Bt_BI",(PetscObject)Bt_BI);CHKERRQ(ierr);
171635509ce9Sstefano_zampini     ierr = PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_Bt_BB",(PetscObject)Bt_BB);CHKERRQ(ierr);
171735509ce9Sstefano_zampini     ierr = MatDestroy(&B_BI);CHKERRQ(ierr);
171835509ce9Sstefano_zampini     ierr = MatDestroy(&B_BB);CHKERRQ(ierr);
171935509ce9Sstefano_zampini     ierr = MatDestroy(&Bt_BI);CHKERRQ(ierr);
172035509ce9Sstefano_zampini     ierr = MatDestroy(&Bt_BB);CHKERRQ(ierr);
172135509ce9Sstefano_zampini   }
172235509ce9Sstefano_zampini 
1723b96c3477SStefano Zampini   /* SetUp coarse and local Neumann solvers */
172499cc7994SStefano Zampini   ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1725b96c3477SStefano Zampini   /* SetUp Scaling operator */
17269d54b7f4SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
1727674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
17280c7d97c5SJed Brown   }
1729c703fcc7SStefano Zampini 
17301dd7afcfSStefano Zampini   /* mark topography as done */
173156282151SStefano Zampini   pcbddc->recompute_topography = PETSC_FALSE;
17320369aaf7SStefano Zampini 
17331dd7afcfSStefano Zampini   /* wrap pcis->A_IB and pcis->A_BI if we did not change explicitly the variables on the pressures */
17341dd7afcfSStefano Zampini   ierr = PCBDDCBenignShellMat(pc,PETSC_FALSE);CHKERRQ(ierr);
17351dd7afcfSStefano Zampini 
173658a03d70SStefano Zampini   if (pcbddc->dbg_flag) {
173758a03d70SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
17382b510759SStefano Zampini   }
17390c7d97c5SJed Brown   PetscFunctionReturn(0);
17400c7d97c5SJed Brown }
17410c7d97c5SJed Brown 
17420c7d97c5SJed Brown /*
174350efa1b5SStefano Zampini    PCApply_BDDC - Applies the BDDC operator to a vector.
17440c7d97c5SJed Brown 
17450c7d97c5SJed Brown    Input Parameters:
17460f202f7eSStefano Zampini +  pc - the preconditioner context
17470f202f7eSStefano Zampini -  r - input vector (global)
17480c7d97c5SJed Brown 
17490c7d97c5SJed Brown    Output Parameter:
17500c7d97c5SJed Brown .  z - output vector (global)
17510c7d97c5SJed Brown 
17520c7d97c5SJed Brown    Application Interface Routine: PCApply()
17530c7d97c5SJed Brown  */
175453cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
17550c7d97c5SJed Brown {
17560c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
17570c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1758b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
17590c7d97c5SJed Brown   PetscErrorCode    ierr;
17603b03a366Sstefano_zampini   const PetscScalar one = 1.0;
17613b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
17622617d88aSStefano Zampini   const PetscScalar zero = 0.0;
17630c7d97c5SJed Brown 
17640c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
17650c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
1766b097fa66SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
17670c7d97c5SJed Brown 
17680c7d97c5SJed Brown   PetscFunctionBegin;
1769f3d41395Sstefano_zampini   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
17701dd7afcfSStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
17711dd7afcfSStefano Zampini     Vec swap;
177227b6a85dSStefano Zampini 
177327b6a85dSStefano Zampini     ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,r,pcbddc->work_change);CHKERRQ(ierr);
17741dd7afcfSStefano Zampini     swap = pcbddc->work_change;
17751dd7afcfSStefano Zampini     pcbddc->work_change = r;
17761dd7afcfSStefano Zampini     r = swap;
17771dd7afcfSStefano Zampini     /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */
17789cc2a9b1Sstefano_zampini     if (pcbddc->benign_apply_coarse_only && pcbddc->use_exact_dirichlet_trick && pcbddc->change_interior) {
17791dd7afcfSStefano Zampini       ierr = VecCopy(r,pcis->vec1_global);CHKERRQ(ierr);
17801dd7afcfSStefano Zampini       ierr = VecLockPush(pcis->vec1_global);CHKERRQ(ierr);
17811dd7afcfSStefano Zampini     }
17821dd7afcfSStefano Zampini   }
178327b6a85dSStefano Zampini   if (pcbddc->benign_have_null) { /* get p0 from r */
1784015636ebSStefano Zampini     ierr = PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);CHKERRQ(ierr);
1785efc2fbd9SStefano Zampini   }
17868ae0ca82SStefano Zampini   if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
1787b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
17880c7d97c5SJed Brown     /* First Dirichlet solve */
17890c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17900c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17910c7d97c5SJed Brown     /*
17920c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1793b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1794674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
17950c7d97c5SJed Brown     */
1796b097fa66SStefano Zampini     if (n_D) {
1797b097fa66SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
17980c7d97c5SJed Brown       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
179916909a7fSStefano Zampini       if (pcbddc->switch_static) {
180016909a7fSStefano Zampini         Mat_IS *matis = (Mat_IS*)(pc->mat->data);
180116909a7fSStefano Zampini 
180216909a7fSStefano Zampini         ierr = VecSet(pcis->vec1_N,0.);CHKERRQ(ierr);
180316909a7fSStefano Zampini         ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
180416909a7fSStefano Zampini         ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
180516909a7fSStefano Zampini         if (!pcbddc->switch_static_change) {
180616909a7fSStefano Zampini           ierr = MatMult(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
180716909a7fSStefano Zampini         } else {
180816909a7fSStefano Zampini           ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
180916909a7fSStefano Zampini           ierr = MatMult(matis->A,pcis->vec2_N,pcis->vec1_N);CHKERRQ(ierr);
181016909a7fSStefano Zampini           ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
181116909a7fSStefano Zampini         }
181216909a7fSStefano Zampini         ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
181316909a7fSStefano Zampini         ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
181416909a7fSStefano Zampini         ierr = VecScatterBegin(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
181516909a7fSStefano Zampini         ierr = VecScatterEnd(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
181616909a7fSStefano Zampini       } else {
1817b097fa66SStefano Zampini         ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
181816909a7fSStefano Zampini       }
1819b097fa66SStefano Zampini     } else {
1820b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1821b097fa66SStefano Zampini     }
18220c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
18230c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1824674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1825b76ba322SStefano Zampini   } else {
18264fee134fSStefano Zampini     if (!pcbddc->benign_apply_coarse_only) {
1827674ae819SStefano Zampini       ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1828b76ba322SStefano Zampini     }
18294fee134fSStefano Zampini   }
1830b76ba322SStefano Zampini 
18312617d88aSStefano Zampini   /* Apply interface preconditioner
18322617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1833dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
18342617d88aSStefano Zampini 
1835674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1836674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
18370c7d97c5SJed Brown 
18383b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
18390c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
18400c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1841b097fa66SStefano Zampini   if (n_B) {
184216909a7fSStefano Zampini     if (pcbddc->switch_static) {
184316909a7fSStefano Zampini       Mat_IS *matis = (Mat_IS*)(pc->mat->data);
184416909a7fSStefano Zampini 
184516909a7fSStefano Zampini       ierr = VecScatterBegin(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
184616909a7fSStefano Zampini       ierr = VecScatterEnd(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
184716909a7fSStefano Zampini       ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
184816909a7fSStefano Zampini       ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
184916909a7fSStefano Zampini       if (!pcbddc->switch_static_change) {
185016909a7fSStefano Zampini         ierr = MatMult(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
185116909a7fSStefano Zampini       } else {
185216909a7fSStefano Zampini         ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
185316909a7fSStefano Zampini         ierr = MatMult(matis->A,pcis->vec2_N,pcis->vec1_N);CHKERRQ(ierr);
185416909a7fSStefano Zampini         ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
185516909a7fSStefano Zampini       }
185616909a7fSStefano Zampini       ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
185716909a7fSStefano Zampini       ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
185816909a7fSStefano Zampini     } else {
18590c7d97c5SJed Brown       ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
186016909a7fSStefano Zampini     }
186116909a7fSStefano Zampini   } else if (pcbddc->switch_static) { /* n_B is zero */
186216909a7fSStefano Zampini     Mat_IS *matis = (Mat_IS*)(pc->mat->data);
186316909a7fSStefano Zampini 
186416909a7fSStefano Zampini     if (!pcbddc->switch_static_change) {
186516909a7fSStefano Zampini       ierr = MatMult(matis->A,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
186616909a7fSStefano Zampini     } else {
186716909a7fSStefano Zampini       ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_D,pcis->vec1_N);CHKERRQ(ierr);
186816909a7fSStefano Zampini       ierr = MatMult(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
186916909a7fSStefano Zampini       ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec2_N,pcis->vec3_D);CHKERRQ(ierr);
187016909a7fSStefano Zampini     }
1871b097fa66SStefano Zampini   }
1872df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1873efc2fbd9SStefano Zampini 
18748ae0ca82SStefano Zampini   if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
1875b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1876b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1877b097fa66SStefano Zampini     } else {
1878b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1879b097fa66SStefano Zampini     }
18800c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
18810c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1882b097fa66SStefano Zampini   } else {
1883b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1884b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1885b097fa66SStefano Zampini     } else {
1886b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1887b097fa66SStefano Zampini     }
1888b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1889b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1890b097fa66SStefano Zampini   }
189127b6a85dSStefano Zampini   if (pcbddc->benign_have_null) { /* set p0 (computed in PCBDDCApplyInterface) */
18921dd7afcfSStefano Zampini     if (pcbddc->benign_apply_coarse_only) {
18931dd7afcfSStefano Zampini       ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
18941dd7afcfSStefano Zampini     }
1895015636ebSStefano Zampini     ierr = PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);CHKERRQ(ierr);
1896efc2fbd9SStefano Zampini   }
18971f4df5f7SStefano Zampini 
18981dd7afcfSStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1899f913dca9SStefano Zampini     pcbddc->work_change = r;
19001dd7afcfSStefano Zampini     ierr = VecCopy(z,pcbddc->work_change);CHKERRQ(ierr);
19011dd7afcfSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,z);CHKERRQ(ierr);
19021dd7afcfSStefano Zampini   }
19030c7d97c5SJed Brown   PetscFunctionReturn(0);
19040c7d97c5SJed Brown }
190550efa1b5SStefano Zampini 
190650efa1b5SStefano Zampini /*
190750efa1b5SStefano Zampini    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
190850efa1b5SStefano Zampini 
190950efa1b5SStefano Zampini    Input Parameters:
19100f202f7eSStefano Zampini +  pc - the preconditioner context
19110f202f7eSStefano Zampini -  r - input vector (global)
191250efa1b5SStefano Zampini 
191350efa1b5SStefano Zampini    Output Parameter:
191450efa1b5SStefano Zampini .  z - output vector (global)
191550efa1b5SStefano Zampini 
191650efa1b5SStefano Zampini    Application Interface Routine: PCApplyTranspose()
191750efa1b5SStefano Zampini  */
191850efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
191950efa1b5SStefano Zampini {
192050efa1b5SStefano Zampini   PC_IS             *pcis = (PC_IS*)(pc->data);
192150efa1b5SStefano Zampini   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1922b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
192350efa1b5SStefano Zampini   PetscErrorCode    ierr;
192450efa1b5SStefano Zampini   const PetscScalar one = 1.0;
192550efa1b5SStefano Zampini   const PetscScalar m_one = -1.0;
192650efa1b5SStefano Zampini   const PetscScalar zero = 0.0;
192750efa1b5SStefano Zampini 
192850efa1b5SStefano Zampini   PetscFunctionBegin;
1929f3d41395Sstefano_zampini   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
19301dd7afcfSStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
19311dd7afcfSStefano Zampini     Vec swap;
193227b6a85dSStefano Zampini 
193327b6a85dSStefano Zampini     ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,r,pcbddc->work_change);CHKERRQ(ierr);
19341dd7afcfSStefano Zampini     swap = pcbddc->work_change;
19351dd7afcfSStefano Zampini     pcbddc->work_change = r;
19361dd7afcfSStefano Zampini     r = swap;
193727b6a85dSStefano Zampini     /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */
19388ae0ca82SStefano Zampini     if (pcbddc->benign_apply_coarse_only && pcbddc->exact_dirichlet_trick_app && pcbddc->change_interior) {
193927b6a85dSStefano Zampini       ierr = VecCopy(r,pcis->vec1_global);CHKERRQ(ierr);
194027b6a85dSStefano Zampini       ierr = VecLockPush(pcis->vec1_global);CHKERRQ(ierr);
19411dd7afcfSStefano Zampini     }
194227b6a85dSStefano Zampini   }
194327b6a85dSStefano Zampini   if (pcbddc->benign_have_null) { /* get p0 from r */
1944537c1cdfSStefano Zampini     ierr = PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);CHKERRQ(ierr);
1945537c1cdfSStefano Zampini   }
19468ae0ca82SStefano Zampini   if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
1947b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
194850efa1b5SStefano Zampini     /* First Dirichlet solve */
194950efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
195050efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
195150efa1b5SStefano Zampini     /*
195250efa1b5SStefano Zampini       Assembling right hand side for BDDC operator
1953b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
195450efa1b5SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
195550efa1b5SStefano Zampini     */
1956b097fa66SStefano Zampini     if (n_D) {
1957b097fa66SStefano Zampini       ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
195850efa1b5SStefano Zampini       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
195916909a7fSStefano Zampini       if (pcbddc->switch_static) {
196016909a7fSStefano Zampini         Mat_IS *matis = (Mat_IS*)(pc->mat->data);
196116909a7fSStefano Zampini 
196216909a7fSStefano Zampini         ierr = VecSet(pcis->vec1_N,0.);CHKERRQ(ierr);
196316909a7fSStefano Zampini         ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
196416909a7fSStefano Zampini         ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
196516909a7fSStefano Zampini         if (!pcbddc->switch_static_change) {
196616909a7fSStefano Zampini           ierr = MatMultTranspose(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
196716909a7fSStefano Zampini         } else {
196816909a7fSStefano Zampini           ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
196916909a7fSStefano Zampini           ierr = MatMultTranspose(matis->A,pcis->vec2_N,pcis->vec1_N);CHKERRQ(ierr);
197016909a7fSStefano Zampini           ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
197116909a7fSStefano Zampini         }
197216909a7fSStefano Zampini         ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
197316909a7fSStefano Zampini         ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
197416909a7fSStefano Zampini         ierr = VecScatterBegin(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
197516909a7fSStefano Zampini         ierr = VecScatterEnd(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
197616909a7fSStefano Zampini       } else {
1977b097fa66SStefano Zampini         ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
197816909a7fSStefano Zampini       }
1979b097fa66SStefano Zampini     } else {
1980b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1981b097fa66SStefano Zampini     }
198250efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
198350efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
198450efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
198550efa1b5SStefano Zampini   } else {
198650efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
198750efa1b5SStefano Zampini   }
198850efa1b5SStefano Zampini 
198950efa1b5SStefano Zampini   /* Apply interface preconditioner
199050efa1b5SStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1991dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
199250efa1b5SStefano Zampini 
199350efa1b5SStefano Zampini   /* Apply transpose of partition of unity operator */
199450efa1b5SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
199550efa1b5SStefano Zampini 
199650efa1b5SStefano Zampini   /* Second Dirichlet solve and assembling of output */
199750efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
199850efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1999b097fa66SStefano Zampini   if (n_B) {
200016909a7fSStefano Zampini     if (pcbddc->switch_static) {
200116909a7fSStefano Zampini       Mat_IS *matis = (Mat_IS*)(pc->mat->data);
200216909a7fSStefano Zampini 
200316909a7fSStefano Zampini       ierr = VecScatterBegin(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
200416909a7fSStefano Zampini       ierr = VecScatterEnd(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
200516909a7fSStefano Zampini       ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
200616909a7fSStefano Zampini       ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
200716909a7fSStefano Zampini       if (!pcbddc->switch_static_change) {
200816909a7fSStefano Zampini         ierr = MatMultTranspose(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
200916909a7fSStefano Zampini       } else {
201016909a7fSStefano Zampini         ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
201116909a7fSStefano Zampini         ierr = MatMultTranspose(matis->A,pcis->vec2_N,pcis->vec1_N);CHKERRQ(ierr);
201216909a7fSStefano Zampini         ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
201316909a7fSStefano Zampini       }
201416909a7fSStefano Zampini       ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
201516909a7fSStefano Zampini       ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
201616909a7fSStefano Zampini     } else {
201750efa1b5SStefano Zampini       ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
201816909a7fSStefano Zampini     }
201916909a7fSStefano Zampini   } else if (pcbddc->switch_static) { /* n_B is zero */
202016909a7fSStefano Zampini     Mat_IS *matis = (Mat_IS*)(pc->mat->data);
202116909a7fSStefano Zampini 
202216909a7fSStefano Zampini     if (!pcbddc->switch_static_change) {
202316909a7fSStefano Zampini       ierr = MatMultTranspose(matis->A,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
202416909a7fSStefano Zampini     } else {
202516909a7fSStefano Zampini       ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_D,pcis->vec1_N);CHKERRQ(ierr);
202616909a7fSStefano Zampini       ierr = MatMultTranspose(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
202716909a7fSStefano Zampini       ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec2_N,pcis->vec3_D);CHKERRQ(ierr);
202816909a7fSStefano Zampini     }
2029b097fa66SStefano Zampini   }
2030b0147a47SStefano Zampini   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
20318ae0ca82SStefano Zampini   if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
2032b097fa66SStefano Zampini     if (pcbddc->switch_static) {
2033b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
2034b097fa66SStefano Zampini     } else {
2035b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
2036b097fa66SStefano Zampini     }
203750efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
203850efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2039b097fa66SStefano Zampini   } else {
2040b097fa66SStefano Zampini     if (pcbddc->switch_static) {
2041b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
2042b097fa66SStefano Zampini     } else {
2043b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
2044b097fa66SStefano Zampini     }
2045b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2046b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2047b097fa66SStefano Zampini   }
204827b6a85dSStefano Zampini   if (pcbddc->benign_have_null) { /* set p0 (computed in PCBDDCApplyInterface) */
2049537c1cdfSStefano Zampini     ierr = PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);CHKERRQ(ierr);
2050537c1cdfSStefano Zampini   }
20511dd7afcfSStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
2052f913dca9SStefano Zampini     pcbddc->work_change = r;
20531dd7afcfSStefano Zampini     ierr = VecCopy(z,pcbddc->work_change);CHKERRQ(ierr);
20541dd7afcfSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,z);CHKERRQ(ierr);
20551dd7afcfSStefano Zampini   }
205650efa1b5SStefano Zampini   PetscFunctionReturn(0);
205750efa1b5SStefano Zampini }
2058674ae819SStefano Zampini 
20599326c5c6Sstefano_zampini PetscErrorCode PCReset_BDDC(PC pc)
2060da1bb401SStefano Zampini {
2061da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
20629326c5c6Sstefano_zampini   PC_IS          *pcis = (PC_IS*)pc->data;
20639326c5c6Sstefano_zampini   KSP            kspD,kspR,kspC;
2064da1bb401SStefano Zampini   PetscErrorCode ierr;
2065da1bb401SStefano Zampini 
2066da1bb401SStefano Zampini   PetscFunctionBegin;
2067674ae819SStefano Zampini   /* free BDDC custom data  */
2068674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
2069674ae819SStefano Zampini   /* destroy objects related to topography */
2070674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
207134a97f8cSStefano Zampini   /* destroy objects for scaling operator */
2072674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
2073674ae819SStefano Zampini   /* free solvers stuff */
2074674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
207562a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
207662a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
207762a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
20781dd7afcfSStefano Zampini   /* free data created by PCIS */
20791dd7afcfSStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
20809326c5c6Sstefano_zampini 
20819326c5c6Sstefano_zampini   /* restore defaults */
20829326c5c6Sstefano_zampini   kspD = pcbddc->ksp_D;
20839326c5c6Sstefano_zampini   kspR = pcbddc->ksp_R;
20849326c5c6Sstefano_zampini   kspC = pcbddc->coarse_ksp;
20859326c5c6Sstefano_zampini   ierr = PetscMemzero(pc->data,sizeof(*pcbddc));CHKERRQ(ierr);
20869326c5c6Sstefano_zampini   pcis->n_neigh                     = -1;
20879326c5c6Sstefano_zampini   pcis->scaling_factor              = 1.0;
20889326c5c6Sstefano_zampini   pcis->reusesubmatrices            = PETSC_TRUE;
20899326c5c6Sstefano_zampini   pcbddc->use_local_adj             = PETSC_TRUE;
20909326c5c6Sstefano_zampini   pcbddc->use_vertices              = PETSC_TRUE;
20919326c5c6Sstefano_zampini   pcbddc->use_edges                 = PETSC_TRUE;
20929326c5c6Sstefano_zampini   pcbddc->symmetric_primal          = PETSC_TRUE;
20939326c5c6Sstefano_zampini   pcbddc->vertex_size               = 1;
20949326c5c6Sstefano_zampini   pcbddc->recompute_topography      = PETSC_TRUE;
20959326c5c6Sstefano_zampini   pcbddc->coarse_size               = -1;
20969326c5c6Sstefano_zampini   pcbddc->use_exact_dirichlet_trick = PETSC_TRUE;
20979326c5c6Sstefano_zampini   pcbddc->coarsening_ratio          = 8;
20989326c5c6Sstefano_zampini   pcbddc->coarse_eqs_per_proc       = 1;
20999326c5c6Sstefano_zampini   pcbddc->benign_compute_correction = PETSC_TRUE;
21009326c5c6Sstefano_zampini   pcbddc->nedfield                  = -1;
21019326c5c6Sstefano_zampini   pcbddc->nedglobal                 = PETSC_TRUE;
21029326c5c6Sstefano_zampini   pcbddc->graphmaxcount             = PETSC_MAX_INT;
21039326c5c6Sstefano_zampini   pcbddc->sub_schurs_layers         = -1;
21049326c5c6Sstefano_zampini   pcbddc->ksp_D                     = kspD;
21059326c5c6Sstefano_zampini   pcbddc->ksp_R                     = kspR;
21069326c5c6Sstefano_zampini   pcbddc->coarse_ksp                = kspC;
21079326c5c6Sstefano_zampini   PetscFunctionReturn(0);
21089326c5c6Sstefano_zampini }
21099326c5c6Sstefano_zampini 
21109326c5c6Sstefano_zampini PetscErrorCode PCDestroy_BDDC(PC pc)
21119326c5c6Sstefano_zampini {
21129326c5c6Sstefano_zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
21139326c5c6Sstefano_zampini   PetscErrorCode ierr;
21149326c5c6Sstefano_zampini 
21159326c5c6Sstefano_zampini   PetscFunctionBegin;
21169326c5c6Sstefano_zampini   ierr = PCReset_BDDC(pc);CHKERRQ(ierr);
21179326c5c6Sstefano_zampini   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
21189326c5c6Sstefano_zampini   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
21199326c5c6Sstefano_zampini   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
2120a13144ffSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDiscreteGradient_C",NULL);CHKERRQ(ierr);
2121a198735bSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDivergenceMat_C",NULL);CHKERRQ(ierr);
2122906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
2123674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
212430368db7SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",NULL);CHKERRQ(ierr);
2125bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
21262b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
2127b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
21282b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
2129bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
213082ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
2131bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
213282ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
2133bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
213482ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
2135bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
2136785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
2137bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
213863602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
2139bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
2140bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
2141bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
2142bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
2143a06fd7f2SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);CHKERRQ(ierr);
2144ab8c8b98SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
2145674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2146da1bb401SStefano Zampini   PetscFunctionReturn(0);
2147da1bb401SStefano Zampini }
21481e6b0712SBarry Smith 
2149ab8c8b98SStefano Zampini static PetscErrorCode PCSetCoordinates_BDDC(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
2150ab8c8b98SStefano Zampini {
2151ab8c8b98SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
2152ab8c8b98SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
2153ab8c8b98SStefano Zampini   PetscErrorCode ierr;
2154ab8c8b98SStefano Zampini 
2155ab8c8b98SStefano Zampini   PetscFunctionBegin;
2156ab8c8b98SStefano Zampini   ierr = PetscFree(mat_graph->coords);CHKERRQ(ierr);
2157ab8c8b98SStefano Zampini   ierr = PetscMalloc1(nloc*dim,&mat_graph->coords);CHKERRQ(ierr);
2158ab8c8b98SStefano Zampini   ierr = PetscMemcpy(mat_graph->coords,coords,nloc*dim*sizeof(PetscReal));CHKERRQ(ierr);
2159ab8c8b98SStefano Zampini   mat_graph->cnloc = nloc;
2160ab8c8b98SStefano Zampini   mat_graph->cdim  = dim;
2161ab8c8b98SStefano Zampini   mat_graph->cloc  = PETSC_FALSE;
2162ab8c8b98SStefano Zampini   PetscFunctionReturn(0);
2163ab8c8b98SStefano Zampini }
2164ab8c8b98SStefano Zampini 
2165a06fd7f2SStefano Zampini static PetscErrorCode PCPreSolveChangeRHS_BDDC(PC pc, PetscBool* change)
2166a06fd7f2SStefano Zampini {
2167a06fd7f2SStefano Zampini   PetscFunctionBegin;
2168a06fd7f2SStefano Zampini   *change = PETSC_TRUE;
2169a06fd7f2SStefano Zampini   PetscFunctionReturn(0);
2170a06fd7f2SStefano Zampini }
2171a06fd7f2SStefano Zampini 
21723425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
21733425bc38SStefano Zampini {
2174674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
2175266e20e9SStefano Zampini   Vec            work;
21763425bc38SStefano Zampini   PC_IS*         pcis;
21773425bc38SStefano Zampini   PC_BDDC*       pcbddc;
21783425bc38SStefano Zampini   PetscErrorCode ierr;
21790c7d97c5SJed Brown 
21803425bc38SStefano Zampini   PetscFunctionBegin;
21813425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
21823425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
21833425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
21843425bc38SStefano Zampini 
2185229984c5Sstefano_zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
2186229984c5Sstefano_zampini   /* copy rhs since we may change it during PCPreSolve_BDDC */
2187229984c5Sstefano_zampini   if (!pcbddc->original_rhs) {
2188229984c5Sstefano_zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
2189229984c5Sstefano_zampini   }
21906cc1294bSstefano_zampini   if (mat_ctx->rhs_flip) {
21916cc1294bSstefano_zampini     ierr = VecPointwiseMult(pcbddc->original_rhs,standard_rhs,mat_ctx->rhs_flip);CHKERRQ(ierr);
21926cc1294bSstefano_zampini   } else {
2193229984c5Sstefano_zampini     ierr = VecCopy(standard_rhs,pcbddc->original_rhs);CHKERRQ(ierr);
21946cc1294bSstefano_zampini   }
2195af140850Sstefano_zampini   if (mat_ctx->g2g_p) {
2196229984c5Sstefano_zampini     /* interface pressure rhs */
2197022d8d2bSstefano_zampini     ierr = VecScatterBegin(mat_ctx->g2g_p,fetidp_flux_rhs,pcbddc->original_rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2198022d8d2bSstefano_zampini     ierr = VecScatterEnd(mat_ctx->g2g_p,fetidp_flux_rhs,pcbddc->original_rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2199229984c5Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->g2g_p,standard_rhs,fetidp_flux_rhs,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2200229984c5Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->g2g_p,standard_rhs,fetidp_flux_rhs,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
22016cc1294bSstefano_zampini     if (!mat_ctx->rhs_flip) {
2202229984c5Sstefano_zampini       ierr = VecScale(fetidp_flux_rhs,-1.);CHKERRQ(ierr);
2203229984c5Sstefano_zampini     }
22046cc1294bSstefano_zampini   }
2205c08af4c6SStefano Zampini   /*
2206c08af4c6SStefano Zampini      change of basis for physical rhs if needed
2207c08af4c6SStefano Zampini      It also changes the rhs in case of dirichlet boundaries
2208c08af4c6SStefano Zampini   */
22093738a8e6SStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,pcbddc->original_rhs,NULL);CHKERRQ(ierr);
2210fc17d649SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
22113738a8e6SStefano Zampini     ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcbddc->original_rhs,pcbddc->work_change);CHKERRQ(ierr);
22123738a8e6SStefano Zampini     work = pcbddc->work_change;
2213fc17d649SStefano Zampini    } else {
22143738a8e6SStefano Zampini     work = pcbddc->original_rhs;
2215fc17d649SStefano Zampini   }
22163425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
2217266e20e9SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,work,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2218266e20e9SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,work,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2219fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
2220fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
2221266e20e9SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2222266e20e9SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2223674ae819SStefano Zampini   /* Apply partition of unity */
22243425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
2225266e20e9SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,work,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
22268eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
22273425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
22283425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
22293425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
22303425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
2231266e20e9SStefano Zampini     ierr = VecSet(work,0.0);CHKERRQ(ierr);
2232266e20e9SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,work,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2233266e20e9SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,work,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2234266e20e9SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,work,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
2235266e20e9SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2236266e20e9SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
22373425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
22383425bc38SStefano Zampini   }
22393425bc38SStefano Zampini   /* BDDC rhs */
22403425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
22418eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
22423425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
22433425bc38SStefano Zampini   }
22443425bc38SStefano Zampini   /* apply BDDC */
2245229984c5Sstefano_zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
2246dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
2247266e20e9SStefano Zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
2248229984c5Sstefano_zampini 
22493425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
22503425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
22513425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
22523425bc38SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2253229984c5Sstefano_zampini   /* Add contribution to interface pressures */
2254229984c5Sstefano_zampini   if (mat_ctx->l2g_p) {
2255229984c5Sstefano_zampini     ierr = MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr);
2256229984c5Sstefano_zampini     if (pcbddc->switch_static) {
2257229984c5Sstefano_zampini       ierr = MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr);
2258229984c5Sstefano_zampini     }
2259229984c5Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2260229984c5Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2261229984c5Sstefano_zampini   }
22623425bc38SStefano Zampini   PetscFunctionReturn(0);
22633425bc38SStefano Zampini }
22641e6b0712SBarry Smith 
22653425bc38SStefano Zampini /*@
22660f202f7eSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side
22673425bc38SStefano Zampini 
22683425bc38SStefano Zampini    Collective
22693425bc38SStefano Zampini 
22703425bc38SStefano Zampini    Input Parameters:
22710f202f7eSStefano Zampini +  fetidp_mat      - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
22720f202f7eSStefano Zampini -  standard_rhs    - the right-hand side of the original linear system
22733425bc38SStefano Zampini 
22743425bc38SStefano Zampini    Output Parameters:
22750f202f7eSStefano Zampini .  fetidp_flux_rhs - the right-hand side for the FETI-DP linear system
22763425bc38SStefano Zampini 
22773425bc38SStefano Zampini    Level: developer
22783425bc38SStefano Zampini 
22793425bc38SStefano Zampini    Notes:
22803425bc38SStefano Zampini 
22810f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
22823425bc38SStefano Zampini @*/
22833425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
22843425bc38SStefano Zampini {
2285674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
22863425bc38SStefano Zampini   PetscErrorCode ierr;
22873425bc38SStefano Zampini 
22883425bc38SStefano Zampini   PetscFunctionBegin;
2289266e20e9SStefano Zampini   PetscValidHeaderSpecific(fetidp_mat,MAT_CLASSID,1);
2290266e20e9SStefano Zampini   PetscValidHeaderSpecific(standard_rhs,VEC_CLASSID,2);
2291266e20e9SStefano Zampini   PetscValidHeaderSpecific(fetidp_flux_rhs,VEC_CLASSID,3);
22923425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
2293163d334eSBarry Smith   ierr = PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
22943425bc38SStefano Zampini   PetscFunctionReturn(0);
22953425bc38SStefano Zampini }
22961e6b0712SBarry Smith 
22973425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
22983425bc38SStefano Zampini {
2299674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
23003425bc38SStefano Zampini   PC_IS*         pcis;
23013425bc38SStefano Zampini   PC_BDDC*       pcbddc;
23023425bc38SStefano Zampini   PetscErrorCode ierr;
2303229984c5Sstefano_zampini   Vec            work;
23043425bc38SStefano Zampini 
23053425bc38SStefano Zampini   PetscFunctionBegin;
23063425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
23073425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
23083425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
23093425bc38SStefano Zampini 
23103425bc38SStefano Zampini   /* apply B_delta^T */
2311af140850Sstefano_zampini   ierr = VecSet(pcis->vec1_B,0.);CHKERRQ(ierr);
23123425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23133425bc38SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23143425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
2315229984c5Sstefano_zampini   if (mat_ctx->l2g_p) {
2316229984c5Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->l2g_p,fetidp_flux_sol,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2317229984c5Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->l2g_p,fetidp_flux_sol,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2318229984c5Sstefano_zampini     ierr = MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
2319229984c5Sstefano_zampini   }
2320229984c5Sstefano_zampini 
23213425bc38SStefano Zampini   /* compute rhs for BDDC application */
23223425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
23238eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
23243425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
2325229984c5Sstefano_zampini     if (mat_ctx->l2g_p) {
2326229984c5Sstefano_zampini       ierr = VecScale(mat_ctx->vP,-1.);CHKERRQ(ierr);
2327229984c5Sstefano_zampini       ierr = MatMultAdd(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr);
23283425bc38SStefano Zampini     }
2329229984c5Sstefano_zampini   }
2330229984c5Sstefano_zampini 
23313425bc38SStefano Zampini   /* apply BDDC */
2332229984c5Sstefano_zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
2333dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
2334229984c5Sstefano_zampini 
2335229984c5Sstefano_zampini   /* put values into global vector */
2336af140850Sstefano_zampini   if (pcbddc->ChangeOfBasisMatrix) work = pcbddc->work_change;
2337af140850Sstefano_zampini   else work = standard_sol;
2338229984c5Sstefano_zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2339229984c5Sstefano_zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23408eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
23413425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
23423425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
234300f6b531SStefano Zampini     ierr = VecAYPX(pcis->vec1_D,-1.0,mat_ctx->temp_solution_D);CHKERRQ(ierr);
234400f6b531SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr);
23453425bc38SStefano Zampini   }
2346229984c5Sstefano_zampini 
2347229984c5Sstefano_zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2348229984c5Sstefano_zampini   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec1_D,work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2349266e20e9SStefano Zampini   /* add p0 solution to final solution */
2350229984c5Sstefano_zampini   ierr = PCBDDCBenignGetOrSetP0(mat_ctx->pc,work,PETSC_FALSE);CHKERRQ(ierr);
2351fc17d649SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
2352af140850Sstefano_zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,work,standard_sol);CHKERRQ(ierr);
2353fc17d649SStefano Zampini   }
2354af140850Sstefano_zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
2355af140850Sstefano_zampini   if (mat_ctx->g2g_p) {
2356229984c5Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->g2g_p,fetidp_flux_sol,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2357229984c5Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->g2g_p,fetidp_flux_sol,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2358229984c5Sstefano_zampini   }
23593425bc38SStefano Zampini   PetscFunctionReturn(0);
23603425bc38SStefano Zampini }
23611e6b0712SBarry Smith 
23623425bc38SStefano Zampini /*@
23630f202f7eSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system
23643425bc38SStefano Zampini 
23653425bc38SStefano Zampini    Collective
23663425bc38SStefano Zampini 
23673425bc38SStefano Zampini    Input Parameters:
23680f202f7eSStefano Zampini +  fetidp_mat      - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
23690f202f7eSStefano Zampini -  fetidp_flux_sol - the solution of the FETI-DP linear system
23703425bc38SStefano Zampini 
23713425bc38SStefano Zampini    Output Parameters:
23720f202f7eSStefano Zampini .  standard_sol    - the solution defined on the physical domain
23733425bc38SStefano Zampini 
23743425bc38SStefano Zampini    Level: developer
23753425bc38SStefano Zampini 
23763425bc38SStefano Zampini    Notes:
23773425bc38SStefano Zampini 
23780f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
23793425bc38SStefano Zampini @*/
23803425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
23813425bc38SStefano Zampini {
2382674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
23833425bc38SStefano Zampini   PetscErrorCode ierr;
23843425bc38SStefano Zampini 
23853425bc38SStefano Zampini   PetscFunctionBegin;
2386266e20e9SStefano Zampini   PetscValidHeaderSpecific(fetidp_mat,MAT_CLASSID,1);
2387266e20e9SStefano Zampini   PetscValidHeaderSpecific(fetidp_flux_sol,VEC_CLASSID,2);
2388266e20e9SStefano Zampini   PetscValidHeaderSpecific(standard_sol,VEC_CLASSID,3);
23893425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
2390163d334eSBarry Smith   ierr = PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
23913425bc38SStefano Zampini   PetscFunctionReturn(0);
23923425bc38SStefano Zampini }
23931e6b0712SBarry Smith 
2394547c9a8eSstefano_zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, PetscBool fully_redundant, const char* prefix, Mat *fetidp_mat, PC *fetidp_pc)
23953425bc38SStefano Zampini {
2396674ae819SStefano Zampini 
2397674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
23983425bc38SStefano Zampini   Mat            newmat;
2399674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
24003425bc38SStefano Zampini   PC             newpc;
2401ce94432eSBarry Smith   MPI_Comm       comm;
24023425bc38SStefano Zampini   PetscErrorCode ierr;
24033425bc38SStefano Zampini 
24043425bc38SStefano Zampini   PetscFunctionBegin;
2405ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
24063425bc38SStefano Zampini   /* FETIDP linear matrix */
24073425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
24081720468bSStefano Zampini   fetidpmat_ctx->fully_redundant = fully_redundant;
24093425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
2410a5bb87b3Sstefano_zampini   ierr = MatCreateShell(comm,fetidpmat_ctx->n,fetidpmat_ctx->n,fetidpmat_ctx->N,fetidpmat_ctx->N,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
24113425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
2412edf7251bSStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
24133425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
2414547c9a8eSstefano_zampini   ierr = MatSetOptionsPrefix(newmat,prefix);CHKERRQ(ierr);
2415547c9a8eSstefano_zampini   ierr = MatAppendOptionsPrefix(newmat,"fetidp_");CHKERRQ(ierr);
24163425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
24173425bc38SStefano Zampini   /* FETIDP preconditioner */
24183425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
24193425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
24203425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
2421e1214c54Sstefano_zampini   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
2422e1214c54Sstefano_zampini   ierr = PCSetOptionsPrefix(newpc,prefix);CHKERRQ(ierr);
2423e1214c54Sstefano_zampini   ierr = PCAppendOptionsPrefix(newpc,"fetidp_");CHKERRQ(ierr);
2424e1214c54Sstefano_zampini   if (!fetidpmat_ctx->l2g_lambda_only) {
24253425bc38SStefano Zampini     ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
24263425bc38SStefano Zampini     ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
24273425bc38SStefano Zampini     ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
2428edf7251bSStefano Zampini     ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
2429c45b8d2dSstefano_zampini     ierr = PCShellSetView(newpc,FETIDPPCView);CHKERRQ(ierr);
24303425bc38SStefano Zampini     ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
2431e1214c54Sstefano_zampini   } else {
2432e1214c54Sstefano_zampini     KSP      *ksps;
2433e1214c54Sstefano_zampini     PC       lagpc;
2434e1214c54Sstefano_zampini     Mat      M,AM,PM;
2435e1214c54Sstefano_zampini     PetscInt nn;
2436e1214c54Sstefano_zampini 
2437e1214c54Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_PPmat",(PetscObject*)&M);CHKERRQ(ierr);
2438e1214c54Sstefano_zampini     ierr = PCSetType(newpc,PCFIELDSPLIT);CHKERRQ(ierr);
2439e1214c54Sstefano_zampini     ierr = PCFieldSplitSetIS(newpc,"lag",fetidpmat_ctx->lagrange);CHKERRQ(ierr);
2440e1214c54Sstefano_zampini     ierr = PCFieldSplitSetIS(newpc,"p",fetidpmat_ctx->pressure);CHKERRQ(ierr);
2441e1214c54Sstefano_zampini     ierr = PCFieldSplitSetType(newpc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
244240c75d76SStefano Zampini     ierr = PCFieldSplitSetSchurFactType(newpc,PC_FIELDSPLIT_SCHUR_FACT_DIAG);CHKERRQ(ierr);
2443e1214c54Sstefano_zampini     ierr = PCFieldSplitSetSchurPre(newpc,PC_FIELDSPLIT_SCHUR_PRE_USER,M);CHKERRQ(ierr);
2444c096484dSStefano Zampini     ierr = PCFieldSplitSetSchurScale(newpc,1.0);CHKERRQ(ierr);
2445e1214c54Sstefano_zampini     ierr = PCSetFromOptions(newpc);CHKERRQ(ierr);
2446e1214c54Sstefano_zampini     ierr = PCSetUp(newpc);CHKERRQ(ierr);
2447e1214c54Sstefano_zampini 
2448e1214c54Sstefano_zampini     /* set the solver for the (0,0) block */
244940c75d76SStefano Zampini     ierr = PCFieldSplitGetSubKSP(newpc,&nn,&ksps);CHKERRQ(ierr);
2450e1214c54Sstefano_zampini     ierr = PCCreate(comm,&lagpc);CHKERRQ(ierr);
2451e1214c54Sstefano_zampini     ierr = PCSetType(lagpc,PCSHELL);CHKERRQ(ierr);
2452e1214c54Sstefano_zampini     ierr = KSPGetOperators(ksps[0],&AM,&PM);CHKERRQ(ierr);
2453e1214c54Sstefano_zampini     ierr = PCSetOperators(lagpc,AM,PM);CHKERRQ(ierr);
2454e1214c54Sstefano_zampini     ierr = PCShellSetContext(lagpc,fetidppc_ctx);CHKERRQ(ierr);
2455e1214c54Sstefano_zampini     ierr = PCShellSetApply(lagpc,FETIDPPCApply);CHKERRQ(ierr);
2456e1214c54Sstefano_zampini     ierr = PCShellSetApplyTranspose(lagpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
2457e1214c54Sstefano_zampini     ierr = PCShellSetView(lagpc,FETIDPPCView);CHKERRQ(ierr);
2458e1214c54Sstefano_zampini     ierr = PCShellSetDestroy(lagpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
2459e1214c54Sstefano_zampini     ierr = KSPSetPC(ksps[0],lagpc);CHKERRQ(ierr);
2460e1214c54Sstefano_zampini     ierr = PCDestroy(&lagpc);CHKERRQ(ierr);
2461e1214c54Sstefano_zampini     ierr = PetscFree(ksps);CHKERRQ(ierr);
2462e1214c54Sstefano_zampini   }
24633425bc38SStefano Zampini   /* return pointers for objects created */
24643425bc38SStefano Zampini   *fetidp_mat = newmat;
24653425bc38SStefano Zampini   *fetidp_pc = newpc;
24663425bc38SStefano Zampini   PetscFunctionReturn(0);
24673425bc38SStefano Zampini }
24681e6b0712SBarry Smith 
246994ef8ddeSSatish Balay /*@C
24700f202f7eSStefano Zampini  PCBDDCCreateFETIDPOperators - Create FETI-DP operators
24713425bc38SStefano Zampini 
24723425bc38SStefano Zampini    Collective
24733425bc38SStefano Zampini 
24743425bc38SStefano Zampini    Input Parameters:
24751720468bSStefano Zampini +  pc - the BDDC preconditioning context (setup should have been called before)
2476547c9a8eSstefano_zampini .  fully_redundant - true for a fully redundant set of Lagrange multipliers
2477547c9a8eSstefano_zampini -  prefix - optional options database prefix for the objects to be created (can be NULL)
247828509bceSStefano Zampini 
247928509bceSStefano Zampini    Output Parameters:
24800f202f7eSStefano Zampini +  fetidp_mat - shell FETI-DP matrix object
24810f202f7eSStefano Zampini -  fetidp_pc  - shell Dirichlet preconditioner for FETI-DP matrix
248228509bceSStefano Zampini 
24833425bc38SStefano Zampini    Level: developer
24843425bc38SStefano Zampini 
24853425bc38SStefano Zampini    Notes:
24860f202f7eSStefano Zampini      Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose
24873425bc38SStefano Zampini 
24880f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
24893425bc38SStefano Zampini @*/
2490547c9a8eSstefano_zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, PetscBool fully_redundant, const char *prefix, Mat *fetidp_mat, PC *fetidp_pc)
24913425bc38SStefano Zampini {
24923425bc38SStefano Zampini   PetscErrorCode ierr;
24933425bc38SStefano Zampini 
24943425bc38SStefano Zampini   PetscFunctionBegin;
24953425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
24963425bc38SStefano Zampini   if (pc->setupcalled) {
2497547c9a8eSstefano_zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,PetscBool,const char*,Mat*,PC*),(pc,fully_redundant,prefix,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
24984247aa23Sstefano_zampini   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
24993425bc38SStefano Zampini   PetscFunctionReturn(0);
25003425bc38SStefano Zampini }
25010c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
2502da1bb401SStefano Zampini /*MC
2503da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
25040c7d97c5SJed Brown 
250528509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
250628509bceSStefano Zampini 
250728509bceSStefano Zampini .vb
250828509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
250928509bceSStefano Zampini    [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
251028509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
25110f202f7eSStefano Zampini    [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
251228509bceSStefano Zampini .ve
251328509bceSStefano Zampini 
251428509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
251528509bceSStefano Zampini 
25160f202f7eSStefano Zampini    Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
251728509bceSStefano Zampini 
251828509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
251928509bceSStefano Zampini 
2520b6fdb6dfSStefano Zampini    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.
2521b6fdb6dfSStefano Zampini 
2522c7017625SStefano Zampini    Approximate local solvers are automatically adapted (see [1]) if the user has attached a nullspace object to the subdomain matrices, and informed BDDC of using approximate solvers (via the command line).
252328509bceSStefano Zampini 
25240f202f7eSStefano Zampini    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()
252530368db7SStefano Zampini    Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesIS() and their local counterparts.
252628509bceSStefano Zampini 
25270f202f7eSStefano Zampini    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.
252828509bceSStefano Zampini 
25290f202f7eSStefano Zampini    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.
25300f202f7eSStefano Zampini    User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()
253128509bceSStefano Zampini 
25320f202f7eSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.
253328509bceSStefano Zampini 
2534df4d28bfSStefano Zampini    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.
253528509bceSStefano Zampini 
25360f202f7eSStefano Zampini    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.
25370f202f7eSStefano Zampini    Deluxe scaling is not supported yet for FETI-DP.
25380f202f7eSStefano Zampini 
25390f202f7eSStefano Zampini    Options Database Keys (some of them, run with -h for a complete list):
25400f202f7eSStefano Zampini 
25410f202f7eSStefano Zampini .    -pc_bddc_use_vertices <true> - use or not vertices in primal space
25420f202f7eSStefano Zampini .    -pc_bddc_use_edges <true> - use or not edges in primal space
25430f202f7eSStefano Zampini .    -pc_bddc_use_faces <false> - use or not faces in primal space
25440f202f7eSStefano Zampini .    -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
25450f202f7eSStefano Zampini .    -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
25460f202f7eSStefano Zampini .    -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
25470f202f7eSStefano Zampini .    -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
254828509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
25490f202f7eSStefano Zampini .    -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)
25505459c157SBarry Smith .    -pc_bddc_coarse_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)
25510f202f7eSStefano Zampini .    -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
25520f202f7eSStefano Zampini .    -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)
2553bd2a564bSStefano Zampini .    -pc_bddc_adaptive_threshold <0.0> - when a value different than zero is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS or MKL_PARDISO installed)
255428509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
255528509bceSStefano Zampini 
255628509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
255728509bceSStefano Zampini .vb
255828509bceSStefano Zampini       -pc_bddc_dirichlet_
255928509bceSStefano Zampini       -pc_bddc_neumann_
256028509bceSStefano Zampini       -pc_bddc_coarse_
256128509bceSStefano Zampini .ve
25620f202f7eSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.
256328509bceSStefano Zampini 
25640f202f7eSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
256528509bceSStefano Zampini .vb
2566312be037SStefano Zampini       -pc_bddc_dirichlet_lN_
2567312be037SStefano Zampini       -pc_bddc_neumann_lN_
2568312be037SStefano Zampini       -pc_bddc_coarse_lN_
256928509bceSStefano Zampini .ve
25700f202f7eSStefano Zampini    Note that level number ranges from the finest (0) to the coarsest (N).
25710f202f7eSStefano Zampini    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.
25720f202f7eSStefano Zampini .vb
25730f202f7eSStefano Zampini      -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
25740f202f7eSStefano Zampini .ve
25750f202f7eSStefano Zampini    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
2576da1bb401SStefano Zampini 
2577da1bb401SStefano Zampini    Level: intermediate
2578da1bb401SStefano Zampini 
2579b6fdb6dfSStefano Zampini    Developer notes:
2580da1bb401SStefano Zampini 
2581da1bb401SStefano Zampini    Contributed by Stefano Zampini
2582da1bb401SStefano Zampini 
2583da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
2584da1bb401SStefano Zampini M*/
2585b2573a8aSBarry Smith 
25868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
2587da1bb401SStefano Zampini {
2588da1bb401SStefano Zampini   PetscErrorCode      ierr;
2589da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
2590da1bb401SStefano Zampini 
2591da1bb401SStefano Zampini   PetscFunctionBegin;
2592b00a9115SJed Brown   ierr     = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
2593da1bb401SStefano Zampini   pc->data = (void*)pcbddc;
2594da1bb401SStefano Zampini 
2595da1bb401SStefano Zampini   /* create PCIS data structure */
2596da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
2597da1bb401SStefano Zampini 
25989326c5c6Sstefano_zampini   /* create local graph structure */
25999326c5c6Sstefano_zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
26009326c5c6Sstefano_zampini 
26019326c5c6Sstefano_zampini   /* BDDC nonzero defaults */
260208a5cf49SStefano Zampini   pcbddc->use_local_adj             = PETSC_TRUE;
260347d04d0dSStefano Zampini   pcbddc->use_vertices              = PETSC_TRUE;
260447d04d0dSStefano Zampini   pcbddc->use_edges                 = PETSC_TRUE;
26053301b35fSStefano Zampini   pcbddc->symmetric_primal          = PETSC_TRUE;
260614f95afaSStefano Zampini   pcbddc->vertex_size               = 1;
2607c703fcc7SStefano Zampini   pcbddc->recompute_topography      = PETSC_TRUE;
260868457ee5SStefano Zampini   pcbddc->coarse_size               = -1;
260985c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = PETSC_TRUE;
261047d04d0dSStefano Zampini   pcbddc->coarsening_ratio          = 8;
261157de7509SStefano Zampini   pcbddc->coarse_eqs_per_proc       = 1;
261227b6a85dSStefano Zampini   pcbddc->benign_compute_correction = PETSC_TRUE;
26131e0482f5SStefano Zampini   pcbddc->nedfield                  = -1;
26141e0482f5SStefano Zampini   pcbddc->nedglobal                 = PETSC_TRUE;
2615be12c134Sstefano_zampini   pcbddc->graphmaxcount             = PETSC_MAX_INT;
2616b96c3477SStefano Zampini   pcbddc->sub_schurs_layers         = -1;
2617bd2a564bSStefano Zampini   pcbddc->adaptive_threshold[0]     = 0.0;
2618bd2a564bSStefano Zampini   pcbddc->adaptive_threshold[1]     = 0.0;
2619b7eb3628SStefano Zampini 
2620da1bb401SStefano Zampini   /* function pointers */
2621da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
262293bd9ae7SStefano Zampini   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2623da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
2624da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
2625da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
26266b78500eSPatrick Sanan   pc->ops->view                = PCView_BDDC;
2627da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
2628da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
2629da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
2630534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
2631534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
26329326c5c6Sstefano_zampini   pc->ops->reset               = PCReset_BDDC;
2633da1bb401SStefano Zampini 
2634da1bb401SStefano Zampini   /* composing function */
2635a13144ffSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDiscreteGradient_C",PCBDDCSetDiscreteGradient_BDDC);CHKERRQ(ierr);
2636a198735bSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDivergenceMat_C",PCBDDCSetDivergenceMat_BDDC);CHKERRQ(ierr);
2637906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
2638674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
263930368db7SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",PCBDDCSetPrimalVerticesIS_BDDC);CHKERRQ(ierr);
2640bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
26412b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
2642b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
26432b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
2644bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
264582ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2646bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
264782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2648bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
264982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2650bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
265182ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2652bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
265363602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
2654bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
2655bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
2656bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
2657bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
2658a06fd7f2SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_BDDC);CHKERRQ(ierr);
2659ab8c8b98SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_BDDC);CHKERRQ(ierr);
2660da1bb401SStefano Zampini   PetscFunctionReturn(0);
2661da1bb401SStefano Zampini }
2662