Actual source code: bddc.c

petsc-master 2016-12-08
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: /* -------------------------------------------------------------------------- */
 29: PetscErrorCode PCSetFromOptions_BDDC(PetscOptionItems *PetscOptionsObject,PC pc)
 30: {
 31:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

 35:   PetscOptionsHead(PetscOptionsObject,"BDDC options");
 36:   /* Verbose debugging */
 37:   PetscOptionsBool("-pc_bddc_dirichlet_approximate","Inform PCBDDC that we are using approximate Dirichlet solvers","none",pcbddc->NullSpace_corr[0],&pcbddc->NullSpace_corr[0],NULL);
 38:   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);
 39:   PetscOptionsBool("-pc_bddc_neumann_approximate","Inform PCBDDC that we are using approximate Neumann solvers","none",pcbddc->NullSpace_corr[2],&pcbddc->NullSpace_corr[2],NULL);
 40:   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);
 41:   PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);
 42:   /* Primal space customization */
 43:   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);
 44:   PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);
 45:   PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);
 46:   PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);
 47:   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);
 48:   PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is used when multiple constraints are present)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);
 49:   /* Change of basis */
 50:   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);
 51:   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);
 52:   if (!pcbddc->use_change_of_basis) {
 53:     pcbddc->use_change_on_faces = PETSC_FALSE;
 54:   }
 55:   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
 56:   PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);
 57:   PetscOptionsInt("-pc_bddc_coarse_redistribute","Number of procs where to redistribute coarse problem","none",pcbddc->redistribute_coarse,&pcbddc->redistribute_coarse,NULL);
 58:   PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);
 59:   PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);
 60:   PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);
 61:   PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);
 62:   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);
 63:   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);
 64:   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);
 65:   PetscOptionsBool("-pc_bddc_deluxe_faster","Faster application of deluxe scaling (requires extra work during setup)","none",pcbddc->faster_deluxe,&pcbddc->faster_deluxe,NULL);
 66:   PetscOptionsReal("-pc_bddc_adaptive_threshold","Threshold to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&pcbddc->adaptive_threshold,NULL);
 67:   PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);
 68:   PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);
 69:   PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);
 70:   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);
 71:   PetscOptionsTail();
 72:   return(0);
 73: }

 75: /* -------------------------------------------------------------------------- */
 78: static PetscErrorCode PCView_BDDC(PC pc,PetscViewer viewer)
 79: {
 80:   PC_BDDC              *pcbddc = (PC_BDDC*)pc->data;
 81:   PetscErrorCode       ierr;
 82:   PetscBool            isascii,isstring;

 85:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 86:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);

 88:   /* In a braindead way, print out anything which the user can control from the command line,
 89:      cribbing from PCSetFromOptions_BDDC */

 91:   /* Nothing printed for the String viewer */

 93:   /* ASCII viewer */
 94:   if (isascii) {
 95:     /* Verbose debugging */
 96:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use verbose output: %d\n",pcbddc->dbg_flag);

 98:     /* Primal space customization */
 99:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use local mat graph: %d\n",pcbddc->use_local_adj);
100:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use vertices: %d\n",pcbddc->use_vertices);
101:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use edges: %d\n",pcbddc->use_edges);
102:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use faces: %d\n",pcbddc->use_faces);
103:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use true near null space: %d\n",pcbddc->use_nnsp_true);
104:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use QR for single constraints on cc: %d\n",pcbddc->use_qr_single);

106:     /* Change of basis */
107:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use change of basis on local edge nodes: %d\n",pcbddc->use_change_of_basis);
108:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use change of basis on local face nodes: %d\n",pcbddc->use_change_on_faces);

110:     /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
111:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Switch on static condensation ops around the interface preconditioner: %d\n",pcbddc->switch_static);
112:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Coarse problem restribute procs: %d\n",pcbddc->redistribute_coarse);
113:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Multilevel coarsening ratio: %d\n",pcbddc->coarsening_ratio);
114:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Multilevel max levels: %d\n",pcbddc->max_levels);
115:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use estimated eigs for coarse problem: %d\n",pcbddc->use_coarse_estimates);
116:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use deluxe scaling: %d\n",pcbddc->use_deluxe_scaling);
117:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Rebuild interface graph for Schur principal minors: %d\n",pcbddc->sub_schurs_rebuild);
118:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Number of dofs' layers for the computation of principal minors: %d\n",pcbddc->sub_schurs_layers);
119:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Use user CSR graph to compute successive layers: %d\n",pcbddc->sub_schurs_use_useradj);
120:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Fast deluxe scaling: %d\n",pcbddc->faster_deluxe);
121:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Adaptive constraint selection threshold: %g\n",pcbddc->adaptive_threshold);
122:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Min constraints / connected component: %d\n",pcbddc->adaptive_nmin);
123:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Max constraints / connected component: %d\n",pcbddc->adaptive_nmax);
124:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Symmetric computation of primal basis functions: %d\n",pcbddc->symmetric_primal);
125:     PetscViewerASCIIPrintf(viewer,    "  BDDC: Num. Procs. to map coarse adjacency list: %d\n",pcbddc->coarse_adj_red);
126:   }
127:   return(0);
128: }

130: /* -------------------------------------------------------------------------- */
133: static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
134: {
135:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

139:   MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);
140:   PetscObjectReference((PetscObject)change);
141:   pcbddc->user_ChangeOfBasisMatrix = change;
142:   return(0);
143: }
146: /*@
147:  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs

149:    Collective on PC

151:    Input Parameters:
152: +  pc - the preconditioning context
153: -  change - the change of basis matrix

155:    Level: intermediate

157:    Notes:

159: .seealso: PCBDDC
160: @*/
161: PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
162: {

169:   if (pc->mat) {
170:     PetscInt rows_c,cols_c,rows,cols;
171:     MatGetSize(pc->mat,&rows,&cols);
172:     MatGetSize(change,&rows_c,&cols_c);
173:     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);
174:     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);
175:     MatGetLocalSize(pc->mat,&rows,&cols);
176:     MatGetLocalSize(change,&rows_c,&cols_c);
177:     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);
178:     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);
179:   }
180:   PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));
181:   return(0);
182: }
183: /* -------------------------------------------------------------------------- */
186: static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
187: {
188:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

192:   ISDestroy(&pcbddc->user_primal_vertices);
193:   PetscObjectReference((PetscObject)PrimalVertices);
194:   pcbddc->user_primal_vertices = PrimalVertices;
195:   return(0);
196: }
199: /*@
200:  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC

202:    Collective

204:    Input Parameters:
205: +  pc - the preconditioning context
206: -  PrimalVertices - index set of primal vertices in local numbering (can be empty)

208:    Level: intermediate

210:    Notes:

212: .seealso: PCBDDC
213: @*/
214: PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
215: {

222:   PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));
223:   return(0);
224: }
225: /* -------------------------------------------------------------------------- */
228: static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
229: {
230:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

233:   pcbddc->coarsening_ratio = k;
234:   return(0);
235: }

239: /*@
240:  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel

242:    Logically collective on PC

244:    Input Parameters:
245: +  pc - the preconditioning context
246: -  k - coarsening ratio (H/h at the coarser level)

248:    Options Database Keys:
249: .    -pc_bddc_coarsening_ratio

251:    Level: intermediate

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

256: .seealso: PCBDDC, PCBDDCSetLevels()
257: @*/
258: PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
259: {

265:   PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));
266:   return(0);
267: }

269: /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
272: static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
273: {
274:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

277:   pcbddc->use_exact_dirichlet_trick = flg;
278:   return(0);
279: }

283: PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
284: {

290:   PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));
291:   return(0);
292: }

296: static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
297: {
298:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

301:   pcbddc->current_level = level;
302:   return(0);
303: }

307: PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
308: {

314:   PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));
315:   return(0);
316: }

320: static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
321: {
322:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

325:   pcbddc->max_levels = levels;
326:   return(0);
327: }

331: /*@
332:  PCBDDCSetLevels - Sets the maximum number of levels for multilevel

334:    Logically collective on PC

336:    Input Parameters:
337: +  pc - the preconditioning context
338: -  levels - the maximum number of levels (max 9)

340:    Options Database Keys:
341: .    -pc_bddc_levels

343:    Level: intermediate

345:    Notes:
346:      Default value is 0, i.e. traditional one-level BDDC

348: .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
349: @*/
350: PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
351: {

357:   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
358:   PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));
359:   return(0);
360: }
361: /* -------------------------------------------------------------------------- */

365: static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
366: {
367:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

371:   /* last user setting takes precendence -> destroy any other customization */
372:   ISDestroy(&pcbddc->DirichletBoundariesLocal);
373:   ISDestroy(&pcbddc->DirichletBoundaries);
374:   PetscObjectReference((PetscObject)DirichletBoundaries);
375:   pcbddc->DirichletBoundaries = DirichletBoundaries;
376:   pcbddc->recompute_topography = PETSC_TRUE;
377:   return(0);
378: }

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

385:    Collective

387:    Input Parameters:
388: +  pc - the preconditioning context
389: -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries

391:    Level: intermediate

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

396: .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
397: @*/
398: PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
399: {

406:   PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));
407:   return(0);
408: }
409: /* -------------------------------------------------------------------------- */

413: static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
414: {
415:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

419:   /* last user setting takes precendence -> destroy any other customization */
420:   ISDestroy(&pcbddc->DirichletBoundariesLocal);
421:   ISDestroy(&pcbddc->DirichletBoundaries);
422:   PetscObjectReference((PetscObject)DirichletBoundaries);
423:   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
424:   pcbddc->recompute_topography = PETSC_TRUE;
425:   return(0);
426: }

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

433:    Collective

435:    Input Parameters:
436: +  pc - the preconditioning context
437: -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)

439:    Level: intermediate

441:    Notes:

443: .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
444: @*/
445: PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
446: {

453:   PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));
454:   return(0);
455: }
456: /* -------------------------------------------------------------------------- */

460: static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
461: {
462:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

466:   /* last user setting takes precendence -> destroy any other customization */
467:   ISDestroy(&pcbddc->NeumannBoundariesLocal);
468:   ISDestroy(&pcbddc->NeumannBoundaries);
469:   PetscObjectReference((PetscObject)NeumannBoundaries);
470:   pcbddc->NeumannBoundaries = NeumannBoundaries;
471:   pcbddc->recompute_topography = PETSC_TRUE;
472:   return(0);
473: }

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

480:    Collective

482:    Input Parameters:
483: +  pc - the preconditioning context
484: -  NeumannBoundaries - parallel IS defining the Neumann boundaries

486:    Level: intermediate

488:    Notes:
489:      Any process can list any global node

491: .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
492: @*/
493: PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
494: {

501:   PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));
502:   return(0);
503: }
504: /* -------------------------------------------------------------------------- */

508: static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
509: {
510:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

514:   /* last user setting takes precendence -> destroy any other customization */
515:   ISDestroy(&pcbddc->NeumannBoundariesLocal);
516:   ISDestroy(&pcbddc->NeumannBoundaries);
517:   PetscObjectReference((PetscObject)NeumannBoundaries);
518:   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
519:   pcbddc->recompute_topography = PETSC_TRUE;
520:   return(0);
521: }

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

528:    Collective

530:    Input Parameters:
531: +  pc - the preconditioning context
532: -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)

534:    Level: intermediate

536:    Notes:

538: .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
539: @*/
540: PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
541: {

548:   PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));
549:   return(0);
550: }
551: /* -------------------------------------------------------------------------- */

555: static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
556: {
557:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

560:   *DirichletBoundaries = pcbddc->DirichletBoundaries;
561:   return(0);
562: }

566: /*@
567:  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries

569:    Collective

571:    Input Parameters:
572: .  pc - the preconditioning context

574:    Output Parameters:
575: .  DirichletBoundaries - index set defining the Dirichlet boundaries

577:    Level: intermediate

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

582: .seealso: PCBDDC
583: @*/
584: PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
585: {

590:   PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));
591:   return(0);
592: }
593: /* -------------------------------------------------------------------------- */

597: static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
598: {
599:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

602:   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
603:   return(0);
604: }

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

611:    Collective

613:    Input Parameters:
614: .  pc - the preconditioning context

616:    Output Parameters:
617: .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries

619:    Level: intermediate

621:    Notes:
622:      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).
623:           In the latter case, the IS will be available after PCSetUp.

625: .seealso: PCBDDC
626: @*/
627: PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
628: {

633:   PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));
634:   return(0);
635: }
636: /* -------------------------------------------------------------------------- */

640: static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
641: {
642:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

645:   *NeumannBoundaries = pcbddc->NeumannBoundaries;
646:   return(0);
647: }

651: /*@
652:  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries

654:    Collective

656:    Input Parameters:
657: .  pc - the preconditioning context

659:    Output Parameters:
660: .  NeumannBoundaries - index set defining the Neumann boundaries

662:    Level: intermediate

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

667: .seealso: PCBDDC
668: @*/
669: PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
670: {

675:   PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));
676:   return(0);
677: }
678: /* -------------------------------------------------------------------------- */

682: static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
683: {
684:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

687:   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
688:   return(0);
689: }

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

696:    Collective

698:    Input Parameters:
699: .  pc - the preconditioning context

701:    Output Parameters:
702: .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries

704:    Level: intermediate

706:    Notes:
707:      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).
708:           In the latter case, the IS will be available after PCSetUp.

710: .seealso: PCBDDC
711: @*/
712: PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
713: {

718:   PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));
719:   return(0);
720: }
721: /* -------------------------------------------------------------------------- */

725: static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
726: {
727:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
728:   PCBDDCGraph    mat_graph = pcbddc->mat_graph;

732:   /* free old CSR */
733:   PCBDDCGraphResetCSR(mat_graph);
734:   /* TODO: PCBDDCGraphSetAdjacency */
735:   /* get CSR into graph structure */
736:   if (copymode == PETSC_COPY_VALUES) {
737:     PetscMalloc1(nvtxs+1,&mat_graph->xadj);
738:     PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);
739:     PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));
740:     PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));
741:   } else if (copymode == PETSC_OWN_POINTER) {
742:     mat_graph->xadj = (PetscInt*)xadj;
743:     mat_graph->adjncy = (PetscInt*)adjncy;
744:   }
745:   mat_graph->nvtxs_csr = nvtxs;
746:   return(0);
747: }

751: /*@
752:  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local matrix

754:    Not collective

756:    Input Parameters:
757: +  pc - the preconditioning context
758: .  nvtxs - number of local vertices of the graph (i.e., the size of the local problem)
759: .  xadj, adjncy - the CSR graph
760: -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.

762:    Level: intermediate

764:    Notes:

766: .seealso: PCBDDC,PetscCopyMode
767: @*/
768: PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
769: {
770:   void (*f)(void) = 0;

777:   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER)  SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d",copymode);
778:   PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));
779:   /* free arrays if PCBDDC is not the PC type */
780:   PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);
781:   if (!f && copymode == PETSC_OWN_POINTER) {
782:     PetscFree(xadj);
783:     PetscFree(adjncy);
784:   }
785:   return(0);
786: }
787: /* -------------------------------------------------------------------------- */

791: static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
792: {
793:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
794:   PetscInt i;

798:   /* Destroy ISes if they were already set */
799:   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
800:     ISDestroy(&pcbddc->ISForDofsLocal[i]);
801:   }
802:   PetscFree(pcbddc->ISForDofsLocal);
803:   /* last user setting takes precendence -> destroy any other customization */
804:   for (i=0;i<pcbddc->n_ISForDofs;i++) {
805:     ISDestroy(&pcbddc->ISForDofs[i]);
806:   }
807:   PetscFree(pcbddc->ISForDofs);
808:   pcbddc->n_ISForDofs = 0;
809:   /* allocate space then set */
810:   if (n_is) {
811:     PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);
812:   }
813:   for (i=0;i<n_is;i++) {
814:     PetscObjectReference((PetscObject)ISForDofs[i]);
815:     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
816:   }
817:   pcbddc->n_ISForDofsLocal=n_is;
818:   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
819:   pcbddc->recompute_topography = PETSC_TRUE;
820:   return(0);
821: }

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

828:    Collective

830:    Input Parameters:
831: +  pc - the preconditioning context
832: .  n_is - number of index sets defining the fields
833: -  ISForDofs - array of IS describing the fields in local ordering

835:    Level: intermediate

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

840: .seealso: PCBDDC
841: @*/
842: PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
843: {
844:   PetscInt       i;

850:   for (i=0;i<n_is;i++) {
853:   }
854:   PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));
855:   return(0);
856: }
857: /* -------------------------------------------------------------------------- */

861: static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
862: {
863:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
864:   PetscInt i;

868:   /* Destroy ISes if they were already set */
869:   for (i=0;i<pcbddc->n_ISForDofs;i++) {
870:     ISDestroy(&pcbddc->ISForDofs[i]);
871:   }
872:   PetscFree(pcbddc->ISForDofs);
873:   /* last user setting takes precendence -> destroy any other customization */
874:   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
875:     ISDestroy(&pcbddc->ISForDofsLocal[i]);
876:   }
877:   PetscFree(pcbddc->ISForDofsLocal);
878:   pcbddc->n_ISForDofsLocal = 0;
879:   /* allocate space then set */
880:   if (n_is) {
881:     PetscMalloc1(n_is,&pcbddc->ISForDofs);
882:   }
883:   for (i=0;i<n_is;i++) {
884:     PetscObjectReference((PetscObject)ISForDofs[i]);
885:     pcbddc->ISForDofs[i]=ISForDofs[i];
886:   }
887:   pcbddc->n_ISForDofs=n_is;
888:   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
889:   pcbddc->recompute_topography = PETSC_TRUE;
890:   return(0);
891: }

895: /*@
896:  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix

898:    Collective

900:    Input Parameters:
901: +  pc - the preconditioning context
902: .  n_is - number of index sets defining the fields
903: -  ISForDofs - array of IS describing the fields in global ordering

905:    Level: intermediate

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

910: .seealso: PCBDDC
911: @*/
912: PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
913: {
914:   PetscInt       i;

920:   for (i=0;i<n_is;i++) {
923:   }
924:   PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));
925:   return(0);
926: }

928: /* -------------------------------------------------------------------------- */
931: /* -------------------------------------------------------------------------- */
932: /*
933:    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
934:                      guess if a transformation of basis approach has been selected.

936:    Input Parameter:
937: +  pc - the preconditioner contex

939:    Application Interface Routine: PCPreSolve()

941:    Notes:
942:      The interface routine PCPreSolve() is not usually called directly by
943:    the user, but instead is called by KSPSolve().
944: */
945: static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
946: {
948:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
949:   PC_IS          *pcis = (PC_IS*)(pc->data);
950:   Vec            used_vec;
951:   PetscBool      copy_rhs = PETSC_TRUE;

954:   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
955:   if (ksp) {
956:     PetscBool iscg;
957:     PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);
958:     if (!iscg) {
959:       PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);
960:     }
961:   }

963:   /* Creates parallel work vectors used in presolve */
964:   if (!pcbddc->original_rhs) {
965:     VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);
966:   }
967:   if (!pcbddc->temp_solution) {
968:     VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);
969:   }

971:   if (x) {
972:     PetscObjectReference((PetscObject)x);
973:     used_vec = x;
974:   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
975:     PetscObjectReference((PetscObject)pcbddc->temp_solution);
976:     used_vec = pcbddc->temp_solution;
977:     VecSet(used_vec,0.0);
978:   }

980:   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
981:   if (ksp) {
982:     /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
983:     KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);
984:     if (!pcbddc->ksp_guess_nonzero) {
985:       VecSet(used_vec,0.0);
986:     }
987:   }

989:   pcbddc->rhs_change = PETSC_FALSE;

991:   /* Take into account zeroed rows -> change rhs and store solution removed */
992:   if (rhs) {
993:     IS dirIS = NULL;

995:     /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
996:     PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);
997:     if (dirIS) {
998:       Mat_IS            *matis = (Mat_IS*)pc->pmat->data;
999:       PetscInt          dirsize,i,*is_indices;
1000:       PetscScalar       *array_x;
1001:       const PetscScalar *array_diagonal;

1003:       MatGetDiagonal(pc->pmat,pcis->vec1_global);
1004:       VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);
1005:       VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);
1006:       VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);
1007:       VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
1008:       VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
1009:       ISGetLocalSize(dirIS,&dirsize);
1010:       VecGetArray(pcis->vec1_N,&array_x);
1011:       VecGetArrayRead(pcis->vec2_N,&array_diagonal);
1012:       ISGetIndices(dirIS,(const PetscInt**)&is_indices);
1013:       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
1014:       ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);
1015:       VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);
1016:       VecRestoreArray(pcis->vec1_N,&array_x);
1017:       VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);
1018:       VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);
1019:       pcbddc->rhs_change = PETSC_TRUE;
1020:       ISDestroy(&dirIS);
1021:     }
1022:   }

1024:   /* remove the computed solution or the initial guess from the rhs */
1025:   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
1026:     /* store the original rhs */
1027:     if (copy_rhs) {
1028:       VecCopy(rhs,pcbddc->original_rhs);
1029:       copy_rhs = PETSC_FALSE;
1030:     }
1031:     pcbddc->rhs_change = PETSC_TRUE;
1032:     VecScale(used_vec,-1.0);
1033:     MatMultAdd(pc->pmat,used_vec,rhs,rhs);
1034:     VecScale(used_vec,-1.0);
1035:     VecCopy(used_vec,pcbddc->temp_solution);
1036:     if (ksp) {
1037:       KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);
1038:     }
1039:   }
1040:   VecDestroy(&used_vec);

1042:   /* store partially computed solution and set initial guess */
1043:   if (x && pcbddc->use_exact_dirichlet_trick) {
1044:     VecSet(x,0.0);
1045:     VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1046:     VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1047:     KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1048:     VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);
1049:     VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);
1050:     if (ksp) {
1051:       KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
1052:     }
1053:   }

1055:   if (pcbddc->ChangeOfBasisMatrix) {
1056:     PCBDDCChange_ctx change_ctx;

1058:     /* get change ctx */
1059:     MatShellGetContext(pcbddc->new_global_mat,&change_ctx);

1061:     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1062:     MatDestroy(&change_ctx->original_mat);
1063:     PetscObjectReference((PetscObject)pc->mat);
1064:     change_ctx->original_mat = pc->mat;

1066:     /* change iteration matrix */
1067:     MatDestroy(&pc->mat);
1068:     PetscObjectReference((PetscObject)pcbddc->new_global_mat);
1069:     pc->mat = pcbddc->new_global_mat;

1071:     /* store the original rhs */
1072:     if (copy_rhs) {
1073:       VecCopy(rhs,pcbddc->original_rhs);
1074:     }

1076:     /* change rhs */
1077:     MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);
1078:     VecCopy(pcis->vec1_global,rhs);
1079:     pcbddc->rhs_change = PETSC_TRUE;
1080:   }
1081:   return(0);
1082: }

1084: /* -------------------------------------------------------------------------- */
1087: /* -------------------------------------------------------------------------- */
1088: /*
1089:    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1090:                      approach has been selected. Also, restores rhs to its original state.

1092:    Input Parameter:
1093: +  pc - the preconditioner contex

1095:    Application Interface Routine: PCPostSolve()

1097:    Notes:
1098:      The interface routine PCPostSolve() is not usually called directly by
1099:      the user, but instead is called by KSPSolve().
1100: */
1101: static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1102: {
1104:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

1107:   if (pcbddc->ChangeOfBasisMatrix) {
1108:     PCBDDCChange_ctx change_ctx;

1110:     /* get change ctx */
1111:     MatShellGetContext(pcbddc->new_global_mat,&change_ctx);

1113:     /* restore iteration matrix */
1114:     MatDestroy(&pc->mat);
1115:     PetscObjectReference((PetscObject)change_ctx->original_mat);
1116:     pc->mat = change_ctx->original_mat;

1118:     /* get solution in original basis */
1119:     if (x) {
1120:       PC_IS *pcis = (PC_IS*)(pc->data);
1121:       MatMult(change_ctx->global_change,x,pcis->vec1_global);
1122:       VecCopy(pcis->vec1_global,x);
1123:     }
1124:   }

1126:   /* add solution removed in presolve */
1127:   if (x && pcbddc->rhs_change) {
1128:     VecAXPY(x,1.0,pcbddc->temp_solution);
1129:   }

1131:   /* restore rhs to its original state */
1132:   if (rhs && pcbddc->rhs_change) {
1133:     VecCopy(pcbddc->original_rhs,rhs);
1134:   }
1135:   pcbddc->rhs_change = PETSC_FALSE;

1137:   /* restore ksp guess state */
1138:   if (ksp) {
1139:     KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);
1140:   }
1141:   return(0);
1142: }
1143: /* -------------------------------------------------------------------------- */
1146: /* -------------------------------------------------------------------------- */
1147: /*
1148:    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1149:                   by setting data structures and options.

1151:    Input Parameter:
1152: +  pc - the preconditioner context

1154:    Application Interface Routine: PCSetUp()

1156:    Notes:
1157:      The interface routine PCSetUp() is not usually called directly by
1158:      the user, but instead is called by PCApply() if necessary.
1159: */
1160: PetscErrorCode PCSetUp_BDDC(PC pc)
1161: {
1163:   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
1164:   Mat_IS*        matis;
1165:   MatNullSpace   nearnullspace;
1166:   PetscInt       nrows,ncols;
1167:   PetscBool      computetopography,computesolvers,computesubschurs;
1168:   PetscBool      computeconstraintsmatrix;
1169:   PetscBool      new_nearnullspace_provided,ismatis;

1172:   PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);
1173:   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
1174:   MatGetSize(pc->pmat,&nrows,&ncols);
1175:   if (nrows != ncols) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
1176:   matis = (Mat_IS*)pc->pmat->data;
1177:   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1178:   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1179:      Also, BDDC directly build the Dirichlet problem */
1180:   /* split work */
1181:   if (pc->setupcalled) {
1182:     if (pc->flag == SAME_NONZERO_PATTERN) {
1183:       computetopography = PETSC_FALSE;
1184:       computesolvers = PETSC_TRUE;
1185:     } else { /* DIFFERENT_NONZERO_PATTERN */
1186:       computetopography = PETSC_TRUE;
1187:       computesolvers = PETSC_TRUE;
1188:     }
1189:   } else {
1190:     computetopography = PETSC_TRUE;
1191:     computesolvers = PETSC_TRUE;
1192:   }
1193:   if (pcbddc->recompute_topography) {
1194:     computetopography = PETSC_TRUE;
1195:   }
1196:   computeconstraintsmatrix = PETSC_FALSE;
1197:   if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling");
1198:   pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling);
1199:   if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;

1201:   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
1202:   if (pcbddc->faster_deluxe && pcbddc->adaptive_selection && pcbddc->use_change_of_basis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute faster deluxe if adaptivity and change of basis are both requested. Rerun with -pc_bddc_deluxe_faster false");

1204:   /* Get stdout for dbg */
1205:   if (pcbddc->dbg_flag) {
1206:     if (!pcbddc->dbg_viewer) {
1207:       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1208:       PetscViewerASCIIPushSynchronized(pcbddc->dbg_viewer);
1209:     }
1210:     PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);
1211:   }

1213:   if (pcbddc->user_ChangeOfBasisMatrix) {
1214:     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1215:     pcbddc->use_change_of_basis = PETSC_FALSE;
1216:     PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);
1217:   } else {
1218:     MatDestroy(&pcbddc->local_mat);
1219:     PetscObjectReference((PetscObject)matis->A);
1220:     pcbddc->local_mat = matis->A;
1221:   }

1223:   /* workaround for reals */
1224: #if !defined(PETSC_USE_COMPLEX)
1225:   if (matis->A->symmetric_set) {
1226:     MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);
1227:   }
1228: #endif

1230:   /* Set up all the "iterative substructuring" common block without computing solvers */
1231:   {
1232:     Mat temp_mat;

1234:     temp_mat = matis->A;
1235:     matis->A = pcbddc->local_mat;
1236:     PCISSetUp(pc,PETSC_FALSE);
1237:     pcbddc->local_mat = matis->A;
1238:     matis->A = temp_mat;
1239:   }

1241:   /* Analyze interface and setup sub_schurs data */
1242:   if (computetopography) {
1243:     PCBDDCAnalyzeInterface(pc);
1244:     computeconstraintsmatrix = PETSC_TRUE;
1245:   }

1247:   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1248:   if (computesolvers) {
1249:     PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs;

1251:     if (computesubschurs && computetopography) {
1252:       PCBDDCInitSubSchurs(pc);
1253:     }
1254:     if (sub_schurs->use_mumps) {
1255:       if (computesubschurs) {
1256:         PCBDDCSetUpSubSchurs(pc);
1257:       }
1258:       PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);
1259:     } else {
1260:       PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);
1261:       if (computesubschurs) {
1262:         PCBDDCSetUpSubSchurs(pc);
1263:       }
1264:     }
1265:     if (pcbddc->adaptive_selection) {
1266:       PCBDDCAdaptiveSelection(pc);
1267:       computeconstraintsmatrix = PETSC_TRUE;
1268:     }
1269:   }

1271:   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1272:   new_nearnullspace_provided = PETSC_FALSE;
1273:   MatGetNearNullSpace(pc->pmat,&nearnullspace);
1274:   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1275:     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1276:       new_nearnullspace_provided = PETSC_TRUE;
1277:     } else {
1278:       /* determine if the two nullspaces are different (should be lightweight) */
1279:       if (nearnullspace != pcbddc->onearnullspace) {
1280:         new_nearnullspace_provided = PETSC_TRUE;
1281:       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1282:         PetscInt         i;
1283:         const Vec        *nearnullvecs;
1284:         PetscObjectState state;
1285:         PetscInt         nnsp_size;
1286:         MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);
1287:         for (i=0;i<nnsp_size;i++) {
1288:           PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);
1289:           if (pcbddc->onearnullvecs_state[i] != state) {
1290:             new_nearnullspace_provided = PETSC_TRUE;
1291:             break;
1292:           }
1293:         }
1294:       }
1295:     }
1296:   } else {
1297:     if (!nearnullspace) { /* both nearnullspaces are null */
1298:       new_nearnullspace_provided = PETSC_FALSE;
1299:     } else { /* nearnullspace attached later */
1300:       new_nearnullspace_provided = PETSC_TRUE;
1301:     }
1302:   }

1304:   /* Setup constraints and related work vectors */
1305:   /* reset primal space flags */
1306:   pcbddc->new_primal_space = PETSC_FALSE;
1307:   pcbddc->new_primal_space_local = PETSC_FALSE;
1308:   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1309:     /* It also sets the primal space flags */
1310:     PCBDDCConstraintsSetUp(pc);
1311:     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1312:     PCBDDCSetUpLocalWorkVectors(pc);
1313:   }

1315:   if (computesolvers || pcbddc->new_primal_space) {
1316:     if (pcbddc->use_change_of_basis) {
1317:       PC_IS *pcis = (PC_IS*)(pc->data);

1319:       PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);
1320:       /* get submatrices */
1321:       MatDestroy(&pcis->A_IB);
1322:       MatDestroy(&pcis->A_BI);
1323:       MatDestroy(&pcis->A_BB);
1324:       MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);
1325:       MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);
1326:       MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);
1327:       /* set flag in pcis to not reuse submatrices during PCISCreate */
1328:       pcis->reusesubmatrices = PETSC_FALSE;
1329:     } else if (!pcbddc->user_ChangeOfBasisMatrix) {
1330:       MatDestroy(&pcbddc->local_mat);
1331:       PetscObjectReference((PetscObject)matis->A);
1332:       pcbddc->local_mat = matis->A;
1333:     }
1334:     /* SetUp coarse and local Neumann solvers */
1335:     PCBDDCSetUpSolvers(pc);
1336:     /* SetUp Scaling operator */
1337:     PCBDDCScalingSetUp(pc);
1338:   }

1340:   if (pcbddc->dbg_flag) {
1341:     PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);
1342:   }
1343:   return(0);
1344: }

1346: /* -------------------------------------------------------------------------- */
1347: /*
1348:    PCApply_BDDC - Applies the BDDC operator to a vector.

1350:    Input Parameters:
1351: +  pc - the preconditioner context
1352: -  r - input vector (global)

1354:    Output Parameter:
1355: .  z - output vector (global)

1357:    Application Interface Routine: PCApply()
1358:  */
1361: PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1362: {
1363:   PC_IS             *pcis = (PC_IS*)(pc->data);
1364:   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1365:   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
1366:   PetscErrorCode    ierr;
1367:   const PetscScalar one = 1.0;
1368:   const PetscScalar m_one = -1.0;
1369:   const PetscScalar zero = 0.0;

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

1376:   if (!pcbddc->use_exact_dirichlet_trick) {
1377:     VecCopy(r,z);
1378:     /* First Dirichlet solve */
1379:     VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1380:     VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1381:     /*
1382:       Assembling right hand side for BDDC operator
1383:       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1384:       - pcis->vec1_B the interface part of the global vector z
1385:     */
1386:     if (n_D) {
1387:       KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1388:       VecScale(pcis->vec2_D,m_one);
1389:       if (pcbddc->switch_static) { MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D); }
1390:       MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);
1391:     } else {
1392:       VecSet(pcis->vec1_B,zero);
1393:     }
1394:     VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1395:     VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1396:     PCBDDCScalingRestriction(pc,z,pcis->vec1_B);
1397:   } else {
1398:     if (pcbddc->switch_static) {
1399:       VecSet(pcis->vec1_D,zero);
1400:     }
1401:     PCBDDCScalingRestriction(pc,r,pcis->vec1_B);
1402:   }

1404:   /* Apply interface preconditioner
1405:      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1406:   PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);

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

1411:   /* Second Dirichlet solve and assembling of output */
1412:   VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1413:   VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1414:   if (n_B) {
1415:     MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);
1416:     if (pcbddc->switch_static) { MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D); }
1417:   } else if (pcbddc->switch_static) {
1418:     MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);
1419:   }
1420:   KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);
1421:   if (!pcbddc->use_exact_dirichlet_trick) {
1422:     if (pcbddc->switch_static) {
1423:       VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);
1424:     } else {
1425:       VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);
1426:     }
1427:     VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1428:     VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1429:   } else {
1430:     if (pcbddc->switch_static) {
1431:       VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);
1432:     } else {
1433:       VecScale(pcis->vec4_D,m_one);
1434:     }
1435:     VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
1436:     VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
1437:   }
1438:   return(0);
1439: }

1441: /* -------------------------------------------------------------------------- */
1442: /*
1443:    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.

1445:    Input Parameters:
1446: +  pc - the preconditioner context
1447: -  r - input vector (global)

1449:    Output Parameter:
1450: .  z - output vector (global)

1452:    Application Interface Routine: PCApplyTranspose()
1453:  */
1456: PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
1457: {
1458:   PC_IS             *pcis = (PC_IS*)(pc->data);
1459:   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1460:   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
1461:   PetscErrorCode    ierr;
1462:   const PetscScalar one = 1.0;
1463:   const PetscScalar m_one = -1.0;
1464:   const PetscScalar zero = 0.0;

1467:   if (!pcbddc->use_exact_dirichlet_trick) {
1468:     VecCopy(r,z);
1469:     /* First Dirichlet solve */
1470:     VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1471:     VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1472:     /*
1473:       Assembling right hand side for BDDC operator
1474:       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1475:       - pcis->vec1_B the interface part of the global vector z
1476:     */
1477:     if (n_D) {
1478:       KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1479:       VecScale(pcis->vec2_D,m_one);
1480:       if (pcbddc->switch_static) { MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D); }
1481:       MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);
1482:     } else {
1483:       VecSet(pcis->vec1_B,zero);
1484:     }
1485:     VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1486:     VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1487:     PCBDDCScalingRestriction(pc,z,pcis->vec1_B);
1488:   } else {
1489:     if (pcbddc->switch_static) {
1490:       VecSet(pcis->vec1_D,zero);
1491:     }
1492:     PCBDDCScalingRestriction(pc,r,pcis->vec1_B);
1493:   }

1495:   /* Apply interface preconditioner
1496:      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1497:   PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);

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

1502:   /* Second Dirichlet solve and assembling of output */
1503:   VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1504:   VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1505:   if (n_B) {
1506:     MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);
1507:     if (pcbddc->switch_static) { MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D); }
1508:   } else if (pcbddc->switch_static) {
1509:     MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);
1510:   }
1511:   KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);
1512:   if (!pcbddc->use_exact_dirichlet_trick) {
1513:     if (pcbddc->switch_static) {
1514:       VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);
1515:     } else {
1516:       VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);
1517:     }
1518:     VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1519:     VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1520:   } else {
1521:     if (pcbddc->switch_static) {
1522:       VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);
1523:     } else {
1524:       VecScale(pcis->vec4_D,m_one);
1525:     }
1526:     VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
1527:     VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
1528:   }
1529:   return(0);
1530: }
1531: /* -------------------------------------------------------------------------- */

1535: PetscErrorCode PCDestroy_BDDC(PC pc)
1536: {
1537:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

1541:   /* free data created by PCIS */
1542:   PCISDestroy(pc);
1543:   /* free BDDC custom data  */
1544:   PCBDDCResetCustomization(pc);
1545:   /* destroy objects related to topography */
1546:   PCBDDCResetTopography(pc);
1547:   /* free allocated graph structure */
1548:   PetscFree(pcbddc->mat_graph);
1549:   /* free allocated sub schurs structure */
1550:   PetscFree(pcbddc->sub_schurs);
1551:   /* destroy objects for scaling operator */
1552:   PCBDDCScalingDestroy(pc);
1553:   PetscFree(pcbddc->deluxe_ctx);
1554:   /* free solvers stuff */
1555:   PCBDDCResetSolvers(pc);
1556:   /* free global vectors needed in presolve */
1557:   VecDestroy(&pcbddc->temp_solution);
1558:   VecDestroy(&pcbddc->original_rhs);
1559:   /* free stuff for change of basis hooks */
1560:   if (pcbddc->new_global_mat) {
1561:     PCBDDCChange_ctx change_ctx;
1562:     MatShellGetContext(pcbddc->new_global_mat,&change_ctx);
1563:     MatDestroy(&change_ctx->original_mat);
1564:     MatDestroy(&change_ctx->global_change);
1565:     VecDestroyVecs(2,&change_ctx->work);
1566:     PetscFree(change_ctx);
1567:   }
1568:   MatDestroy(&pcbddc->new_global_mat);
1569:   /* remove functions */
1570:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);
1571:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);
1572:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);
1573:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);
1574:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);
1575:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);
1576:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);
1577:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);
1578:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);
1579:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);
1580:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);
1581:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);
1582:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);
1583:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);
1584:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);
1585:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);
1586:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);
1587:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);
1588:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);
1589:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);
1590:   PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);
1591:   /* Free the private data structure */
1592:   PetscFree(pc->data);
1593:   return(0);
1594: }
1595: /* -------------------------------------------------------------------------- */

1599: static PetscErrorCode PCPreSolveChangeRHS_BDDC(PC pc, PetscBool* change)
1600: {
1602:   *change = PETSC_TRUE;
1603:   return(0);
1604: }

1608: static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1609: {
1610:   FETIDPMat_ctx  mat_ctx;
1611:   Vec            copy_standard_rhs;
1612:   PC_IS*         pcis;
1613:   PC_BDDC*       pcbddc;

1617:   MatShellGetContext(fetidp_mat,&mat_ctx);
1618:   pcis = (PC_IS*)mat_ctx->pc->data;
1619:   pcbddc = (PC_BDDC*)mat_ctx->pc->data;

1621:   /*
1622:      change of basis for physical rhs if needed
1623:      It also changes the rhs in case of dirichlet boundaries
1624:      TODO: better management when FETIDP will have its own class
1625:   */
1626:   VecDuplicate(standard_rhs,&copy_standard_rhs);
1627:   VecCopy(standard_rhs,copy_standard_rhs);
1628:   PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);
1629:   /* store vectors for computation of fetidp final solution */
1630:   VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);
1631:   VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);
1632:   /* scale rhs since it should be unassembled */
1633:   /* TODO use counter scaling? (also below) */
1634:   VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1635:   VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1636:   /* Apply partition of unity */
1637:   VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);
1638:   /* PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B); */
1639:   if (!pcbddc->switch_static) {
1640:     /* compute partially subassembled Schur complement right-hand side */
1641:     KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);
1642:     MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);
1643:     VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);
1644:     VecSet(copy_standard_rhs,0.0);
1645:     VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);
1646:     VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);
1647:     /* PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B); */
1648:     VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1649:     VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1650:     VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);
1651:   }
1652:   VecDestroy(&copy_standard_rhs);
1653:   /* BDDC rhs */
1654:   VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);
1655:   if (pcbddc->switch_static) {
1656:     VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);
1657:   }
1658:   /* apply BDDC */
1659:   PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
1660:   /* Application of B_delta and assembling of rhs for fetidp fluxes */
1661:   VecSet(fetidp_flux_rhs,0.0);
1662:   MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);
1663:   VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
1664:   VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
1665:   return(0);
1666: }

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

1673:    Collective

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

1679:    Output Parameters:
1680: .  fetidp_flux_rhs - the right-hand side for the FETI-DP linear system

1682:    Level: developer

1684:    Notes:

1686: .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
1687: @*/
1688: PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1689: {
1690:   FETIDPMat_ctx  mat_ctx;

1694:   MatShellGetContext(fetidp_mat,&mat_ctx);
1695:   PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));
1696:   return(0);
1697: }
1698: /* -------------------------------------------------------------------------- */

1702: static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1703: {
1704:   FETIDPMat_ctx  mat_ctx;
1705:   PC_IS*         pcis;
1706:   PC_BDDC*       pcbddc;

1710:   MatShellGetContext(fetidp_mat,&mat_ctx);
1711:   pcis = (PC_IS*)mat_ctx->pc->data;
1712:   pcbddc = (PC_BDDC*)mat_ctx->pc->data;

1714:   /* apply B_delta^T */
1715:   VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
1716:   VecScatterEnd(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
1717:   MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);
1718:   /* compute rhs for BDDC application */
1719:   VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);
1720:   if (pcbddc->switch_static) {
1721:     VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);
1722:   }
1723:   /* apply BDDC */
1724:   PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
1725:   /* put values into standard global vector */
1726:   VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1727:   VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1728:   if (!pcbddc->switch_static) {
1729:     /* compute values into the interior if solved for the partially subassembled Schur complement */
1730:     MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);
1731:     VecAYPX(pcis->vec1_D,-1.0,mat_ctx->temp_solution_D);
1732:     KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec1_D);
1733:   }
1734:   VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1735:   VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1736:   /* final change of basis if needed
1737:      Is also sums the dirichlet part removed during RHS assembling */
1738:   PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);
1739:   return(0);
1740: }

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

1747:    Collective

1749:    Input Parameters:
1750: +  fetidp_mat      - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
1751: -  fetidp_flux_sol - the solution of the FETI-DP linear system

1753:    Output Parameters:
1754: .  standard_sol    - the solution defined on the physical domain

1756:    Level: developer

1758:    Notes:

1760: .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
1761: @*/
1762: PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1763: {
1764:   FETIDPMat_ctx  mat_ctx;

1768:   MatShellGetContext(fetidp_mat,&mat_ctx);
1769:   PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));
1770:   return(0);
1771: }
1772: /* -------------------------------------------------------------------------- */

1774: extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1775: extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1776: extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1777: extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1778: extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1779: extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);

1783: static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, PetscBool fully_redundant, Mat *fetidp_mat, PC *fetidp_pc)
1784: {

1786:   FETIDPMat_ctx  fetidpmat_ctx;
1787:   Mat            newmat;
1788:   FETIDPPC_ctx   fetidppc_ctx;
1789:   PC             newpc;
1790:   MPI_Comm       comm;

1794:   PetscObjectGetComm((PetscObject)pc,&comm);
1795:   /* FETIDP linear matrix */
1796:   PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);
1797:   fetidpmat_ctx->fully_redundant = fully_redundant;
1798:   PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);
1799:   MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);
1800:   MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);
1801:   MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);
1802:   MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);
1803:   MatSetUp(newmat);
1804:   /* FETIDP preconditioner */
1805:   PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);
1806:   PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);
1807:   PCCreate(comm,&newpc);
1808:   PCSetType(newpc,PCSHELL);
1809:   PCShellSetContext(newpc,fetidppc_ctx);
1810:   PCShellSetApply(newpc,FETIDPPCApply);
1811:   PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);
1812:   PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);
1813:   PCSetOperators(newpc,newmat,newmat);
1814:   PCSetUp(newpc);
1815:   /* return pointers for objects created */
1816:   *fetidp_mat=newmat;
1817:   *fetidp_pc=newpc;
1818:   return(0);
1819: }

1823: /*@
1824:  PCBDDCCreateFETIDPOperators - Create FETI-DP operators

1826:    Collective

1828:    Input Parameters:
1829: +  pc - the BDDC preconditioning context (setup should have been called before)
1830: -  fully_redundant - true for a fully redundant set of Lagrange multipliers

1832:    Output Parameters:
1833: +  fetidp_mat - shell FETI-DP matrix object
1834: -  fetidp_pc  - shell Dirichlet preconditioner for FETI-DP matrix

1836:    Level: developer

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

1841: .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
1842: @*/
1843: PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, PetscBool fully_redundant, Mat *fetidp_mat, PC *fetidp_pc)
1844: {

1849:   if (pc->setupcalled) {
1850:     PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,PetscBool,Mat*,PC*),(pc,fully_redundant,fetidp_mat,fetidp_pc));
1851:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1852:   return(0);
1853: }
1854: /* -------------------------------------------------------------------------- */
1855: /*MC
1856:    PCBDDC - Balancing Domain Decomposition by Constraints.

1858:    An implementation of the BDDC preconditioner based on

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

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

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

1871:    It also works with unsymmetric and indefinite problems.

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

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

1877:    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()
1878:    Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesLocalIS()

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

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

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

1887:    Adaptive selection of primal constraints [4] is supported for SPD systems with high-contrast in the coefficients if MUMPS is present. Future versions of the code will also consider using MKL_PARDISO or PASTIX.

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

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

1894: .    -pc_bddc_use_vertices <true> - use or not vertices in primal space
1895: .    -pc_bddc_use_edges <true> - use or not edges in primal space
1896: .    -pc_bddc_use_faces <false> - use or not faces in primal space
1897: .    -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
1898: .    -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
1899: .    -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
1900: .    -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
1901: .    -pc_bddc_levels <0> - maximum number of levels for multilevel
1902: .    -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)
1903: .    -pc_bddc_redistribute <0> - size of a subset of processors where the coarse problem will be remapped (the value is ignored if not at the coarsest level)
1904: .    -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
1905: .    -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)
1906: .    -pc_bddc_adaptive_threshold <0.0> - when a value greater than one is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS installed)
1907: -    -pc_bddc_check_level <0> - set verbosity level of debugging output

1909:    Options for Dirichlet, Neumann or coarse solver can be set with
1910: .vb
1911:       -pc_bddc_dirichlet_
1912:       -pc_bddc_neumann_
1913:       -pc_bddc_coarse_
1914: .ve
1915:    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.

1917:    When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
1918: .vb
1919:       -pc_bddc_dirichlet_lN_
1920:       -pc_bddc_neumann_lN_
1921:       -pc_bddc_coarse_lN_
1922: .ve
1923:    Note that level number ranges from the finest (0) to the coarsest (N).
1924:    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.
1925: .vb
1926:      -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
1927: .ve
1928:    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

1930:    Level: intermediate

1932:    Developer notes:

1934:    Contributed by Stefano Zampini

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

1941: PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1942: {
1943:   PetscErrorCode      ierr;
1944:   PC_BDDC             *pcbddc;

1947:   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1948:   PetscNewLog(pc,&pcbddc);
1949:   pc->data  = (void*)pcbddc;

1951:   /* create PCIS data structure */
1952:   PCISCreate(pc);

1954:   /* BDDC customization */
1955:   pcbddc->use_local_adj       = PETSC_TRUE;
1956:   pcbddc->use_vertices        = PETSC_TRUE;
1957:   pcbddc->use_edges           = PETSC_TRUE;
1958:   pcbddc->use_faces           = PETSC_FALSE;
1959:   pcbddc->use_change_of_basis = PETSC_FALSE;
1960:   pcbddc->use_change_on_faces = PETSC_FALSE;
1961:   pcbddc->switch_static       = PETSC_FALSE;
1962:   pcbddc->use_nnsp_true       = PETSC_FALSE;
1963:   pcbddc->use_qr_single       = PETSC_FALSE;
1964:   pcbddc->symmetric_primal    = PETSC_TRUE;
1965:   pcbddc->dbg_flag            = 0;
1966:   /* private */
1967:   pcbddc->local_primal_size          = 0;
1968:   pcbddc->local_primal_size_cc       = 0;
1969:   pcbddc->local_primal_ref_node      = 0;
1970:   pcbddc->local_primal_ref_mult      = 0;
1971:   pcbddc->n_vertices                 = 0;
1972:   pcbddc->primal_indices_local_idxs  = 0;
1973:   pcbddc->recompute_topography       = PETSC_FALSE;
1974:   pcbddc->coarse_size                = -1;
1975:   pcbddc->new_primal_space           = PETSC_FALSE;
1976:   pcbddc->new_primal_space_local     = PETSC_FALSE;
1977:   pcbddc->global_primal_indices      = 0;
1978:   pcbddc->onearnullspace             = 0;
1979:   pcbddc->onearnullvecs_state        = 0;
1980:   pcbddc->user_primal_vertices       = 0;
1981:   pcbddc->temp_solution              = 0;
1982:   pcbddc->original_rhs               = 0;
1983:   pcbddc->local_mat                  = 0;
1984:   pcbddc->ChangeOfBasisMatrix        = 0;
1985:   pcbddc->user_ChangeOfBasisMatrix   = 0;
1986:   pcbddc->new_global_mat             = 0;
1987:   pcbddc->coarse_vec                 = 0;
1988:   pcbddc->coarse_ksp                 = 0;
1989:   pcbddc->coarse_phi_B               = 0;
1990:   pcbddc->coarse_phi_D               = 0;
1991:   pcbddc->coarse_psi_B               = 0;
1992:   pcbddc->coarse_psi_D               = 0;
1993:   pcbddc->vec1_P                     = 0;
1994:   pcbddc->vec1_R                     = 0;
1995:   pcbddc->vec2_R                     = 0;
1996:   pcbddc->local_auxmat1              = 0;
1997:   pcbddc->local_auxmat2              = 0;
1998:   pcbddc->R_to_B                     = 0;
1999:   pcbddc->R_to_D                     = 0;
2000:   pcbddc->ksp_D                      = 0;
2001:   pcbddc->ksp_R                      = 0;
2002:   pcbddc->NeumannBoundaries          = 0;
2003:   pcbddc->NeumannBoundariesLocal     = 0;
2004:   pcbddc->DirichletBoundaries        = 0;
2005:   pcbddc->DirichletBoundariesLocal   = 0;
2006:   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
2007:   pcbddc->n_ISForDofs                = 0;
2008:   pcbddc->n_ISForDofsLocal           = 0;
2009:   pcbddc->ISForDofs                  = 0;
2010:   pcbddc->ISForDofsLocal             = 0;
2011:   pcbddc->ConstraintMatrix           = 0;
2012:   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
2013:   pcbddc->coarse_loc_to_glob         = 0;
2014:   pcbddc->coarsening_ratio           = 8;
2015:   pcbddc->coarse_adj_red             = 0;
2016:   pcbddc->current_level              = 0;
2017:   pcbddc->max_levels                 = 0;
2018:   pcbddc->use_coarse_estimates       = PETSC_FALSE;
2019:   pcbddc->redistribute_coarse        = 0;
2020:   pcbddc->coarse_subassembling       = 0;
2021:   pcbddc->coarse_subassembling_init  = 0;

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

2026:   /* scaling */
2027:   pcbddc->work_scaling          = 0;
2028:   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
2029:   pcbddc->faster_deluxe         = PETSC_FALSE;

2031:   /* create sub schurs structure */
2032:   PCBDDCSubSchursCreate(&pcbddc->sub_schurs);
2033:   pcbddc->sub_schurs_rebuild     = PETSC_FALSE;
2034:   pcbddc->sub_schurs_layers      = -1;
2035:   pcbddc->sub_schurs_use_useradj = PETSC_FALSE;

2037:   pcbddc->computed_rowadj = PETSC_FALSE;

2039:   /* adaptivity */
2040:   pcbddc->adaptive_threshold      = 0.0;
2041:   pcbddc->adaptive_nmax           = 0;
2042:   pcbddc->adaptive_nmin           = 0;

2044:   /* function pointers */
2045:   pc->ops->apply               = PCApply_BDDC;
2046:   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2047:   pc->ops->setup               = PCSetUp_BDDC;
2048:   pc->ops->destroy             = PCDestroy_BDDC;
2049:   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
2050:   pc->ops->view                = PCView_BDDC;
2051:   pc->ops->applyrichardson     = 0;
2052:   pc->ops->applysymmetricleft  = 0;
2053:   pc->ops->applysymmetricright = 0;
2054:   pc->ops->presolve            = PCPreSolve_BDDC;
2055:   pc->ops->postsolve           = PCPostSolve_BDDC;

2057:   /* composing function */
2058:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);
2059:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);
2060:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);
2061:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);
2062:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);
2063:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);
2064:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);
2065:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);
2066:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);
2067:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);
2068:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);
2069:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);
2070:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);
2071:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);
2072:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);
2073:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);
2074:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);
2075:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);
2076:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);
2077:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);
2078:   PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_BDDC);
2079:   return(0);
2080: }