Actual source code: bddc.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: /* TODOLIST

  3:    Solvers
  4:    - Add support for cholesky for coarse solver (similar to local solvers)
  5:    - Propagate ksp prefixes for solvers to mat objects?

  7:    User interface
  8:    - ** DM attached to pc?

 10:    Debugging output
 11:    - * Better management of verbosity levels of debugging output

 13:    Extra
 14:    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
 15:    - BDDC with MG framework?

 17:    MATIS related operations contained in BDDC code
 18:    - Provide general case for subassembling

 20: */

 22:  #include <../src/ksp/pc/impls/bddc/bddc.h>
 23:  #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
 24:  #include <petscblaslapack.h>

 26: static PetscBool  cited = PETSC_FALSE;
 27: static const char citation[] =
 28: "@article{ZampiniPCBDDC,\n"
 29: "author = {Stefano Zampini},\n"
 30: "title = {{PCBDDC}: A Class of Robust Dual-Primal Methods in {PETS}c},\n"
 31: "journal = {SIAM Journal on Scientific Computing},\n"
 32: "volume = {38},\n"
 33: "number = {5},\n"
 34: "pages = {S282-S306},\n"
 35: "year = {2016},\n"
 36: "doi = {10.1137/15M1025785},\n"
 37: "URL = {http://dx.doi.org/10.1137/15M1025785},\n"
 38: "eprint = {http://dx.doi.org/10.1137/15M1025785}\n"
 39: "}\n";

 41: /* temporarily declare it */
 42: PetscErrorCode PCApply_BDDC(PC,Vec,Vec);

 44: PetscErrorCode PCSetFromOptions_BDDC(PetscOptionItems *PetscOptionsObject,PC pc)
 45: {
 46:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
 47:   PetscInt       nt;

 51:   PetscOptionsHead(PetscOptionsObject,"BDDC options");
 52:   /* Verbose debugging */
 53:   PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);
 54:   /* Approximate solvers */
 55:   PetscOptionsBool("-pc_bddc_dirichlet_approximate","Inform PCBDDC that we are using approximate Dirichlet solvers","none",pcbddc->NullSpace_corr[0],&pcbddc->NullSpace_corr[0],NULL);
 56:   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);
 57:   PetscOptionsBool("-pc_bddc_neumann_approximate","Inform PCBDDC that we are using approximate Neumann solvers","none",pcbddc->NullSpace_corr[2],&pcbddc->NullSpace_corr[2],NULL);
 58:   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);
 59:   /* Primal space customization */
 60:   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);
 61:   PetscOptionsInt("-pc_bddc_graph_maxcount","Maximum number of shared subdomains for a connected component","none",pcbddc->graphmaxcount,&pcbddc->graphmaxcount,NULL);
 62:   PetscOptionsBool("-pc_bddc_corner_selection","Activates face-based corner selection","none",pcbddc->corner_selection,&pcbddc->corner_selection,NULL);
 63:   PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);
 64:   PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);
 65:   PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);
 66:   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);
 67:   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);
 68:   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);
 69:   /* Change of basis */
 70:   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);
 71:   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);
 72:   if (!pcbddc->use_change_of_basis) {
 73:     pcbddc->use_change_on_faces = PETSC_FALSE;
 74:   }
 75:   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
 76:   PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);
 77:   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);
 78:   PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);
 79:   PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);
 80:   PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);
 81:   PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);
 82:   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);
 83:   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);
 84:   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);
 85:   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);
 86:   PetscOptionsBool("-pc_bddc_deluxe_zerorows","Zero rows and columns of deluxe operators associated with primal dofs","none",pcbddc->deluxe_zerorows,&pcbddc->deluxe_zerorows,NULL);
 87:   PetscOptionsBool("-pc_bddc_deluxe_singlemat","Collapse deluxe operators","none",pcbddc->deluxe_singlemat,&pcbddc->deluxe_singlemat,NULL);
 88:   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);
 89:   nt   = 2;
 90:   PetscOptionsRealArray("-pc_bddc_adaptive_threshold","Thresholds to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&nt,NULL);
 91:   if (nt == 1) pcbddc->adaptive_threshold[1] = pcbddc->adaptive_threshold[0];
 92:   PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);
 93:   PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);
 94:   PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);
 95:   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);
 96:   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);
 97:   PetscOptionsBool("-pc_bddc_benign_change","Compute the pressure change of basis explicitly","none",pcbddc->benign_change_explicit,&pcbddc->benign_change_explicit,NULL);
 98:   PetscOptionsBool("-pc_bddc_benign_compute_correction","Compute the benign correction during PreSolve","none",pcbddc->benign_compute_correction,&pcbddc->benign_compute_correction,NULL);
 99:   PetscOptionsBool("-pc_bddc_nonetflux","Automatic computation of no-net-flux quadrature weights","none",pcbddc->compute_nonetflux,&pcbddc->compute_nonetflux,NULL);
100:   PetscOptionsBool("-pc_bddc_detect_disconnected","Detects disconnected subdomains","none",pcbddc->detect_disconnected,&pcbddc->detect_disconnected,NULL);
101:   PetscOptionsBool("-pc_bddc_eliminate_dirichlet","Whether or not we want to eliminate dirichlet dofs during presolve","none",pcbddc->eliminate_dirdofs,&pcbddc->eliminate_dirdofs,NULL);
102:   PetscOptionsTail();
103:   return(0);
104: }

106: static PetscErrorCode PCView_BDDC(PC pc,PetscViewer viewer)
107: {
108:   PC_BDDC              *pcbddc = (PC_BDDC*)pc->data;
109:   PC_IS                *pcis = (PC_IS*)pc->data;
110:   PetscErrorCode       ierr;
111:   PetscBool            isascii;
112:   PetscSubcomm         subcomm;
113:   PetscViewer          subviewer;

116:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
117:   /* ASCII viewer */
118:   if (isascii) {
119:     PetscMPIInt   color,rank,size;
120:     PetscInt64    loc[7],gsum[6],gmax[6],gmin[6],totbenign;
121:     PetscScalar   interface_size;
122:     PetscReal     ratio1=0.,ratio2=0.;
123:     Vec           counter;

125:     if (!pc->setupcalled) {
126:       PetscViewerASCIIPrintf(viewer,"  Partial information available: preconditioner has not been setup yet\n");
127:     }
128:     PetscViewerASCIIPrintf(viewer,"  Use verbose output: %D\n",pcbddc->dbg_flag);
129:     PetscViewerASCIIPrintf(viewer,"  Use user-defined CSR: %d\n",!!pcbddc->mat_graph->nvtxs_csr);
130:     PetscViewerASCIIPrintf(viewer,"  Use local mat graph: %d\n",pcbddc->use_local_adj && !pcbddc->mat_graph->nvtxs_csr);
131:     if (pcbddc->mat_graph->twodim) {
132:       PetscViewerASCIIPrintf(viewer,"  Connectivity graph topological dimension: 2\n");
133:     } else {
134:       PetscViewerASCIIPrintf(viewer,"  Connectivity graph topological dimension: 3\n");
135:     }
136:     if (pcbddc->graphmaxcount != PETSC_MAX_INT) {
137:       PetscViewerASCIIPrintf(viewer,"  Graph max count: %D\n",pcbddc->graphmaxcount);
138:     }
139:     PetscViewerASCIIPrintf(viewer,"  Use vertices: %d (vertex size %D)\n",pcbddc->use_vertices,pcbddc->vertex_size);
140:     PetscViewerASCIIPrintf(viewer,"  Use edges: %d\n",pcbddc->use_edges);
141:     PetscViewerASCIIPrintf(viewer,"  Use faces: %d\n",pcbddc->use_faces);
142:     PetscViewerASCIIPrintf(viewer,"  Use true near null space: %d\n",pcbddc->use_nnsp_true);
143:     PetscViewerASCIIPrintf(viewer,"  Use QR for single constraints on cc: %d\n",pcbddc->use_qr_single);
144:     PetscViewerASCIIPrintf(viewer,"  Use change of basis on local edge nodes: %d\n",pcbddc->use_change_of_basis);
145:     PetscViewerASCIIPrintf(viewer,"  Use change of basis on local face nodes: %d\n",pcbddc->use_change_on_faces);
146:     PetscViewerASCIIPrintf(viewer,"  User defined change of basis matrix: %d\n",!!pcbddc->user_ChangeOfBasisMatrix);
147:     PetscViewerASCIIPrintf(viewer,"  Has change of basis matrix: %d\n",!!pcbddc->ChangeOfBasisMatrix);
148:     PetscViewerASCIIPrintf(viewer,"  Eliminate dirichlet boundary dofs: %d\n",pcbddc->eliminate_dirdofs);
149:     PetscViewerASCIIPrintf(viewer,"  Switch on static condensation ops around the interface preconditioner: %d\n",pcbddc->switch_static);
150:     PetscViewerASCIIPrintf(viewer,"  Use exact dirichlet trick: %d\n",pcbddc->use_exact_dirichlet_trick);
151:     PetscViewerASCIIPrintf(viewer,"  Multilevel max levels: %D\n",pcbddc->max_levels);
152:     PetscViewerASCIIPrintf(viewer,"  Multilevel coarsening ratio: %D\n",pcbddc->coarsening_ratio);
153:     PetscViewerASCIIPrintf(viewer,"  Use estimated eigs for coarse problem: %d\n",pcbddc->use_coarse_estimates);
154:     PetscViewerASCIIPrintf(viewer,"  Use deluxe scaling: %d\n",pcbddc->use_deluxe_scaling);
155:     PetscViewerASCIIPrintf(viewer,"  Use deluxe zerorows: %d\n",pcbddc->deluxe_zerorows);
156:     PetscViewerASCIIPrintf(viewer,"  Use deluxe singlemat: %d\n",pcbddc->deluxe_singlemat);
157:     PetscViewerASCIIPrintf(viewer,"  Rebuild interface graph for Schur principal minors: %d\n",pcbddc->sub_schurs_rebuild);
158:     PetscViewerASCIIPrintf(viewer,"  Number of dofs' layers for the computation of principal minors: %D\n",pcbddc->sub_schurs_layers);
159:     PetscViewerASCIIPrintf(viewer,"  Use user CSR graph to compute successive layers: %d\n",pcbddc->sub_schurs_use_useradj);
160:     if (pcbddc->adaptive_threshold[1] != pcbddc->adaptive_threshold[0]) {
161:       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]);
162:     } else {
163:       PetscViewerASCIIPrintf(viewer,"  Adaptive constraint selection threshold (active %d, userdefined %d): %g\n",pcbddc->adaptive_selection,pcbddc->adaptive_userdefined,pcbddc->adaptive_threshold[0]);
164:     }
165:     PetscViewerASCIIPrintf(viewer,"  Min constraints / connected component: %D\n",pcbddc->adaptive_nmin);
166:     PetscViewerASCIIPrintf(viewer,"  Max constraints / connected component: %D\n",pcbddc->adaptive_nmax);
167:     PetscViewerASCIIPrintf(viewer,"  Invert exact Schur complement for adaptive selection: %d\n",pcbddc->sub_schurs_exact_schur);
168:     PetscViewerASCIIPrintf(viewer,"  Symmetric computation of primal basis functions: %d\n",pcbddc->symmetric_primal);
169:     PetscViewerASCIIPrintf(viewer,"  Num. Procs. to map coarse adjacency list: %D\n",pcbddc->coarse_adj_red);
170:     PetscViewerASCIIPrintf(viewer,"  Coarse eqs per proc (significant at the coarsest level): %D\n",pcbddc->coarse_eqs_per_proc);
171:     PetscViewerASCIIPrintf(viewer,"  Detect disconnected: %d\n",pcbddc->detect_disconnected);
172:     PetscViewerASCIIPrintf(viewer,"  Benign subspace trick: %d (change explicit %d)\n",pcbddc->benign_saddle_point,pcbddc->benign_change_explicit);
173:     PetscViewerASCIIPrintf(viewer,"  Benign subspace trick is active: %d\n",pcbddc->benign_have_null);
174:     PetscViewerASCIIPrintf(viewer,"  Algebraic computation of no-net-flux %d\n",pcbddc->compute_nonetflux);
175:     if (!pc->setupcalled) return(0);

177:     /* compute interface size */
178:     VecSet(pcis->vec1_B,1.0);
179:     MatCreateVecs(pc->pmat,&counter,0);
180:     VecSet(counter,0.0);
181:     VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,INSERT_VALUES,SCATTER_REVERSE);
182:     VecScatterEnd(pcis->global_to_B,pcis->vec1_B,counter,INSERT_VALUES,SCATTER_REVERSE);
183:     VecSum(counter,&interface_size);
184:     VecDestroy(&counter);

186:     /* compute some statistics on the domain decomposition */
187:     gsum[0] = 1;
188:     gsum[1] = gsum[2] = gsum[3] = gsum[4] = gsum[5] = 0;
189:     loc[0]  = !!pcis->n;
190:     loc[1]  = pcis->n - pcis->n_B;
191:     loc[2]  = pcis->n_B;
192:     loc[3]  = pcbddc->local_primal_size;
193:     loc[4]  = pcis->n;
194:     loc[5]  = pcbddc->n_local_subs > 0 ? pcbddc->n_local_subs : (pcis->n ? 1 : 0);
195:     loc[6]  = pcbddc->benign_n;
196:     MPI_Reduce(loc,gsum,6,MPIU_INT64,MPI_SUM,0,PetscObjectComm((PetscObject)pc));
197:     if (!loc[0]) loc[1] = loc[2] = loc[3] = loc[4] = loc[5] = -1;
198:     MPI_Reduce(loc,gmax,6,MPIU_INT64,MPI_MAX,0,PetscObjectComm((PetscObject)pc));
199:     if (!loc[0]) loc[1] = loc[2] = loc[3] = loc[4] = loc[5] = PETSC_MAX_INT;
200:     MPI_Reduce(loc,gmin,6,MPIU_INT64,MPI_MIN,0,PetscObjectComm((PetscObject)pc));
201:     MPI_Reduce(&loc[6],&totbenign,1,MPIU_INT64,MPI_SUM,0,PetscObjectComm((PetscObject)pc));
202:     if (pcbddc->coarse_size) {
203:       ratio1 = pc->pmat->rmap->N/(1.*pcbddc->coarse_size);
204:       ratio2 = PetscRealPart(interface_size)/pcbddc->coarse_size;
205:     }
206:     PetscViewerASCIIPrintf(viewer,"  ********************************** STATISTICS AT LEVEL %d **********************************\n",pcbddc->current_level);
207:     PetscViewerASCIIPrintf(viewer,"  Global dofs sizes: all %D interface %D coarse %D\n",pc->pmat->rmap->N,(PetscInt)PetscRealPart(interface_size),pcbddc->coarse_size);
208:     PetscViewerASCIIPrintf(viewer,"  Coarsening ratios: all/coarse %D interface/coarse %D\n",(PetscInt)ratio1,(PetscInt)ratio2);
209:     PetscViewerASCIIPrintf(viewer,"  Active processes : %D\n",(PetscInt)gsum[0]);
210:     PetscViewerASCIIPrintf(viewer,"  Total subdomains : %D\n",(PetscInt)gsum[5]);
211:     if (pcbddc->benign_have_null) {
212:       PetscViewerASCIIPrintf(viewer,"  Benign subs      : %D\n",(PetscInt)totbenign);
213:     }
214:     PetscViewerASCIIPrintf(viewer,"  Dofs type        :\tMIN\tMAX\tMEAN\n");
215:     PetscViewerASCIIPrintf(viewer,"  Interior  dofs   :\t%D\t%D\t%D\n",(PetscInt)gmin[1],(PetscInt)gmax[1],(PetscInt)(gsum[1]/gsum[0]));
216:     PetscViewerASCIIPrintf(viewer,"  Interface dofs   :\t%D\t%D\t%D\n",(PetscInt)gmin[2],(PetscInt)gmax[2],(PetscInt)(gsum[2]/gsum[0]));
217:     PetscViewerASCIIPrintf(viewer,"  Primal    dofs   :\t%D\t%D\t%D\n",(PetscInt)gmin[3],(PetscInt)gmax[3],(PetscInt)(gsum[3]/gsum[0]));
218:     PetscViewerASCIIPrintf(viewer,"  Local     dofs   :\t%D\t%D\t%D\n",(PetscInt)gmin[4],(PetscInt)gmax[4],(PetscInt)(gsum[4]/gsum[0]));
219:     PetscViewerASCIIPrintf(viewer,"  Local     subs   :\t%D\t%D\n"    ,(PetscInt)gmin[5],(PetscInt)gmax[5]);
220:     PetscViewerASCIIPrintf(viewer,"  ********************************** COARSE PROBLEM DETAILS *********************************\n");
221:     PetscViewerFlush(viewer);

223:     /* the coarse problem can be handled by a different communicator */
224:     if (pcbddc->coarse_ksp) color = 1;
225:     else color = 0;
226:     MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
227:     MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
228:     PetscSubcommCreate(PetscObjectComm((PetscObject)pc),&subcomm);
229:     PetscSubcommSetNumber(subcomm,PetscMin(size,2));
230:     PetscSubcommSetTypeGeneral(subcomm,color,rank);
231:     PetscViewerGetSubViewer(viewer,PetscSubcommChild(subcomm),&subviewer);
232:     if (color == 1) {
233:       KSPView(pcbddc->coarse_ksp,subviewer);
234:       PetscViewerFlush(subviewer);
235:     }
236:     PetscViewerRestoreSubViewer(viewer,PetscSubcommChild(subcomm),&subviewer);
237:     PetscSubcommDestroy(&subcomm);
238:     PetscViewerFlush(viewer);
239:   }
240:   return(0);
241: }

243: static PetscErrorCode PCBDDCSetDiscreteGradient_BDDC(PC pc, Mat G, PetscInt order, PetscInt field, PetscBool global, PetscBool conforming)
244: {
245:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

249:   PetscObjectReference((PetscObject)G);
250:   MatDestroy(&pcbddc->discretegradient);
251:   pcbddc->discretegradient = G;
252:   pcbddc->nedorder         = order > 0 ? order : -order;
253:   pcbddc->nedfield         = field;
254:   pcbddc->nedglobal        = global;
255:   pcbddc->conforming       = conforming;
256:   return(0);
257: }

259: /*@
260:  PCBDDCSetDiscreteGradient - Sets the discrete gradient

262:    Collective on PC

264:    Input Parameters:
265: +  pc         - the preconditioning context
266: .  G          - the discrete gradient matrix (should be in AIJ format)
267: .  order      - the order of the Nedelec space (1 for the lowest order)
268: .  field      - the field id of the Nedelec dofs (not used if the fields have not been specified)
269: .  global     - the type of global ordering for the rows of G
270: -  conforming - whether the mesh is conforming or not

272:    Level: advanced

274:    Notes: The discrete gradient matrix G is used to analyze the subdomain edges, and it should not contain any zero entry.
275:           For variable order spaces, the order should be set to zero.
276:           If global is true, the rows of G should be given in global ordering for the whole dofs;
277:           if false, the ordering should be global for the Nedelec field.
278:           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
279:           and geid the one for the Nedelec field.

281: .seealso: PCBDDC,PCBDDCSetDofsSplitting(),PCBDDCSetDofsSplittingLocal()
282: @*/
283: PetscErrorCode PCBDDCSetDiscreteGradient(PC pc, Mat G, PetscInt order, PetscInt field, PetscBool global, PetscBool conforming)
284: {

295:   PetscTryMethod(pc,"PCBDDCSetDiscreteGradient_C",(PC,Mat,PetscInt,PetscInt,PetscBool,PetscBool),(pc,G,order,field,global,conforming));
296:   return(0);
297: }

299: static PetscErrorCode PCBDDCSetDivergenceMat_BDDC(PC pc, Mat divudotp, PetscBool trans, IS vl2l)
300: {
301:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

305:   PetscObjectReference((PetscObject)divudotp);
306:   MatDestroy(&pcbddc->divudotp);
307:   pcbddc->divudotp = divudotp;
308:   pcbddc->divudotp_trans = trans;
309:   pcbddc->compute_nonetflux = PETSC_TRUE;
310:   if (vl2l) {
311:     PetscObjectReference((PetscObject)vl2l);
312:     ISDestroy(&pcbddc->divudotp_vl2l);
313:     pcbddc->divudotp_vl2l = vl2l;
314:   }
315:   return(0);
316: }

318: /*@
319:  PCBDDCSetDivergenceMat - Sets the linear operator representing \int_\Omega \div {\bf u} \cdot p dx

321:    Collective on PC

323:    Input Parameters:
324: +  pc - the preconditioning context
325: .  divudotp - the matrix (must be of type MATIS)
326: .  trans - if trans if false (resp. true), then pressures are in the test (trial) space and velocities are in the trial (test) space.
327: -  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

329:    Level: advanced

331:    Notes: This auxiliary matrix is used to compute quadrature weights representing the net-flux across subdomain boundaries
332:           If vl2l is NULL, the local ordering for velocities in divudotp should match that of the preconditioning matrix

334: .seealso: PCBDDC
335: @*/
336: PetscErrorCode PCBDDCSetDivergenceMat(PC pc, Mat divudotp, PetscBool trans, IS vl2l)
337: {
338:   PetscBool      ismatis;

347:   PetscObjectTypeCompare((PetscObject)divudotp,MATIS,&ismatis);
348:   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Divergence matrix needs to be of type MATIS");
349:   PetscTryMethod(pc,"PCBDDCSetDivergenceMat_C",(PC,Mat,PetscBool,IS),(pc,divudotp,trans,vl2l));
350:   return(0);
351: }

353: static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change, PetscBool interior)
354: {
355:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

359:   PetscObjectReference((PetscObject)change);
360:   MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);
361:   pcbddc->user_ChangeOfBasisMatrix = change;
362:   pcbddc->change_interior = interior;
363:   return(0);
364: }
365: /*@
366:  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs

368:    Collective on PC

370:    Input Parameters:
371: +  pc - the preconditioning context
372: .  change - the change of basis matrix
373: -  interior - whether or not the change of basis modifies interior dofs

375:    Level: intermediate

377:    Notes:

379: .seealso: PCBDDC
380: @*/
381: PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change, PetscBool interior)
382: {

389:   if (pc->mat) {
390:     PetscInt rows_c,cols_c,rows,cols;
391:     MatGetSize(pc->mat,&rows,&cols);
392:     MatGetSize(change,&rows_c,&cols_c);
393:     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);
394:     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);
395:     MatGetLocalSize(pc->mat,&rows,&cols);
396:     MatGetLocalSize(change,&rows_c,&cols_c);
397:     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);
398:     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);
399:   }
400:   PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat,PetscBool),(pc,change,interior));
401:   return(0);
402: }

404: static PetscErrorCode PCBDDCSetPrimalVerticesIS_BDDC(PC pc, IS PrimalVertices)
405: {
406:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
407:   PetscBool      isequal = PETSC_FALSE;

411:   PetscObjectReference((PetscObject)PrimalVertices);
412:   if (pcbddc->user_primal_vertices) {
413:     ISEqual(PrimalVertices,pcbddc->user_primal_vertices,&isequal);
414:   }
415:   ISDestroy(&pcbddc->user_primal_vertices);
416:   ISDestroy(&pcbddc->user_primal_vertices_local);
417:   pcbddc->user_primal_vertices = PrimalVertices;
418:   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
419:   return(0);
420: }

422: /*@
423:  PCBDDCSetPrimalVerticesIS - Set additional user defined primal vertices in PCBDDC

425:    Collective

427:    Input Parameters:
428: +  pc - the preconditioning context
429: -  PrimalVertices - index set of primal vertices in global numbering (can be empty)

431:    Level: intermediate

433:    Notes:
434:      Any process can list any global node

436: .seealso: PCBDDC
437: @*/
438: PetscErrorCode PCBDDCSetPrimalVerticesIS(PC pc, IS PrimalVertices)
439: {

446:   PetscTryMethod(pc,"PCBDDCSetPrimalVerticesIS_C",(PC,IS),(pc,PrimalVertices));
447:   return(0);
448: }

450: static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
451: {
452:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
453:   PetscBool      isequal = PETSC_FALSE;

457:   PetscObjectReference((PetscObject)PrimalVertices);
458:   if (pcbddc->user_primal_vertices_local) {
459:     ISEqual(PrimalVertices,pcbddc->user_primal_vertices_local,&isequal);
460:   }
461:   ISDestroy(&pcbddc->user_primal_vertices);
462:   ISDestroy(&pcbddc->user_primal_vertices_local);
463:   pcbddc->user_primal_vertices_local = PrimalVertices;
464:   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
465:   return(0);
466: }
467: /*@
468:  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC

470:    Collective

472:    Input Parameters:
473: +  pc - the preconditioning context
474: -  PrimalVertices - index set of primal vertices in local numbering (can be empty)

476:    Level: intermediate

478:    Notes:

480: .seealso: PCBDDC
481: @*/
482: PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
483: {

490:   PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));
491:   return(0);
492: }

494: static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
495: {
496:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

499:   pcbddc->coarsening_ratio = k;
500:   return(0);
501: }

503: /*@
504:  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel

506:    Logically collective on PC

508:    Input Parameters:
509: +  pc - the preconditioning context
510: -  k - coarsening ratio (H/h at the coarser level)

512:    Options Database Keys:
513: .    -pc_bddc_coarsening_ratio

515:    Level: intermediate

517:    Notes:
518:      Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level

520: .seealso: PCBDDC, PCBDDCSetLevels()
521: @*/
522: PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
523: {

529:   PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));
530:   return(0);
531: }

533: /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
534: static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
535: {
536:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

539:   pcbddc->use_exact_dirichlet_trick = flg;
540:   return(0);
541: }

543: PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
544: {

550:   PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));
551:   return(0);
552: }

554: static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
555: {
556:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

559:   pcbddc->current_level = level;
560:   return(0);
561: }

563: PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
564: {

570:   PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));
571:   return(0);
572: }

574: static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
575: {
576:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

579:   pcbddc->max_levels = levels;
580:   return(0);
581: }

583: /*@
584:  PCBDDCSetLevels - Sets the maximum number of levels for multilevel

586:    Logically collective on PC

588:    Input Parameters:
589: +  pc - the preconditioning context
590: -  levels - the maximum number of levels (max 9)

592:    Options Database Keys:
593: .    -pc_bddc_levels

595:    Level: intermediate

597:    Notes:
598:      Default value is 0, i.e. traditional one-level BDDC

600: .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
601: @*/
602: PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
603: {

609:   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
610:   PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));
611:   return(0);
612: }

614: static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
615: {
616:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
617:   PetscBool      isequal = PETSC_FALSE;

621:   PetscObjectReference((PetscObject)DirichletBoundaries);
622:   if (pcbddc->DirichletBoundaries) {
623:     ISEqual(DirichletBoundaries,pcbddc->DirichletBoundaries,&isequal);
624:   }
625:   /* last user setting takes precendence -> destroy any other customization */
626:   ISDestroy(&pcbddc->DirichletBoundariesLocal);
627:   ISDestroy(&pcbddc->DirichletBoundaries);
628:   pcbddc->DirichletBoundaries = DirichletBoundaries;
629:   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
630:   return(0);
631: }

633: /*@
634:  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.

636:    Collective

638:    Input Parameters:
639: +  pc - the preconditioning context
640: -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries

642:    Level: intermediate

644:    Notes:
645:      Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node

647: .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
648: @*/
649: PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
650: {

657:   PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));
658:   return(0);
659: }

661: static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
662: {
663:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
664:   PetscBool      isequal = PETSC_FALSE;

668:   PetscObjectReference((PetscObject)DirichletBoundaries);
669:   if (pcbddc->DirichletBoundariesLocal) {
670:     ISEqual(DirichletBoundaries,pcbddc->DirichletBoundariesLocal,&isequal);
671:   }
672:   /* last user setting takes precendence -> destroy any other customization */
673:   ISDestroy(&pcbddc->DirichletBoundariesLocal);
674:   ISDestroy(&pcbddc->DirichletBoundaries);
675:   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
676:   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
677:   return(0);
678: }

680: /*@
681:  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.

683:    Collective

685:    Input Parameters:
686: +  pc - the preconditioning context
687: -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)

689:    Level: intermediate

691:    Notes:

693: .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
694: @*/
695: PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
696: {

703:   PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));
704:   return(0);
705: }

707: static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
708: {
709:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
710:   PetscBool      isequal = PETSC_FALSE;

714:   PetscObjectReference((PetscObject)NeumannBoundaries);
715:   if (pcbddc->NeumannBoundaries) {
716:     ISEqual(NeumannBoundaries,pcbddc->NeumannBoundaries,&isequal);
717:   }
718:   /* last user setting takes precendence -> destroy any other customization */
719:   ISDestroy(&pcbddc->NeumannBoundariesLocal);
720:   ISDestroy(&pcbddc->NeumannBoundaries);
721:   pcbddc->NeumannBoundaries = NeumannBoundaries;
722:   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
723:   return(0);
724: }

726: /*@
727:  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.

729:    Collective

731:    Input Parameters:
732: +  pc - the preconditioning context
733: -  NeumannBoundaries - parallel IS defining the Neumann boundaries

735:    Level: intermediate

737:    Notes:
738:      Any process can list any global node

740: .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
741: @*/
742: PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
743: {

750:   PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));
751:   return(0);
752: }

754: static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
755: {
756:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
757:   PetscBool      isequal = PETSC_FALSE;

761:   PetscObjectReference((PetscObject)NeumannBoundaries);
762:   if (pcbddc->NeumannBoundariesLocal) {
763:     ISEqual(NeumannBoundaries,pcbddc->NeumannBoundariesLocal,&isequal);
764:   }
765:   /* last user setting takes precendence -> destroy any other customization */
766:   ISDestroy(&pcbddc->NeumannBoundariesLocal);
767:   ISDestroy(&pcbddc->NeumannBoundaries);
768:   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
769:   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
770:   return(0);
771: }

773: /*@
774:  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.

776:    Collective

778:    Input Parameters:
779: +  pc - the preconditioning context
780: -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)

782:    Level: intermediate

784:    Notes:

786: .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
787: @*/
788: PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
789: {

796:   PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));
797:   return(0);
798: }

800: static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
801: {
802:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

805:   *DirichletBoundaries = pcbddc->DirichletBoundaries;
806:   return(0);
807: }

809: /*@
810:  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries

812:    Collective

814:    Input Parameters:
815: .  pc - the preconditioning context

817:    Output Parameters:
818: .  DirichletBoundaries - index set defining the Dirichlet boundaries

820:    Level: intermediate

822:    Notes:
823:      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries

825: .seealso: PCBDDC
826: @*/
827: PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
828: {

833:   PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));
834:   return(0);
835: }

837: static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
838: {
839:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

842:   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
843:   return(0);
844: }

846: /*@
847:  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)

849:    Collective

851:    Input Parameters:
852: .  pc - the preconditioning context

854:    Output Parameters:
855: .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries

857:    Level: intermediate

859:    Notes:
860:      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).
861:           In the latter case, the IS will be available after PCSetUp.

863: .seealso: PCBDDC
864: @*/
865: PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
866: {

871:   PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));
872:   return(0);
873: }

875: static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
876: {
877:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

880:   *NeumannBoundaries = pcbddc->NeumannBoundaries;
881:   return(0);
882: }

884: /*@
885:  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries

887:    Collective

889:    Input Parameters:
890: .  pc - the preconditioning context

892:    Output Parameters:
893: .  NeumannBoundaries - index set defining the Neumann boundaries

895:    Level: intermediate

897:    Notes:
898:      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries

900: .seealso: PCBDDC
901: @*/
902: PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
903: {

908:   PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));
909:   return(0);
910: }

912: static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
913: {
914:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

917:   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
918:   return(0);
919: }

921: /*@
922:  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)

924:    Collective

926:    Input Parameters:
927: .  pc - the preconditioning context

929:    Output Parameters:
930: .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries

932:    Level: intermediate

934:    Notes:
935:      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).
936:           In the latter case, the IS will be available after PCSetUp.

938: .seealso: PCBDDC
939: @*/
940: PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
941: {

946:   PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));
947:   return(0);
948: }

950: static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
951: {
952:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
953:   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
954:   PetscBool      same_data = PETSC_FALSE;

958:   if (!nvtxs) {
959:     if (copymode == PETSC_OWN_POINTER) {
960:       PetscFree(xadj);
961:       PetscFree(adjncy);
962:     }
963:     PCBDDCGraphResetCSR(mat_graph);
964:     return(0);
965:   }
966:   if (mat_graph->nvtxs == nvtxs && mat_graph->freecsr) { /* we own the data */
967:     if (mat_graph->xadj == xadj && mat_graph->adjncy == adjncy) same_data = PETSC_TRUE;
968:     if (!same_data && mat_graph->xadj[nvtxs] == xadj[nvtxs]) {
969:       PetscMemcmp(xadj,mat_graph->xadj,(nvtxs+1)*sizeof(PetscInt),&same_data);
970:       if (same_data) {
971:         PetscMemcmp(adjncy,mat_graph->adjncy,xadj[nvtxs]*sizeof(PetscInt),&same_data);
972:       }
973:     }
974:   }
975:   if (!same_data) {
976:     /* free old CSR */
977:     PCBDDCGraphResetCSR(mat_graph);
978:     /* get CSR into graph structure */
979:     if (copymode == PETSC_COPY_VALUES) {
980:       PetscMalloc1(nvtxs+1,&mat_graph->xadj);
981:       PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);
982:       PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));
983:       PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));
984:       mat_graph->freecsr = PETSC_TRUE;
985:     } else if (copymode == PETSC_OWN_POINTER) {
986:       mat_graph->xadj    = (PetscInt*)xadj;
987:       mat_graph->adjncy  = (PetscInt*)adjncy;
988:       mat_graph->freecsr = PETSC_TRUE;
989:     } else if (copymode == PETSC_USE_POINTER) {
990:       mat_graph->xadj    = (PetscInt*)xadj;
991:       mat_graph->adjncy  = (PetscInt*)adjncy;
992:       mat_graph->freecsr = PETSC_FALSE;
993:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %D",copymode);
994:     mat_graph->nvtxs_csr = nvtxs;
995:     pcbddc->recompute_topography = PETSC_TRUE;
996:   }
997:   return(0);
998: }

1000: /*@
1001:  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local degrees of freedom.

1003:    Not collective

1005:    Input Parameters:
1006: +  pc - the preconditioning context.
1007: .  nvtxs - number of local vertices of the graph (i.e., the number of local dofs).
1008: .  xadj, adjncy - the connectivity of the dofs in CSR format.
1009: -  copymode - supported modes are PETSC_COPY_VALUES, PETSC_USE_POINTER or PETSC_OWN_POINTER.

1011:    Level: intermediate

1013:    Notes: A dof is considered connected with all local dofs if xadj[dof+1]-xadj[dof] == 1 and adjncy[xadj[dof]] is negative.

1015: .seealso: PCBDDC,PetscCopyMode
1016: @*/
1017: PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
1018: {
1019:   void (*f)(void) = 0;

1024:   if (nvtxs) {
1027:   }
1028:   PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));
1029:   /* free arrays if PCBDDC is not the PC type */
1030:   PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);
1031:   if (!f && copymode == PETSC_OWN_POINTER) {
1032:     PetscFree(xadj);
1033:     PetscFree(adjncy);
1034:   }
1035:   return(0);
1036: }

1038: static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
1039: {
1040:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1041:   PetscInt       i;
1042:   PetscBool      isequal = PETSC_FALSE;

1046:   if (pcbddc->n_ISForDofsLocal == n_is) {
1047:     for (i=0;i<n_is;i++) {
1048:       PetscBool isequalt;
1049:       ISEqual(ISForDofs[i],pcbddc->ISForDofsLocal[i],&isequalt);
1050:       if (!isequalt) break;
1051:     }
1052:     if (i == n_is) isequal = PETSC_TRUE;
1053:   }
1054:   for (i=0;i<n_is;i++) {
1055:     PetscObjectReference((PetscObject)ISForDofs[i]);
1056:   }
1057:   /* Destroy ISes if they were already set */
1058:   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
1059:     ISDestroy(&pcbddc->ISForDofsLocal[i]);
1060:   }
1061:   PetscFree(pcbddc->ISForDofsLocal);
1062:   /* last user setting takes precendence -> destroy any other customization */
1063:   for (i=0;i<pcbddc->n_ISForDofs;i++) {
1064:     ISDestroy(&pcbddc->ISForDofs[i]);
1065:   }
1066:   PetscFree(pcbddc->ISForDofs);
1067:   pcbddc->n_ISForDofs = 0;
1068:   /* allocate space then set */
1069:   if (n_is) {
1070:     PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);
1071:   }
1072:   for (i=0;i<n_is;i++) {
1073:     pcbddc->ISForDofsLocal[i] = ISForDofs[i];
1074:   }
1075:   pcbddc->n_ISForDofsLocal = n_is;
1076:   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
1077:   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
1078:   return(0);
1079: }

1081: /*@
1082:  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix

1084:    Collective

1086:    Input Parameters:
1087: +  pc - the preconditioning context
1088: .  n_is - number of index sets defining the fields
1089: -  ISForDofs - array of IS describing the fields in local ordering

1091:    Level: intermediate

1093:    Notes:
1094:      n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field.

1096: .seealso: PCBDDC
1097: @*/
1098: PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
1099: {
1100:   PetscInt       i;

1106:   for (i=0;i<n_is;i++) {
1109:   }
1110:   PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));
1111:   return(0);
1112: }

1114: static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
1115: {
1116:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1117:   PetscInt       i;
1118:   PetscBool      isequal = PETSC_FALSE;

1122:   if (pcbddc->n_ISForDofs == n_is) {
1123:     for (i=0;i<n_is;i++) {
1124:       PetscBool isequalt;
1125:       ISEqual(ISForDofs[i],pcbddc->ISForDofs[i],&isequalt);
1126:       if (!isequalt) break;
1127:     }
1128:     if (i == n_is) isequal = PETSC_TRUE;
1129:   }
1130:   for (i=0;i<n_is;i++) {
1131:     PetscObjectReference((PetscObject)ISForDofs[i]);
1132:   }
1133:   /* Destroy ISes if they were already set */
1134:   for (i=0;i<pcbddc->n_ISForDofs;i++) {
1135:     ISDestroy(&pcbddc->ISForDofs[i]);
1136:   }
1137:   PetscFree(pcbddc->ISForDofs);
1138:   /* last user setting takes precendence -> destroy any other customization */
1139:   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
1140:     ISDestroy(&pcbddc->ISForDofsLocal[i]);
1141:   }
1142:   PetscFree(pcbddc->ISForDofsLocal);
1143:   pcbddc->n_ISForDofsLocal = 0;
1144:   /* allocate space then set */
1145:   if (n_is) {
1146:     PetscMalloc1(n_is,&pcbddc->ISForDofs);
1147:   }
1148:   for (i=0;i<n_is;i++) {
1149:     pcbddc->ISForDofs[i] = ISForDofs[i];
1150:   }
1151:   pcbddc->n_ISForDofs = n_is;
1152:   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
1153:   if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
1154:   return(0);
1155: }

1157: /*@
1158:  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix

1160:    Collective

1162:    Input Parameters:
1163: +  pc - the preconditioning context
1164: .  n_is - number of index sets defining the fields
1165: -  ISForDofs - array of IS describing the fields in global ordering

1167:    Level: intermediate

1169:    Notes:
1170:      Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field.

1172: .seealso: PCBDDC
1173: @*/
1174: PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
1175: {
1176:   PetscInt       i;

1182:   for (i=0;i<n_is;i++) {
1185:   }
1186:   PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));
1187:   return(0);
1188: }

1190: /*
1191:    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
1192:                      guess if a transformation of basis approach has been selected.

1194:    Input Parameter:
1195: +  pc - the preconditioner contex

1197:    Application Interface Routine: PCPreSolve()

1199:    Notes:
1200:      The interface routine PCPreSolve() is not usually called directly by
1201:    the user, but instead is called by KSPSolve().
1202: */
1203: static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1204: {
1206:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1207:   PC_IS          *pcis = (PC_IS*)(pc->data);
1208:   Vec            used_vec;
1209:   PetscBool      save_rhs = PETSC_TRUE, benign_correction_computed;

1212:   /* if we are working with CG, one dirichlet solve can be avoided during Krylov iterations */
1213:   if (ksp) {
1214:     PetscBool iscg, isgroppcg, ispipecg, ispipecgrr;
1215:     PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);
1216:     PetscObjectTypeCompare((PetscObject)ksp,KSPGROPPCG,&isgroppcg);
1217:     PetscObjectTypeCompare((PetscObject)ksp,KSPPIPECG,&ispipecg);
1218:     PetscObjectTypeCompare((PetscObject)ksp,KSPPIPECGRR,&ispipecgrr);
1219:     if (pcbddc->benign_apply_coarse_only || pcbddc->switch_static || (!iscg && !isgroppcg && !ispipecg && !ispipecgrr)) {
1220:       PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);
1221:     }
1222:   }
1223:   if (pcbddc->benign_apply_coarse_only || pcbddc->switch_static) {
1224:     PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);
1225:   }

1227:   /* Creates parallel work vectors used in presolve */
1228:   if (!pcbddc->original_rhs) {
1229:     VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);
1230:   }
1231:   if (!pcbddc->temp_solution) {
1232:     VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);
1233:   }

1235:   pcbddc->temp_solution_used = PETSC_FALSE;
1236:   if (x) {
1237:     PetscObjectReference((PetscObject)x);
1238:     used_vec = x;
1239:   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
1240:     PetscObjectReference((PetscObject)pcbddc->temp_solution);
1241:     used_vec = pcbddc->temp_solution;
1242:     VecSet(used_vec,0.0);
1243:     pcbddc->temp_solution_used = PETSC_TRUE;
1244:     VecCopy(rhs,pcbddc->original_rhs);
1245:     save_rhs = PETSC_FALSE;
1246:     pcbddc->eliminate_dirdofs = PETSC_TRUE;
1247:   }

1249:   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
1250:   if (ksp) {
1251:     /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
1252:     KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);
1253:     if (!pcbddc->ksp_guess_nonzero) {
1254:       VecSet(used_vec,0.0);
1255:     }
1256:   }

1258:   pcbddc->rhs_change = PETSC_FALSE;
1259:   /* Take into account zeroed rows -> change rhs and store solution removed */
1260:   if (rhs && pcbddc->eliminate_dirdofs) {
1261:     IS dirIS = NULL;

1263:     /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
1264:     PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);
1265:     if (dirIS) {
1266:       Mat_IS            *matis = (Mat_IS*)pc->pmat->data;
1267:       PetscInt          dirsize,i,*is_indices;
1268:       PetscScalar       *array_x;
1269:       const PetscScalar *array_diagonal;

1271:       MatGetDiagonal(pc->pmat,pcis->vec1_global);
1272:       VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);
1273:       VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);
1274:       VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);
1275:       VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
1276:       VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
1277:       ISGetLocalSize(dirIS,&dirsize);
1278:       VecGetArray(pcis->vec1_N,&array_x);
1279:       VecGetArrayRead(pcis->vec2_N,&array_diagonal);
1280:       ISGetIndices(dirIS,(const PetscInt**)&is_indices);
1281:       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
1282:       ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);
1283:       VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);
1284:       VecRestoreArray(pcis->vec1_N,&array_x);
1285:       VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);
1286:       VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);
1287:       pcbddc->rhs_change = PETSC_TRUE;
1288:       ISDestroy(&dirIS);
1289:     }
1290:   }

1292:   /* remove the computed solution or the initial guess from the rhs */
1293:   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
1294:     /* save the original rhs */
1295:     if (save_rhs) {
1296:       VecSwap(rhs,pcbddc->original_rhs);
1297:       save_rhs = PETSC_FALSE;
1298:     }
1299:     pcbddc->rhs_change = PETSC_TRUE;
1300:     VecScale(used_vec,-1.0);
1301:     MatMultAdd(pc->mat,used_vec,pcbddc->original_rhs,rhs);
1302:     VecScale(used_vec,-1.0);
1303:     VecCopy(used_vec,pcbddc->temp_solution);
1304:     pcbddc->temp_solution_used = PETSC_TRUE;
1305:     if (ksp) {
1306:       KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);
1307:     }
1308:   }
1309:   VecDestroy(&used_vec);

1311:   /* compute initial vector in benign space if needed
1312:      and remove non-benign solution from the rhs */
1313:   benign_correction_computed = PETSC_FALSE;
1314:   if (rhs && pcbddc->benign_compute_correction && (pcbddc->benign_have_null || pcbddc->benign_apply_coarse_only)) {
1315:     /* compute u^*_h using ideas similar to those in Xuemin Tu's PhD thesis (see Section 4.8.1)
1316:        Recursively apply BDDC in the multilevel case */
1317:     if (!pcbddc->benign_vec) {
1318:       VecDuplicate(rhs,&pcbddc->benign_vec);
1319:     }
1320:     /* keep applying coarse solver unless we no longer have benign subdomains */
1321:     pcbddc->benign_apply_coarse_only = pcbddc->benign_have_null ? PETSC_TRUE : PETSC_FALSE;
1322:     if (!pcbddc->benign_skip_correction) {
1323:       PCApply_BDDC(pc,rhs,pcbddc->benign_vec);
1324:       benign_correction_computed = PETSC_TRUE;
1325:       if (pcbddc->temp_solution_used) {
1326:         VecAXPY(pcbddc->temp_solution,1.0,pcbddc->benign_vec);
1327:       }
1328:       VecScale(pcbddc->benign_vec,-1.0);
1329:       /* store the original rhs if not done earlier */
1330:       if (save_rhs) {
1331:         VecSwap(rhs,pcbddc->original_rhs);
1332:       }
1333:       if (pcbddc->rhs_change) {
1334:         MatMultAdd(pc->mat,pcbddc->benign_vec,rhs,rhs);
1335:       } else {
1336:         MatMultAdd(pc->mat,pcbddc->benign_vec,pcbddc->original_rhs,rhs);
1337:       }
1338:       pcbddc->rhs_change = PETSC_TRUE;
1339:     }
1340:     pcbddc->benign_apply_coarse_only = PETSC_FALSE;
1341:   }

1343:   /* dbg output */
1344:   if (pcbddc->dbg_flag && benign_correction_computed) {
1345:     Vec v;

1347:     VecDuplicate(pcis->vec1_global,&v);
1348:     if (pcbddc->ChangeOfBasisMatrix) {
1349:       MatMultTranspose(pcbddc->ChangeOfBasisMatrix,rhs,v);
1350:     } else {
1351:       VecCopy(rhs,v);
1352:     }
1353:     PCBDDCBenignGetOrSetP0(pc,v,PETSC_TRUE);
1354:     PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"LEVEL %D: is the correction benign?\n",pcbddc->current_level);
1355:     PetscScalarView(pcbddc->benign_n,pcbddc->benign_p0,pcbddc->dbg_viewer);
1356:     PetscViewerFlush(pcbddc->dbg_viewer);
1357:     VecDestroy(&v);
1358:   }

1360:   /* set initial guess if using PCG */
1361:   pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
1362:   if (x && pcbddc->use_exact_dirichlet_trick) {
1363:     VecSet(x,0.0);
1364:     if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) {
1365:       if (benign_correction_computed) { /* we have already saved the changed rhs */
1366:         VecLockPop(pcis->vec1_global);
1367:       } else {
1368:         MatMultTranspose(pcbddc->ChangeOfBasisMatrix,rhs,pcis->vec1_global);
1369:       }
1370:       VecScatterBegin(pcis->global_to_D,pcis->vec1_global,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1371:       VecScatterEnd(pcis->global_to_D,pcis->vec1_global,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1372:     } else {
1373:       VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1374:       VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1375:     }
1376:     KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1377:     if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) {
1378:       VecSet(pcis->vec1_global,0.);
1379:       VecScatterBegin(pcis->global_to_D,pcis->vec2_D,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);
1380:       VecScatterEnd(pcis->global_to_D,pcis->vec2_D,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);
1381:       MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_global,x);
1382:     } else {
1383:       VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);
1384:       VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);
1385:     }
1386:     if (ksp) {
1387:       KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
1388:     }
1389:     pcbddc->exact_dirichlet_trick_app = PETSC_TRUE;
1390:   } else if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior && benign_correction_computed && pcbddc->use_exact_dirichlet_trick) {
1391:     VecLockPop(pcis->vec1_global);
1392:   }
1393:   return(0);
1394: }

1396: /*
1397:    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1398:                      approach has been selected. Also, restores rhs to its original state.

1400:    Input Parameter:
1401: +  pc - the preconditioner contex

1403:    Application Interface Routine: PCPostSolve()

1405:    Notes:
1406:      The interface routine PCPostSolve() is not usually called directly by
1407:      the user, but instead is called by KSPSolve().
1408: */
1409: static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1410: {
1412:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

1415:   /* add solution removed in presolve */
1416:   if (x && pcbddc->rhs_change) {
1417:     if (pcbddc->temp_solution_used) {
1418:       VecAXPY(x,1.0,pcbddc->temp_solution);
1419:     } else if (pcbddc->benign_compute_correction && pcbddc->benign_vec) {
1420:       VecAXPY(x,-1.0,pcbddc->benign_vec);
1421:     }
1422:     /* restore to original state (not for FETI-DP) */
1423:     if (ksp) pcbddc->temp_solution_used = PETSC_FALSE;
1424:   }

1426:   /* restore rhs to its original state (not needed for FETI-DP) */
1427:   if (rhs && pcbddc->rhs_change) {
1428:     VecSwap(rhs,pcbddc->original_rhs);
1429:     pcbddc->rhs_change = PETSC_FALSE;
1430:   }
1431:   /* restore ksp guess state */
1432:   if (ksp) {
1433:     KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);
1434:     /* reset flag for exact dirichlet trick */
1435:     pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
1436:   }
1437:   return(0);
1438: }

1440: /*
1441:    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1442:                   by setting data structures and options.

1444:    Input Parameter:
1445: +  pc - the preconditioner context

1447:    Application Interface Routine: PCSetUp()

1449:    Notes:
1450:      The interface routine PCSetUp() is not usually called directly by
1451:      the user, but instead is called by PCApply() if necessary.
1452: */
1453: PetscErrorCode PCSetUp_BDDC(PC pc)
1454: {
1455:   PC_BDDC*        pcbddc = (PC_BDDC*)pc->data;
1456:   PCBDDCSubSchurs sub_schurs;
1457:   Mat_IS*         matis;
1458:   MatNullSpace    nearnullspace;
1459:   Mat             lA;
1460:   IS              lP,zerodiag = NULL;
1461:   PetscInt        nrows,ncols;
1462:   PetscMPIInt     size;
1463:   PetscBool       computesubschurs;
1464:   PetscBool       computeconstraintsmatrix;
1465:   PetscBool       new_nearnullspace_provided,ismatis;
1466:   PetscErrorCode  ierr;

1469:   PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);
1470:   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
1471:   MatGetSize(pc->pmat,&nrows,&ncols);
1472:   if (nrows != ncols) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
1473:   MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);

1475:   matis = (Mat_IS*)pc->pmat->data;
1476:   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1477:   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1478:      Also, BDDC builds its own KSP for the Dirichlet problem */
1479:   if (!pc->setupcalled || pc->flag == DIFFERENT_NONZERO_PATTERN) pcbddc->recompute_topography = PETSC_TRUE;
1480:   if (pcbddc->recompute_topography) {
1481:     pcbddc->graphanalyzed    = PETSC_FALSE;
1482:     computeconstraintsmatrix = PETSC_TRUE;
1483:   } else {
1484:     computeconstraintsmatrix = PETSC_FALSE;
1485:   }

1487:   /* check parameters' compatibility */
1488:   if (!pcbddc->use_deluxe_scaling) pcbddc->deluxe_zerorows = PETSC_FALSE;
1489:   pcbddc->adaptive_selection   = (PetscBool)(pcbddc->adaptive_threshold[0] != 0.0 || pcbddc->adaptive_threshold[1] != 0.0);
1490:   pcbddc->use_deluxe_scaling   = (PetscBool)(pcbddc->use_deluxe_scaling && size > 1);
1491:   pcbddc->adaptive_selection   = (PetscBool)(pcbddc->adaptive_selection && size > 1);
1492:   pcbddc->adaptive_userdefined = (PetscBool)(pcbddc->adaptive_selection && pcbddc->adaptive_userdefined);
1493:   if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;

1495:   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
1496:   if (pcbddc->switch_static) {
1497:     PetscBool ismatis;
1498:     PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&ismatis);
1499:     if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"When the static switch is one, the iteration matrix should be of type MATIS");
1500:   }

1502:   /* activate all connected components if the netflux has been requested */
1503:   if (pcbddc->compute_nonetflux) {
1504:     pcbddc->use_vertices = PETSC_TRUE;
1505:     pcbddc->use_edges = PETSC_TRUE;
1506:     pcbddc->use_faces = PETSC_TRUE;
1507:   }

1509:   /* Get stdout for dbg */
1510:   if (pcbddc->dbg_flag) {
1511:     if (!pcbddc->dbg_viewer) {
1512:       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1513:       PetscViewerASCIIPushSynchronized(pcbddc->dbg_viewer);
1514:     }
1515:     PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);
1516:   }

1518:   /* process topology information */
1519:   if (pcbddc->recompute_topography) {
1520:     PCBDDCComputeLocalTopologyInfo(pc);
1521:     if (pcbddc->discretegradient) {
1522:       PCBDDCNedelecSupport(pc);
1523:     }
1524:   }

1526:   /* change basis if requested by the user */
1527:   if (pcbddc->user_ChangeOfBasisMatrix) {
1528:     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1529:     pcbddc->use_change_of_basis = PETSC_FALSE;
1530:     PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);
1531:   } else {
1532:     MatDestroy(&pcbddc->local_mat);
1533:     PetscObjectReference((PetscObject)matis->A);
1534:     pcbddc->local_mat = matis->A;
1535:   }

1537:   /*
1538:      Compute change of basis on local pressures (aka zerodiag dofs) with the benign trick
1539:      This should come earlier then PCISSetUp for extracting the correct subdomain matrices
1540:   */
1541:   PCBDDCBenignShellMat(pc,PETSC_TRUE);
1542:   if (pcbddc->benign_saddle_point) {
1543:     PC_IS* pcis = (PC_IS*)pc->data;

1545:     if (pcbddc->user_ChangeOfBasisMatrix || pcbddc->use_change_of_basis || !computesubschurs) pcbddc->benign_change_explicit = PETSC_TRUE;
1546:     /* detect local saddle point and change the basis in pcbddc->local_mat (TODO: reuse case) */
1547:     PCBDDCBenignDetectSaddlePoint(pc,&zerodiag);
1548:     /* pop B0 mat from local mat */
1549:     PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);
1550:     /* give pcis a hint to not reuse submatrices during PCISCreate */
1551:     if (pc->flag == SAME_NONZERO_PATTERN && pcis->reusesubmatrices == PETSC_TRUE) {
1552:       if (pcbddc->benign_n && (pcbddc->benign_change_explicit || pcbddc->dbg_flag)) {
1553:         pcis->reusesubmatrices = PETSC_FALSE;
1554:       } else {
1555:         pcis->reusesubmatrices = PETSC_TRUE;
1556:       }
1557:     } else {
1558:       pcis->reusesubmatrices = PETSC_FALSE;
1559:     }
1560:   }

1562:   /* propagate relevant information */
1563:   if (matis->A->symmetric_set) {
1564:     MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);
1565:   }
1566:   if (matis->A->spd_set) {
1567:     MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);
1568:   }

1570:   /* Set up all the "iterative substructuring" common block without computing solvers */
1571:   {
1572:     Mat temp_mat;

1574:     temp_mat = matis->A;
1575:     matis->A = pcbddc->local_mat;
1576:     PCISSetUp(pc,PETSC_FALSE);
1577:     pcbddc->local_mat = matis->A;
1578:     matis->A = temp_mat;
1579:   }

1581:   /* Analyze interface */
1582:   if (!pcbddc->graphanalyzed) {
1583:     PCBDDCAnalyzeInterface(pc);
1584:     computeconstraintsmatrix = PETSC_TRUE;
1585:     if (pcbddc->adaptive_selection && !pcbddc->use_deluxe_scaling && !pcbddc->mat_graph->twodim) {
1586:       SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute the adaptive primal space for a problem with 3D edges without deluxe scaling");
1587:     }
1588:     if (pcbddc->compute_nonetflux) {
1589:       MatNullSpace nnfnnsp;

1591:       if (!pcbddc->divudotp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Missing divudotp operator");
1592:       PCBDDCComputeNoNetFlux(pc->pmat,pcbddc->divudotp,pcbddc->divudotp_trans,pcbddc->divudotp_vl2l,pcbddc->mat_graph,&nnfnnsp);
1593:       /* TODO what if a nearnullspace is already attached? */
1594:       if (nnfnnsp) {
1595:         MatSetNearNullSpace(pc->pmat,nnfnnsp);
1596:         MatNullSpaceDestroy(&nnfnnsp);
1597:       }
1598:     }
1599:   }

1601:   /* check existence of a divergence free extension, i.e.
1602:      b(v_I,p_0) = 0 for all v_I (raise error if not).
1603:      Also, check that PCBDDCBenignGetOrSetP0 works */
1604:   if (pcbddc->benign_saddle_point && pcbddc->dbg_flag > 1) {
1605:     PCBDDCBenignCheck(pc,zerodiag);
1606:   }
1607:   ISDestroy(&zerodiag);

1609:   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1610:   if (computesubschurs && pcbddc->recompute_topography) {
1611:     PCBDDCInitSubSchurs(pc);
1612:   }
1613:   /* SetUp Scaling operator (scaling matrices could be needed in SubSchursSetUp)*/
1614:   if (!pcbddc->use_deluxe_scaling) {
1615:     PCBDDCScalingSetUp(pc);
1616:   }

1618:   /* finish setup solvers and do adaptive selection of constraints */
1619:   sub_schurs = pcbddc->sub_schurs;
1620:   if (sub_schurs && sub_schurs->schur_explicit) {
1621:     if (computesubschurs) {
1622:       PCBDDCSetUpSubSchurs(pc);
1623:     }
1624:     PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);
1625:   } else {
1626:     PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);
1627:     if (computesubschurs) {
1628:       PCBDDCSetUpSubSchurs(pc);
1629:     }
1630:   }
1631:   if (pcbddc->adaptive_selection) {
1632:     PCBDDCAdaptiveSelection(pc);
1633:     computeconstraintsmatrix = PETSC_TRUE;
1634:   }

1636:   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1637:   new_nearnullspace_provided = PETSC_FALSE;
1638:   MatGetNearNullSpace(pc->pmat,&nearnullspace);
1639:   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1640:     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1641:       new_nearnullspace_provided = PETSC_TRUE;
1642:     } else {
1643:       /* determine if the two nullspaces are different (should be lightweight) */
1644:       if (nearnullspace != pcbddc->onearnullspace) {
1645:         new_nearnullspace_provided = PETSC_TRUE;
1646:       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1647:         PetscInt         i;
1648:         const Vec        *nearnullvecs;
1649:         PetscObjectState state;
1650:         PetscInt         nnsp_size;
1651:         MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);
1652:         for (i=0;i<nnsp_size;i++) {
1653:           PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);
1654:           if (pcbddc->onearnullvecs_state[i] != state) {
1655:             new_nearnullspace_provided = PETSC_TRUE;
1656:             break;
1657:           }
1658:         }
1659:       }
1660:     }
1661:   } else {
1662:     if (!nearnullspace) { /* both nearnullspaces are null */
1663:       new_nearnullspace_provided = PETSC_FALSE;
1664:     } else { /* nearnullspace attached later */
1665:       new_nearnullspace_provided = PETSC_TRUE;
1666:     }
1667:   }

1669:   /* Setup constraints and related work vectors */
1670:   /* reset primal space flags */
1671:   pcbddc->new_primal_space = PETSC_FALSE;
1672:   pcbddc->new_primal_space_local = PETSC_FALSE;
1673:   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1674:     /* It also sets the primal space flags */
1675:     PCBDDCConstraintsSetUp(pc);
1676:   }
1677:   /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1678:   PCBDDCSetUpLocalWorkVectors(pc);

1680:   if (pcbddc->use_change_of_basis) {
1681:     PC_IS *pcis = (PC_IS*)(pc->data);

1683:     PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);
1684:     if (pcbddc->benign_change) {
1685:       MatDestroy(&pcbddc->benign_B0);
1686:       /* pop B0 from pcbddc->local_mat */
1687:       PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);
1688:     }
1689:     /* get submatrices */
1690:     MatDestroy(&pcis->A_IB);
1691:     MatDestroy(&pcis->A_BI);
1692:     MatDestroy(&pcis->A_BB);
1693:     MatCreateSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);
1694:     MatCreateSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);
1695:     MatCreateSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);
1696:     /* set flag in pcis to not reuse submatrices during PCISCreate */
1697:     pcis->reusesubmatrices = PETSC_FALSE;
1698:   } else if (!pcbddc->user_ChangeOfBasisMatrix && !pcbddc->benign_change) {
1699:     MatDestroy(&pcbddc->local_mat);
1700:     PetscObjectReference((PetscObject)matis->A);
1701:     pcbddc->local_mat = matis->A;
1702:   }

1704:   /* interface pressure block row for B_C */
1705:   PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_lP" ,(PetscObject*)&lP);
1706:   PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_lA" ,(PetscObject*)&lA);
1707:   if (lA && lP) {
1708:     PC_IS*    pcis = (PC_IS*)pc->data;
1709:     Mat       B_BI,B_BB,Bt_BI,Bt_BB;
1710:     PetscBool issym;
1711:     MatIsSymmetric(lA,PETSC_SMALL,&issym);
1712:     if (issym) {
1713:       MatCreateSubMatrix(lA,lP,pcis->is_I_local,MAT_INITIAL_MATRIX,&B_BI);
1714:       MatCreateSubMatrix(lA,lP,pcis->is_B_local,MAT_INITIAL_MATRIX,&B_BB);
1715:       MatCreateTranspose(B_BI,&Bt_BI);
1716:       MatCreateTranspose(B_BB,&Bt_BB);
1717:     } else {
1718:       MatCreateSubMatrix(lA,lP,pcis->is_I_local,MAT_INITIAL_MATRIX,&B_BI);
1719:       MatCreateSubMatrix(lA,lP,pcis->is_B_local,MAT_INITIAL_MATRIX,&B_BB);
1720:       MatCreateSubMatrix(lA,pcis->is_I_local,lP,MAT_INITIAL_MATRIX,&Bt_BI);
1721:       MatCreateSubMatrix(lA,pcis->is_B_local,lP,MAT_INITIAL_MATRIX,&Bt_BB);
1722:     }
1723:     PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_B_BI",(PetscObject)B_BI);
1724:     PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_B_BB",(PetscObject)B_BB);
1725:     PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_Bt_BI",(PetscObject)Bt_BI);
1726:     PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_Bt_BB",(PetscObject)Bt_BB);
1727:     MatDestroy(&B_BI);
1728:     MatDestroy(&B_BB);
1729:     MatDestroy(&Bt_BI);
1730:     MatDestroy(&Bt_BB);
1731:   }

1733:   /* SetUp coarse and local Neumann solvers */
1734:   PCBDDCSetUpSolvers(pc);
1735:   /* SetUp Scaling operator */
1736:   if (pcbddc->use_deluxe_scaling) {
1737:     PCBDDCScalingSetUp(pc);
1738:   }

1740:   /* mark topography as done */
1741:   pcbddc->recompute_topography = PETSC_FALSE;

1743:   /* wrap pcis->A_IB and pcis->A_BI if we did not change explicitly the variables on the pressures */
1744:   PCBDDCBenignShellMat(pc,PETSC_FALSE);

1746:   if (pcbddc->dbg_flag) {
1747:     PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);
1748:   }
1749:   return(0);
1750: }

1752: /*
1753:    PCApply_BDDC - Applies the BDDC operator to a vector.

1755:    Input Parameters:
1756: +  pc - the preconditioner context
1757: -  r - input vector (global)

1759:    Output Parameter:
1760: .  z - output vector (global)

1762:    Application Interface Routine: PCApply()
1763:  */
1764: PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1765: {
1766:   PC_IS             *pcis = (PC_IS*)(pc->data);
1767:   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1768:   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
1769:   PetscErrorCode    ierr;
1770:   const PetscScalar one = 1.0;
1771:   const PetscScalar m_one = -1.0;
1772:   const PetscScalar zero = 0.0;

1774: /* This code is similar to that provided in nn.c for PCNN
1775:    NN interface preconditioner changed to BDDC
1776:    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */

1779:   PetscCitationsRegister(citation,&cited);
1780:   if (pcbddc->ChangeOfBasisMatrix) {
1781:     Vec swap;

1783:     MatMultTranspose(pcbddc->ChangeOfBasisMatrix,r,pcbddc->work_change);
1784:     swap = pcbddc->work_change;
1785:     pcbddc->work_change = r;
1786:     r = swap;
1787:     /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */
1788:     if (pcbddc->benign_apply_coarse_only && pcbddc->use_exact_dirichlet_trick && pcbddc->change_interior) {
1789:       VecCopy(r,pcis->vec1_global);
1790:       VecLockPush(pcis->vec1_global);
1791:     }
1792:   }
1793:   if (pcbddc->benign_have_null) { /* get p0 from r */
1794:     PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);
1795:   }
1796:   if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
1797:     VecCopy(r,z);
1798:     /* First Dirichlet solve */
1799:     VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1800:     VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1801:     /*
1802:       Assembling right hand side for BDDC operator
1803:       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1804:       - pcis->vec1_B the interface part of the global vector z
1805:     */
1806:     if (n_D) {
1807:       KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1808:       VecScale(pcis->vec2_D,m_one);
1809:       if (pcbddc->switch_static) {
1810:         Mat_IS *matis = (Mat_IS*)(pc->mat->data);

1812:         VecSet(pcis->vec1_N,0.);
1813:         VecScatterBegin(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1814:         VecScatterEnd(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1815:         if (!pcbddc->switch_static_change) {
1816:           MatMult(matis->A,pcis->vec1_N,pcis->vec2_N);
1817:         } else {
1818:           MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
1819:           MatMult(matis->A,pcis->vec2_N,pcis->vec1_N);
1820:           MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
1821:         }
1822:         VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);
1823:         VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);
1824:         VecScatterBegin(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1825:         VecScatterEnd(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1826:       } else {
1827:         MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);
1828:       }
1829:     } else {
1830:       VecSet(pcis->vec1_B,zero);
1831:     }
1832:     VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1833:     VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1834:     PCBDDCScalingRestriction(pc,z,pcis->vec1_B);
1835:   } else {
1836:     if (!pcbddc->benign_apply_coarse_only) {
1837:       PCBDDCScalingRestriction(pc,r,pcis->vec1_B);
1838:     }
1839:   }

1841:   /* Apply interface preconditioner
1842:      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1843:   PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);

1845:   /* Apply transpose of partition of unity operator */
1846:   PCBDDCScalingExtension(pc,pcis->vec1_B,z);

1848:   /* Second Dirichlet solve and assembling of output */
1849:   VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1850:   VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1851:   if (n_B) {
1852:     if (pcbddc->switch_static) {
1853:       Mat_IS *matis = (Mat_IS*)(pc->mat->data);

1855:       VecScatterBegin(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1856:       VecScatterEnd(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1857:       VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1858:       VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1859:       if (!pcbddc->switch_static_change) {
1860:         MatMult(matis->A,pcis->vec1_N,pcis->vec2_N);
1861:       } else {
1862:         MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
1863:         MatMult(matis->A,pcis->vec2_N,pcis->vec1_N);
1864:         MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
1865:       }
1866:       VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);
1867:       VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);
1868:     } else {
1869:       MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);
1870:     }
1871:   } else if (pcbddc->switch_static) { /* n_B is zero */
1872:     Mat_IS *matis = (Mat_IS*)(pc->mat->data);

1874:     if (!pcbddc->switch_static_change) {
1875:       MatMult(matis->A,pcis->vec1_D,pcis->vec3_D);
1876:     } else {
1877:       MatMult(pcbddc->switch_static_change,pcis->vec1_D,pcis->vec1_N);
1878:       MatMult(matis->A,pcis->vec1_N,pcis->vec2_N);
1879:       MatMultTranspose(pcbddc->switch_static_change,pcis->vec2_N,pcis->vec3_D);
1880:     }
1881:   }
1882:   KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);

1884:   if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
1885:     if (pcbddc->switch_static) {
1886:       VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);
1887:     } else {
1888:       VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);
1889:     }
1890:     VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1891:     VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1892:   } else {
1893:     if (pcbddc->switch_static) {
1894:       VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);
1895:     } else {
1896:       VecScale(pcis->vec4_D,m_one);
1897:     }
1898:     VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
1899:     VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
1900:   }
1901:   if (pcbddc->benign_have_null) { /* set p0 (computed in PCBDDCApplyInterface) */
1902:     if (pcbddc->benign_apply_coarse_only) {
1903:       PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));
1904:     }
1905:     PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);
1906:   }

1908:   if (pcbddc->ChangeOfBasisMatrix) {
1909:     pcbddc->work_change = r;
1910:     VecCopy(z,pcbddc->work_change);
1911:     MatMult(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,z);
1912:   }
1913:   return(0);
1914: }

1916: /*
1917:    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.

1919:    Input Parameters:
1920: +  pc - the preconditioner context
1921: -  r - input vector (global)

1923:    Output Parameter:
1924: .  z - output vector (global)

1926:    Application Interface Routine: PCApplyTranspose()
1927:  */
1928: PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
1929: {
1930:   PC_IS             *pcis = (PC_IS*)(pc->data);
1931:   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1932:   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
1933:   PetscErrorCode    ierr;
1934:   const PetscScalar one = 1.0;
1935:   const PetscScalar m_one = -1.0;
1936:   const PetscScalar zero = 0.0;

1939:   PetscCitationsRegister(citation,&cited);
1940:   if (pcbddc->ChangeOfBasisMatrix) {
1941:     Vec swap;

1943:     MatMultTranspose(pcbddc->ChangeOfBasisMatrix,r,pcbddc->work_change);
1944:     swap = pcbddc->work_change;
1945:     pcbddc->work_change = r;
1946:     r = swap;
1947:     /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */
1948:     if (pcbddc->benign_apply_coarse_only && pcbddc->exact_dirichlet_trick_app && pcbddc->change_interior) {
1949:       VecCopy(r,pcis->vec1_global);
1950:       VecLockPush(pcis->vec1_global);
1951:     }
1952:   }
1953:   if (pcbddc->benign_have_null) { /* get p0 from r */
1954:     PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);
1955:   }
1956:   if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
1957:     VecCopy(r,z);
1958:     /* First Dirichlet solve */
1959:     VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1960:     VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1961:     /*
1962:       Assembling right hand side for BDDC operator
1963:       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1964:       - pcis->vec1_B the interface part of the global vector z
1965:     */
1966:     if (n_D) {
1967:       KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1968:       VecScale(pcis->vec2_D,m_one);
1969:       if (pcbddc->switch_static) {
1970:         Mat_IS *matis = (Mat_IS*)(pc->mat->data);

1972:         VecSet(pcis->vec1_N,0.);
1973:         VecScatterBegin(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1974:         VecScatterEnd(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1975:         if (!pcbddc->switch_static_change) {
1976:           MatMultTranspose(matis->A,pcis->vec1_N,pcis->vec2_N);
1977:         } else {
1978:           MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
1979:           MatMultTranspose(matis->A,pcis->vec2_N,pcis->vec1_N);
1980:           MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
1981:         }
1982:         VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);
1983:         VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);
1984:         VecScatterBegin(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1985:         VecScatterEnd(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1986:       } else {
1987:         MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);
1988:       }
1989:     } else {
1990:       VecSet(pcis->vec1_B,zero);
1991:     }
1992:     VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1993:     VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1994:     PCBDDCScalingRestriction(pc,z,pcis->vec1_B);
1995:   } else {
1996:     PCBDDCScalingRestriction(pc,r,pcis->vec1_B);
1997:   }

1999:   /* Apply interface preconditioner
2000:      input/output vecs: pcis->vec1_B and pcis->vec1_D */
2001:   PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);

2003:   /* Apply transpose of partition of unity operator */
2004:   PCBDDCScalingExtension(pc,pcis->vec1_B,z);

2006:   /* Second Dirichlet solve and assembling of output */
2007:   VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
2008:   VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
2009:   if (n_B) {
2010:     if (pcbddc->switch_static) {
2011:       Mat_IS *matis = (Mat_IS*)(pc->mat->data);

2013:       VecScatterBegin(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
2014:       VecScatterEnd(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
2015:       VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
2016:       VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
2017:       if (!pcbddc->switch_static_change) {
2018:         MatMultTranspose(matis->A,pcis->vec1_N,pcis->vec2_N);
2019:       } else {
2020:         MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
2021:         MatMultTranspose(matis->A,pcis->vec2_N,pcis->vec1_N);
2022:         MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
2023:       }
2024:       VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);
2025:       VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);
2026:     } else {
2027:       MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);
2028:     }
2029:   } else if (pcbddc->switch_static) { /* n_B is zero */
2030:     Mat_IS *matis = (Mat_IS*)(pc->mat->data);

2032:     if (!pcbddc->switch_static_change) {
2033:       MatMultTranspose(matis->A,pcis->vec1_D,pcis->vec3_D);
2034:     } else {
2035:       MatMult(pcbddc->switch_static_change,pcis->vec1_D,pcis->vec1_N);
2036:       MatMultTranspose(matis->A,pcis->vec1_N,pcis->vec2_N);
2037:       MatMultTranspose(pcbddc->switch_static_change,pcis->vec2_N,pcis->vec3_D);
2038:     }
2039:   }
2040:   KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);
2041:   if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
2042:     if (pcbddc->switch_static) {
2043:       VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);
2044:     } else {
2045:       VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);
2046:     }
2047:     VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
2048:     VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
2049:   } else {
2050:     if (pcbddc->switch_static) {
2051:       VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);
2052:     } else {
2053:       VecScale(pcis->vec4_D,m_one);
2054:     }
2055:     VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
2056:     VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
2057:   }
2058:   if (pcbddc->benign_have_null) { /* set p0 (computed in PCBDDCApplyInterface) */
2059:     PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);
2060:   }
2061:   if (pcbddc->ChangeOfBasisMatrix) {
2062:     pcbddc->work_change = r;
2063:     VecCopy(z,pcbddc->work_change);
2064:     MatMult(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,z);
2065:   }
2066:   return(0);
2067: }

2069: PetscErrorCode PCReset_BDDC(PC pc)
2070: {
2071:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
2072:   PC_IS          *pcis = (PC_IS*)pc->data;
2073:   KSP            kspD,kspR,kspC;

2077:   /* free BDDC custom data  */
2078:   PCBDDCResetCustomization(pc);
2079:   /* destroy objects related to topography */
2080:   PCBDDCResetTopography(pc);
2081:   /* destroy objects for scaling operator */
2082:   PCBDDCScalingDestroy(pc);
2083:   /* free solvers stuff */
2084:   PCBDDCResetSolvers(pc);
2085:   /* free global vectors needed in presolve */
2086:   VecDestroy(&pcbddc->temp_solution);
2087:   VecDestroy(&pcbddc->original_rhs);
2088:   /* free data created by PCIS */
2089:   PCISDestroy(pc);

2091:   /* restore defaults */
2092:   kspD = pcbddc->ksp_D;
2093:   kspR = pcbddc->ksp_R;
2094:   kspC = pcbddc->coarse_ksp;
2095:   PetscMemzero(pc->data,sizeof(*pcbddc));
2096:   pcis->n_neigh                     = -1;
2097:   pcis->scaling_factor              = 1.0;
2098:   pcis->reusesubmatrices            = PETSC_TRUE;
2099:   pcbddc->use_local_adj             = PETSC_TRUE;
2100:   pcbddc->use_vertices              = PETSC_TRUE;
2101:   pcbddc->use_edges                 = PETSC_TRUE;
2102:   pcbddc->symmetric_primal          = PETSC_TRUE;
2103:   pcbddc->vertex_size               = 1;
2104:   pcbddc->recompute_topography      = PETSC_TRUE;
2105:   pcbddc->coarse_size               = -1;
2106:   pcbddc->use_exact_dirichlet_trick = PETSC_TRUE;
2107:   pcbddc->coarsening_ratio          = 8;
2108:   pcbddc->coarse_eqs_per_proc       = 1;
2109:   pcbddc->benign_compute_correction = PETSC_TRUE;
2110:   pcbddc->nedfield                  = -1;
2111:   pcbddc->nedglobal                 = PETSC_TRUE;
2112:   pcbddc->graphmaxcount             = PETSC_MAX_INT;
2113:   pcbddc->sub_schurs_layers         = -1;
2114:   pcbddc->ksp_D                     = kspD;
2115:   pcbddc->ksp_R                     = kspR;
2116:   pcbddc->coarse_ksp                = kspC;
2117:   return(0);
2118: }

2120: PetscErrorCode PCDestroy_BDDC(PC pc)
2121: {
2122:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

2126:   PCReset_BDDC(pc);
2127:   KSPDestroy(&pcbddc->ksp_D);
2128:   KSPDestroy(&pcbddc->ksp_R);
2129:   KSPDestroy(&pcbddc->coarse_ksp);
2130:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDiscreteGradient_C",NULL);
2131:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDivergenceMat_C",NULL);
2132:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);
2133:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);
2134:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",NULL);
2135:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);
2136:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);
2137:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);
2138:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);
2139:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);
2140:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);
2141:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);
2142:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);
2143:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);
2144:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);
2145:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);
2146:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);
2147:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);
2148:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);
2149:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);
2150:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);
2151:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);
2152:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);
2153:   PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);
2154:   PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);
2155:   PetscFree(pc->data);
2156:   return(0);
2157: }

2159: static PetscErrorCode PCSetCoordinates_BDDC(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
2160: {
2161:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
2162:   PCBDDCGraph    mat_graph = pcbddc->mat_graph;

2166:   PetscFree(mat_graph->coords);
2167:   PetscMalloc1(nloc*dim,&mat_graph->coords);
2168:   PetscMemcpy(mat_graph->coords,coords,nloc*dim*sizeof(PetscReal));
2169:   mat_graph->cnloc = nloc;
2170:   mat_graph->cdim  = dim;
2171:   mat_graph->cloc  = PETSC_FALSE;
2172:   return(0);
2173: }

2175: static PetscErrorCode PCPreSolveChangeRHS_BDDC(PC pc, PetscBool* change)
2176: {
2178:   *change = PETSC_TRUE;
2179:   return(0);
2180: }

2182: static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
2183: {
2184:   FETIDPMat_ctx  mat_ctx;
2185:   Vec            work;
2186:   PC_IS*         pcis;
2187:   PC_BDDC*       pcbddc;

2191:   MatShellGetContext(fetidp_mat,&mat_ctx);
2192:   pcis = (PC_IS*)mat_ctx->pc->data;
2193:   pcbddc = (PC_BDDC*)mat_ctx->pc->data;

2195:   VecSet(fetidp_flux_rhs,0.0);
2196:   /* copy rhs since we may change it during PCPreSolve_BDDC */
2197:   if (!pcbddc->original_rhs) {
2198:     VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);
2199:   }
2200:   if (mat_ctx->rhs_flip) {
2201:     VecPointwiseMult(pcbddc->original_rhs,standard_rhs,mat_ctx->rhs_flip);
2202:   } else {
2203:     VecCopy(standard_rhs,pcbddc->original_rhs);
2204:   }
2205:   if (mat_ctx->g2g_p) {
2206:     /* interface pressure rhs */
2207:     VecScatterBegin(mat_ctx->g2g_p,fetidp_flux_rhs,pcbddc->original_rhs,INSERT_VALUES,SCATTER_REVERSE);
2208:     VecScatterEnd(mat_ctx->g2g_p,fetidp_flux_rhs,pcbddc->original_rhs,INSERT_VALUES,SCATTER_REVERSE);
2209:     VecScatterBegin(mat_ctx->g2g_p,standard_rhs,fetidp_flux_rhs,INSERT_VALUES,SCATTER_FORWARD);
2210:     VecScatterEnd(mat_ctx->g2g_p,standard_rhs,fetidp_flux_rhs,INSERT_VALUES,SCATTER_FORWARD);
2211:     if (!mat_ctx->rhs_flip) {
2212:       VecScale(fetidp_flux_rhs,-1.);
2213:     }
2214:   }
2215:   /*
2216:      change of basis for physical rhs if needed
2217:      It also changes the rhs in case of dirichlet boundaries
2218:   */
2219:   PCPreSolve_BDDC(mat_ctx->pc,NULL,pcbddc->original_rhs,NULL);
2220:   if (pcbddc->ChangeOfBasisMatrix) {
2221:     MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcbddc->original_rhs,pcbddc->work_change);
2222:     work = pcbddc->work_change;
2223:    } else {
2224:     work = pcbddc->original_rhs;
2225:   }
2226:   /* store vectors for computation of fetidp final solution */
2227:   VecScatterBegin(pcis->global_to_D,work,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);
2228:   VecScatterEnd(pcis->global_to_D,work,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);
2229:   /* scale rhs since it should be unassembled */
2230:   /* TODO use counter scaling? (also below) */
2231:   VecScatterBegin(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
2232:   VecScatterEnd(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
2233:   /* Apply partition of unity */
2234:   VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);
2235:   /* PCBDDCScalingRestriction(mat_ctx->pc,work,mat_ctx->temp_solution_B); */
2236:   if (!pcbddc->switch_static) {
2237:     /* compute partially subassembled Schur complement right-hand side */
2238:     KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);
2239:     MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);
2240:     VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);
2241:     VecSet(work,0.0);
2242:     VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,work,ADD_VALUES,SCATTER_REVERSE);
2243:     VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,work,ADD_VALUES,SCATTER_REVERSE);
2244:     /* PCBDDCScalingRestriction(mat_ctx->pc,work,mat_ctx->temp_solution_B); */
2245:     VecScatterBegin(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
2246:     VecScatterEnd(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
2247:     VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);
2248:   }
2249:   /* BDDC rhs */
2250:   VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);
2251:   if (pcbddc->switch_static) {
2252:     VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);
2253:   }
2254:   /* apply BDDC */
2255:   PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));
2256:   PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
2257:   PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));

2259:   /* Application of B_delta and assembling of rhs for fetidp fluxes */
2260:   MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);
2261:   VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
2262:   VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
2263:   /* Add contribution to interface pressures */
2264:   if (mat_ctx->l2g_p) {
2265:     MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);
2266:     if (pcbddc->switch_static) {
2267:       MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);
2268:     }
2269:     VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
2270:     VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
2271:   }
2272:   return(0);
2273: }

2275: /*@
2276:  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side

2278:    Collective

2280:    Input Parameters:
2281: +  fetidp_mat      - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
2282: -  standard_rhs    - the right-hand side of the original linear system

2284:    Output Parameters:
2285: .  fetidp_flux_rhs - the right-hand side for the FETI-DP linear system

2287:    Level: developer

2289:    Notes:

2291: .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
2292: @*/
2293: PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
2294: {
2295:   FETIDPMat_ctx  mat_ctx;

2302:   MatShellGetContext(fetidp_mat,&mat_ctx);
2303:   PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));
2304:   return(0);
2305: }

2307: static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
2308: {
2309:   FETIDPMat_ctx  mat_ctx;
2310:   PC_IS*         pcis;
2311:   PC_BDDC*       pcbddc;
2313:   Vec            work;

2316:   MatShellGetContext(fetidp_mat,&mat_ctx);
2317:   pcis = (PC_IS*)mat_ctx->pc->data;
2318:   pcbddc = (PC_BDDC*)mat_ctx->pc->data;

2320:   /* apply B_delta^T */
2321:   VecSet(pcis->vec1_B,0.);
2322:   VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
2323:   VecScatterEnd(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
2324:   MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);
2325:   if (mat_ctx->l2g_p) {
2326:     VecScatterBegin(mat_ctx->l2g_p,fetidp_flux_sol,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);
2327:     VecScatterEnd(mat_ctx->l2g_p,fetidp_flux_sol,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);
2328:     MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);
2329:   }

2331:   /* compute rhs for BDDC application */
2332:   VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);
2333:   if (pcbddc->switch_static) {
2334:     VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);
2335:     if (mat_ctx->l2g_p) {
2336:       VecScale(mat_ctx->vP,-1.);
2337:       MatMultAdd(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D,pcis->vec1_D);
2338:     }
2339:   }

2341:   /* apply BDDC */
2342:   PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));
2343:   PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);

2345:   /* put values into global vector */
2346:   if (pcbddc->ChangeOfBasisMatrix) work = pcbddc->work_change;
2347:   else work = standard_sol;
2348:   VecScatterBegin(pcis->global_to_B,pcis->vec1_B,work,INSERT_VALUES,SCATTER_REVERSE);
2349:   VecScatterEnd(pcis->global_to_B,pcis->vec1_B,work,INSERT_VALUES,SCATTER_REVERSE);
2350:   if (!pcbddc->switch_static) {
2351:     /* compute values into the interior if solved for the partially subassembled Schur complement */
2352:     MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);
2353:     VecAYPX(pcis->vec1_D,-1.0,mat_ctx->temp_solution_D);
2354:     KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec1_D);
2355:   }

2357:   VecScatterBegin(pcis->global_to_D,pcis->vec1_D,work,INSERT_VALUES,SCATTER_REVERSE);
2358:   VecScatterEnd(pcis->global_to_D,pcis->vec1_D,work,INSERT_VALUES,SCATTER_REVERSE);
2359:   /* add p0 solution to final solution */
2360:   PCBDDCBenignGetOrSetP0(mat_ctx->pc,work,PETSC_FALSE);
2361:   if (pcbddc->ChangeOfBasisMatrix) {
2362:     MatMult(pcbddc->ChangeOfBasisMatrix,work,standard_sol);
2363:   }
2364:   PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);
2365:   if (mat_ctx->g2g_p) {
2366:     VecScatterBegin(mat_ctx->g2g_p,fetidp_flux_sol,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
2367:     VecScatterEnd(mat_ctx->g2g_p,fetidp_flux_sol,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
2368:   }
2369:   return(0);
2370: }

2372: /*@
2373:  PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system

2375:    Collective

2377:    Input Parameters:
2378: +  fetidp_mat      - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
2379: -  fetidp_flux_sol - the solution of the FETI-DP linear system

2381:    Output Parameters:
2382: .  standard_sol    - the solution defined on the physical domain

2384:    Level: developer

2386:    Notes:

2388: .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
2389: @*/
2390: PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
2391: {
2392:   FETIDPMat_ctx  mat_ctx;

2399:   MatShellGetContext(fetidp_mat,&mat_ctx);
2400:   PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));
2401:   return(0);
2402: }

2404: static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, PetscBool fully_redundant, const char* prefix, Mat *fetidp_mat, PC *fetidp_pc)
2405: {

2407:   FETIDPMat_ctx  fetidpmat_ctx;
2408:   Mat            newmat;
2409:   FETIDPPC_ctx   fetidppc_ctx;
2410:   PC             newpc;
2411:   MPI_Comm       comm;

2415:   PetscObjectGetComm((PetscObject)pc,&comm);
2416:   /* FETIDP linear matrix */
2417:   PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);
2418:   fetidpmat_ctx->fully_redundant = fully_redundant;
2419:   PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);
2420:   MatCreateShell(comm,fetidpmat_ctx->n,fetidpmat_ctx->n,fetidpmat_ctx->N,fetidpmat_ctx->N,fetidpmat_ctx,&newmat);
2421:   MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);
2422:   MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);
2423:   MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);
2424:   MatSetOptionsPrefix(newmat,prefix);
2425:   MatAppendOptionsPrefix(newmat,"fetidp_");
2426:   MatSetUp(newmat);
2427:   /* FETIDP preconditioner */
2428:   PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);
2429:   PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);
2430:   PCCreate(comm,&newpc);
2431:   PCSetOperators(newpc,newmat,newmat);
2432:   PCSetOptionsPrefix(newpc,prefix);
2433:   PCAppendOptionsPrefix(newpc,"fetidp_");
2434:   if (!fetidpmat_ctx->l2g_lambda_only) {
2435:     PCSetType(newpc,PCSHELL);
2436:     PCShellSetContext(newpc,fetidppc_ctx);
2437:     PCShellSetApply(newpc,FETIDPPCApply);
2438:     PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);
2439:     PCShellSetView(newpc,FETIDPPCView);
2440:     PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);
2441:   } else {
2442:     KSP      *ksps;
2443:     PC       lagpc;
2444:     Mat      M,AM,PM;
2445:     PetscInt nn;

2447:     PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_PPmat",(PetscObject*)&M);
2448:     PCSetType(newpc,PCFIELDSPLIT);
2449:     PCFieldSplitSetIS(newpc,"lag",fetidpmat_ctx->lagrange);
2450:     PCFieldSplitSetIS(newpc,"p",fetidpmat_ctx->pressure);
2451:     PCFieldSplitSetType(newpc,PC_COMPOSITE_SCHUR);
2452:     PCFieldSplitSetSchurFactType(newpc,PC_FIELDSPLIT_SCHUR_FACT_DIAG);
2453:     PCFieldSplitSetSchurPre(newpc,PC_FIELDSPLIT_SCHUR_PRE_USER,M);
2454:     PCFieldSplitSetSchurScale(newpc,1.0);
2455:     PCSetFromOptions(newpc);
2456:     PCSetUp(newpc);

2458:     /* set the solver for the (0,0) block */
2459:     PCFieldSplitGetSubKSP(newpc,&nn,&ksps);
2460:     PCCreate(comm,&lagpc);
2461:     PCSetType(lagpc,PCSHELL);
2462:     KSPGetOperators(ksps[0],&AM,&PM);
2463:     PCSetOperators(lagpc,AM,PM);
2464:     PCShellSetContext(lagpc,fetidppc_ctx);
2465:     PCShellSetApply(lagpc,FETIDPPCApply);
2466:     PCShellSetApplyTranspose(lagpc,FETIDPPCApplyTranspose);
2467:     PCShellSetView(lagpc,FETIDPPCView);
2468:     PCShellSetDestroy(lagpc,PCBDDCDestroyFETIDPPC);
2469:     KSPSetPC(ksps[0],lagpc);
2470:     PCDestroy(&lagpc);
2471:     PetscFree(ksps);
2472:   }
2473:   /* return pointers for objects created */
2474:   *fetidp_mat = newmat;
2475:   *fetidp_pc = newpc;
2476:   return(0);
2477: }

2479: /*@C
2480:  PCBDDCCreateFETIDPOperators - Create FETI-DP operators

2482:    Collective

2484:    Input Parameters:
2485: +  pc - the BDDC preconditioning context (setup should have been called before)
2486: .  fully_redundant - true for a fully redundant set of Lagrange multipliers
2487: -  prefix - optional options database prefix for the objects to be created (can be NULL)

2489:    Output Parameters:
2490: +  fetidp_mat - shell FETI-DP matrix object
2491: -  fetidp_pc  - shell Dirichlet preconditioner for FETI-DP matrix

2493:    Level: developer

2495:    Notes:
2496:      Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose

2498: .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
2499: @*/
2500: PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, PetscBool fully_redundant, const char *prefix, Mat *fetidp_mat, PC *fetidp_pc)
2501: {

2506:   if (pc->setupcalled) {
2507:     PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,PetscBool,const char*,Mat*,PC*),(pc,fully_redundant,prefix,fetidp_mat,fetidp_pc));
2508:   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
2509:   return(0);
2510: }
2511: /* -------------------------------------------------------------------------- */
2512: /*MC
2513:    PCBDDC - Balancing Domain Decomposition by Constraints.

2515:    An implementation of the BDDC preconditioner based on

2517: .vb
2518:    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
2519:    [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
2520:    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
2521:    [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
2522: .ve

2524:    The matrix to be preconditioned (Pmat) must be of type MATIS.

2526:    Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.

2528:    It also works with unsymmetric and indefinite problems.

2530:    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.

2532:    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).

2534:    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()
2535:    Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesIS() and their local counterparts.

2537:    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.

2539:    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.
2540:    User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()

2542:    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.

2544:    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.

2546:    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.
2547:    Deluxe scaling is not supported yet for FETI-DP.

2549:    Options Database Keys (some of them, run with -h for a complete list):

2551: .    -pc_bddc_use_vertices <true> - use or not vertices in primal space
2552: .    -pc_bddc_use_edges <true> - use or not edges in primal space
2553: .    -pc_bddc_use_faces <false> - use or not faces in primal space
2554: .    -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
2555: .    -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
2556: .    -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
2557: .    -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
2558: .    -pc_bddc_levels <0> - maximum number of levels for multilevel
2559: .    -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)
2560: .    -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)
2561: .    -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
2562: .    -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)
2563: .    -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)
2564: -    -pc_bddc_check_level <0> - set verbosity level of debugging output

2566:    Options for Dirichlet, Neumann or coarse solver can be set with
2567: .vb
2568:       -pc_bddc_dirichlet_
2569:       -pc_bddc_neumann_
2570:       -pc_bddc_coarse_
2571: .ve
2572:    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.

2574:    When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
2575: .vb
2576:       -pc_bddc_dirichlet_lN_
2577:       -pc_bddc_neumann_lN_
2578:       -pc_bddc_coarse_lN_
2579: .ve
2580:    Note that level number ranges from the finest (0) to the coarsest (N).
2581:    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.
2582: .vb
2583:      -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
2584: .ve
2585:    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

2587:    Level: intermediate

2589:    Developer notes:

2591:    Contributed by Stefano Zampini

2593: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
2594: M*/

2596: PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
2597: {
2598:   PetscErrorCode      ierr;
2599:   PC_BDDC             *pcbddc;

2602:   PetscNewLog(pc,&pcbddc);
2603:   pc->data = (void*)pcbddc;

2605:   /* create PCIS data structure */
2606:   PCISCreate(pc);

2608:   /* create local graph structure */
2609:   PCBDDCGraphCreate(&pcbddc->mat_graph);

2611:   /* BDDC nonzero defaults */
2612:   pcbddc->use_local_adj             = PETSC_TRUE;
2613:   pcbddc->use_vertices              = PETSC_TRUE;
2614:   pcbddc->use_edges                 = PETSC_TRUE;
2615:   pcbddc->symmetric_primal          = PETSC_TRUE;
2616:   pcbddc->vertex_size               = 1;
2617:   pcbddc->recompute_topography      = PETSC_TRUE;
2618:   pcbddc->coarse_size               = -1;
2619:   pcbddc->use_exact_dirichlet_trick = PETSC_TRUE;
2620:   pcbddc->coarsening_ratio          = 8;
2621:   pcbddc->coarse_eqs_per_proc       = 1;
2622:   pcbddc->benign_compute_correction = PETSC_TRUE;
2623:   pcbddc->nedfield                  = -1;
2624:   pcbddc->nedglobal                 = PETSC_TRUE;
2625:   pcbddc->graphmaxcount             = PETSC_MAX_INT;
2626:   pcbddc->sub_schurs_layers         = -1;
2627:   pcbddc->adaptive_threshold[0]     = 0.0;
2628:   pcbddc->adaptive_threshold[1]     = 0.0;

2630:   /* function pointers */
2631:   pc->ops->apply               = PCApply_BDDC;
2632:   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2633:   pc->ops->setup               = PCSetUp_BDDC;
2634:   pc->ops->destroy             = PCDestroy_BDDC;
2635:   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
2636:   pc->ops->view                = PCView_BDDC;
2637:   pc->ops->applyrichardson     = 0;
2638:   pc->ops->applysymmetricleft  = 0;
2639:   pc->ops->applysymmetricright = 0;
2640:   pc->ops->presolve            = PCPreSolve_BDDC;
2641:   pc->ops->postsolve           = PCPostSolve_BDDC;
2642:   pc->ops->reset               = PCReset_BDDC;

2644:   /* composing function */
2645:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDiscreteGradient_C",PCBDDCSetDiscreteGradient_BDDC);
2646:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDivergenceMat_C",PCBDDCSetDivergenceMat_BDDC);
2647:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);
2648:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);
2649:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",PCBDDCSetPrimalVerticesIS_BDDC);
2650:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);
2651:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);
2652:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);
2653:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);
2654:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);
2655:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);
2656:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);
2657:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);
2658:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);
2659:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);
2660:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);
2661:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);
2662:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);
2663:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);
2664:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);
2665:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);
2666:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);
2667:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);
2668:   PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_BDDC);
2669:   PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_BDDC);
2670:   return(0);
2671: }