Actual source code: bddc.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: /* TODOLIST

  3:    ConstraintsSetup
  4:    - tolerances for constraints as an option (take care of single precision!)

  6:    Solvers
  7:    - Add support for cholesky for coarse solver (similar to local solvers)
  8:    - Propagate ksp prefixes for solvers to mat objects?
  9:    - Propagate nearnullspace info among levels

 11:    User interface
 12:    - ** DofSplitting and DM attached to pc?

 14:    Debugging output
 15:    - * Better management of verbosity levels of debugging output

 17:    Build
 18:    - make runexe59

 20:    Extra
 21:    - ** GetRid of PCBDDCApplySchur, use MatSchur instead
 22:    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
 23:    - add support for computing h,H and related using coordinates?
 24:    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
 25:    - Better management in PCIS code
 26:    - BDDC with MG framework?

 28:    FETIDP
 29:    - Move FETIDP code to its own classes

 31:    MATIS related operations contained in BDDC code
 32:    - Provide general case for subassembling
 33:    - *** Preallocation routines in MatISGetMPIAXAIJ

 35: */

 37: /* ----------------------------------------------------------------------------------------------------------------------------------------------
 38:    Implementation of BDDC preconditioner based on:
 39:    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
 40:    ---------------------------------------------------------------------------------------------------------------------------------------------- */

 42: #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
 43: #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
 44: #include <petscblaslapack.h>

 46: /* -------------------------------------------------------------------------- */
 49: PetscErrorCode PCSetFromOptions_BDDC(PC pc)
 50: {
 51:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

 55:   PetscOptionsHead("BDDC options");
 56:   /* Verbose debugging */
 57:   PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);
 58:   /* Primal space cumstomization */
 59:   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);
 60:   PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);
 61:   PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);
 62:   PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);
 63:   /* Change of basis */
 64:   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);
 65:   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);
 66:   if (!pcbddc->use_change_of_basis) {
 67:     pcbddc->use_change_on_faces = PETSC_FALSE;
 68:   }
 69:   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
 70:   PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);
 71:   PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);
 72:   PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);
 73:   PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);
 74:   PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);
 75:   PetscOptionsTail();
 76:   return(0);
 77: }
 78: /* -------------------------------------------------------------------------- */
 81: static PetscErrorCode PCBDDCSetChangeOfBasisLocalMat_BDDC(PC pc, Mat change)
 82: {
 83:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

 87:   MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);
 88:   PetscObjectReference((PetscObject)change);
 89:   pcbddc->user_ChangeOfBasisMatrix = change;
 90:   return(0);
 91: }
 94: /*@
 95:  PCBDDCSetChangeOfBasisLocalMat - Set user defined change of basis for local boundary dofs

 97:    Collective on PC

 99:    Input Parameters:
100: +  pc - the preconditioning context
101: -  change - the local change of basis matrix, either in local (internal + boundary) or in local boundary numbering

103:    Level: intermediate

105:    Notes:

107: .seealso: PCBDDC
108: @*/
109: PetscErrorCode PCBDDCSetChangeOfBasisLocalMat(PC pc, Mat change)
110: {

116:   PetscTryMethod(pc,"PCBDDCSetChangeOfBasisLocalMat_C",(PC,Mat),(pc,change));
117:   return(0);
118: }
119: /* -------------------------------------------------------------------------- */
122: static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
123: {
124:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

128:   ISDestroy(&pcbddc->user_primal_vertices);
129:   PetscObjectReference((PetscObject)PrimalVertices);
130:   pcbddc->user_primal_vertices = PrimalVertices;
131:   return(0);
132: }
135: /*@
136:  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC

138:    Not collective

140:    Input Parameters:
141: +  pc - the preconditioning context
142: -  PrimalVertices - index set of primal vertices in local numbering

144:    Level: intermediate

146:    Notes:

148: .seealso: PCBDDC
149: @*/
150: PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
151: {

157:   PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));
158:   return(0);
159: }
160: /* -------------------------------------------------------------------------- */
163: static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
164: {
165:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

168:   pcbddc->coarsening_ratio = k;
169:   return(0);
170: }

174: /*@
175:  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel

177:    Logically collective on PC

179:    Input Parameters:
180: +  pc - the preconditioning context
181: -  k - coarsening ratio (H/h at the coarser level)

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

185:    Level: intermediate

187:    Notes:

189: .seealso: PCBDDC
190: @*/
191: PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
192: {

198:   PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));
199:   return(0);
200: }

202: /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
205: static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
206: {
207:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

210:   pcbddc->use_exact_dirichlet_trick = flg;
211:   return(0);
212: }

216: PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
217: {

223:   PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));
224:   return(0);
225: }

229: static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
230: {
231:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

234:   pcbddc->current_level = level;
235:   return(0);
236: }

240: PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
241: {

247:   PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));
248:   return(0);
249: }

253: static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
254: {
255:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

258:   pcbddc->max_levels = levels;
259:   return(0);
260: }

264: /*@
265:  PCBDDCSetLevels - Sets the maximum number of levels for multilevel

267:    Logically collective on PC

269:    Input Parameters:
270: +  pc - the preconditioning context
271: -  levels - the maximum number of levels (max 9)

273:    Default value is 0, i.e. traditional one-level BDDC

275:    Level: intermediate

277:    Notes:

279: .seealso: PCBDDC
280: @*/
281: PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
282: {

288:   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
289:   PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));
290:   return(0);
291: }
292: /* -------------------------------------------------------------------------- */

296: static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
297: {
298:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

302:   PetscObjectReference((PetscObject)NullSpace);
303:   MatNullSpaceDestroy(&pcbddc->NullSpace);
304:   pcbddc->NullSpace = NullSpace;
305:   return(0);
306: }

310: /*@
311:  PCBDDCSetNullSpace - Set nullspace for BDDC operator

313:    Logically collective on PC and MatNullSpace

315:    Input Parameters:
316: +  pc - the preconditioning context
317: -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)

319:    Level: intermediate

321:    Notes:

323: .seealso: PCBDDC
324: @*/
325: PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
326: {

333:   PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));
334:   return(0);
335: }
336: /* -------------------------------------------------------------------------- */

340: static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
341: {
342:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

346:   /* last user setting takes precendence -> destroy any other customization */
347:   ISDestroy(&pcbddc->DirichletBoundariesLocal);
348:   ISDestroy(&pcbddc->DirichletBoundaries);
349:   PetscObjectReference((PetscObject)DirichletBoundaries);
350:   pcbddc->DirichletBoundaries = DirichletBoundaries;
351:   pcbddc->recompute_topography = PETSC_TRUE;
352:   return(0);
353: }

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

360:    Collective

362:    Input Parameters:
363: +  pc - the preconditioning context
364: -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries

366:    Level: intermediate

368:    Notes: Any process can list any global node

370: .seealso: PCBDDC
371: @*/
372: PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
373: {

380:   PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));
381:   return(0);
382: }
383: /* -------------------------------------------------------------------------- */

387: static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
388: {
389:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

393:   /* last user setting takes precendence -> destroy any other customization */
394:   ISDestroy(&pcbddc->DirichletBoundariesLocal);
395:   ISDestroy(&pcbddc->DirichletBoundaries);
396:   PetscObjectReference((PetscObject)DirichletBoundaries);
397:   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
398:   pcbddc->recompute_topography = PETSC_TRUE;
399:   return(0);
400: }

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

407:    Collective

409:    Input Parameters:
410: +  pc - the preconditioning context
411: -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)

413:    Level: intermediate

415:    Notes:

417: .seealso: PCBDDC
418: @*/
419: PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
420: {

427:   PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));
428:   return(0);
429: }
430: /* -------------------------------------------------------------------------- */

434: static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
435: {
436:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

440:   /* last user setting takes precendence -> destroy any other customization */
441:   ISDestroy(&pcbddc->NeumannBoundariesLocal);
442:   ISDestroy(&pcbddc->NeumannBoundaries);
443:   PetscObjectReference((PetscObject)NeumannBoundaries);
444:   pcbddc->NeumannBoundaries = NeumannBoundaries;
445:   pcbddc->recompute_topography = PETSC_TRUE;
446:   return(0);
447: }

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

454:    Collective

456:    Input Parameters:
457: +  pc - the preconditioning context
458: -  NeumannBoundaries - parallel IS defining the Neumann boundaries

460:    Level: intermediate

462:    Notes: Any process can list any global node

464: .seealso: PCBDDC
465: @*/
466: PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
467: {

474:   PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));
475:   return(0);
476: }
477: /* -------------------------------------------------------------------------- */

481: static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
482: {
483:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

487:   /* last user setting takes precendence -> destroy any other customization */
488:   ISDestroy(&pcbddc->NeumannBoundariesLocal);
489:   ISDestroy(&pcbddc->NeumannBoundaries);
490:   PetscObjectReference((PetscObject)NeumannBoundaries);
491:   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
492:   pcbddc->recompute_topography = PETSC_TRUE;
493:   return(0);
494: }

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

501:    Collective

503:    Input Parameters:
504: +  pc - the preconditioning context
505: -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)

507:    Level: intermediate

509:    Notes:

511: .seealso: PCBDDC
512: @*/
513: PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
514: {

521:   PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));
522:   return(0);
523: }
524: /* -------------------------------------------------------------------------- */

528: static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
529: {
530:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

533:   *DirichletBoundaries = pcbddc->DirichletBoundaries;
534:   return(0);
535: }

539: /*@
540:  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries

542:    Collective

544:    Input Parameters:
545: .  pc - the preconditioning context

547:    Output Parameters:
548: .  DirichletBoundaries - index set defining the Dirichlet boundaries

550:    Level: intermediate

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

554: .seealso: PCBDDC
555: @*/
556: PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
557: {

562:   PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));
563:   return(0);
564: }
565: /* -------------------------------------------------------------------------- */

569: static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
570: {
571:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

574:   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
575:   return(0);
576: }

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

583:    Collective

585:    Input Parameters:
586: .  pc - the preconditioning context

588:    Output Parameters:
589: .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries

591:    Level: intermediate

593:    Notes:

595: .seealso: PCBDDC
596: @*/
597: PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
598: {

603:   PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));
604:   return(0);
605: }
606: /* -------------------------------------------------------------------------- */

610: static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
611: {
612:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

615:   *NeumannBoundaries = pcbddc->NeumannBoundaries;
616:   return(0);
617: }

621: /*@
622:  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries

624:    Collective

626:    Input Parameters:
627: .  pc - the preconditioning context

629:    Output Parameters:
630: .  NeumannBoundaries - index set defining the Neumann boundaries

632:    Level: intermediate

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

636: .seealso: PCBDDC
637: @*/
638: PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
639: {

644:   PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));
645:   return(0);
646: }
647: /* -------------------------------------------------------------------------- */

651: static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
652: {
653:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;

656:   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
657:   return(0);
658: }

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

665:    Collective

667:    Input Parameters:
668: .  pc - the preconditioning context

670:    Output Parameters:
671: .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries

673:    Level: intermediate

675:    Notes:

677: .seealso: PCBDDC
678: @*/
679: PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
680: {

685:   PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));
686:   return(0);
687: }
688: /* -------------------------------------------------------------------------- */

692: static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
693: {
694:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
695:   PCBDDCGraph    mat_graph = pcbddc->mat_graph;

699:   /* free old CSR */
700:   PCBDDCGraphResetCSR(mat_graph);
701:   /* TODO: PCBDDCGraphSetAdjacency */
702:   /* get CSR into graph structure */
703:   if (copymode == PETSC_COPY_VALUES) {
704:     PetscMalloc1((nvtxs+1),&mat_graph->xadj);
705:     PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);
706:     PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));
707:     PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));
708:   } else if (copymode == PETSC_OWN_POINTER) {
709:     mat_graph->xadj = (PetscInt*)xadj;
710:     mat_graph->adjncy = (PetscInt*)adjncy;
711:   }
712:   mat_graph->nvtxs_csr = nvtxs;
713:   return(0);
714: }

718: /*@
719:  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix

721:    Not collective

723:    Input Parameters:
724: +  pc - the preconditioning context
725: .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
726: .  xadj, adjncy - the CSR graph
727: -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.

729:    Level: intermediate

731:    Notes:

733: .seealso: PCBDDC,PetscCopyMode
734: @*/
735: PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
736: {
737:   void (*f)(void) = 0;

744:   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
745:     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
746:   }
747:   PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));
748:   /* free arrays if PCBDDC is not the PC type */
749:   PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);
750:   if (!f && copymode == PETSC_OWN_POINTER) {
751:     PetscFree(xadj);
752:     PetscFree(adjncy);
753:   }
754:   return(0);
755: }
756: /* -------------------------------------------------------------------------- */

760: static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
761: {
762:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
763:   PetscInt i;

767:   /* Destroy ISes if they were already set */
768:   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
769:     ISDestroy(&pcbddc->ISForDofsLocal[i]);
770:   }
771:   PetscFree(pcbddc->ISForDofsLocal);
772:   /* last user setting takes precendence -> destroy any other customization */
773:   for (i=0;i<pcbddc->n_ISForDofs;i++) {
774:     ISDestroy(&pcbddc->ISForDofs[i]);
775:   }
776:   PetscFree(pcbddc->ISForDofs);
777:   pcbddc->n_ISForDofs = 0;
778:   /* allocate space then set */
779:   PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);
780:   for (i=0;i<n_is;i++) {
781:     PetscObjectReference((PetscObject)ISForDofs[i]);
782:     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
783:   }
784:   pcbddc->n_ISForDofsLocal=n_is;
785:   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
786:   pcbddc->recompute_topography = PETSC_TRUE;
787:   return(0);
788: }

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

795:    Collective

797:    Input Parameters:
798: +  pc - the preconditioning context
799: -  n_is - number of index sets defining the fields
800: .  ISForDofs - array of IS describing the fields in local ordering

802:    Level: intermediate

804:    Notes: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to a different field.

806: .seealso: PCBDDC
807: @*/
808: PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
809: {
810:   PetscInt       i;

816:   for (i=0;i<n_is;i++) {
819:   }
820:   PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));
821:   return(0);
822: }
823: /* -------------------------------------------------------------------------- */

827: static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
828: {
829:   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
830:   PetscInt i;

834:   /* Destroy ISes if they were already set */
835:   for (i=0;i<pcbddc->n_ISForDofs;i++) {
836:     ISDestroy(&pcbddc->ISForDofs[i]);
837:   }
838:   PetscFree(pcbddc->ISForDofs);
839:   /* last user setting takes precendence -> destroy any other customization */
840:   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
841:     ISDestroy(&pcbddc->ISForDofsLocal[i]);
842:   }
843:   PetscFree(pcbddc->ISForDofsLocal);
844:   pcbddc->n_ISForDofsLocal = 0;
845:   /* allocate space then set */
846:   PetscMalloc1(n_is,&pcbddc->ISForDofs);
847:   for (i=0;i<n_is;i++) {
848:     PetscObjectReference((PetscObject)ISForDofs[i]);
849:     pcbddc->ISForDofs[i]=ISForDofs[i];
850:   }
851:   pcbddc->n_ISForDofs=n_is;
852:   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
853:   pcbddc->recompute_topography = PETSC_TRUE;
854:   return(0);
855: }

859: /*@
860:  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix

862:    Collective

864:    Input Parameters:
865: +  pc - the preconditioning context
866: -  n_is - number of index sets defining the fields
867: .  ISForDofs - array of IS describing the fields in global ordering

869:    Level: intermediate

871:    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.

873: .seealso: PCBDDC
874: @*/
875: PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
876: {
877:   PetscInt       i;

883:   for (i=0;i<n_is;i++) {
886:   }
887:   PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));
888:   return(0);
889: }
890: /* -------------------------------------------------------------------------- */
893: /* -------------------------------------------------------------------------- */
894: /*
895:    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
896:                      guess if a transformation of basis approach has been selected.

898:    Input Parameter:
899: +  pc - the preconditioner contex

901:    Application Interface Routine: PCPreSolve()

903:    Notes:
904:    The interface routine PCPreSolve() is not usually called directly by
905:    the user, but instead is called by KSPSolve().
906: */
907: static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
908: {
910:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
911:   PC_IS          *pcis = (PC_IS*)(pc->data);
912:   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
913:   Mat            temp_mat;
914:   IS             dirIS;
915:   Vec            used_vec;
916:   PetscBool      guess_nonzero;

919:   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
920:   if (ksp) {
921:     PetscBool iscg;
922:     PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);
923:     if (!iscg) {
924:       PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);
925:     }
926:   }
927:   /* Creates parallel work vectors used in presolve */
928:   if (!pcbddc->original_rhs) {
929:     VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);
930:   }
931:   if (!pcbddc->temp_solution) {
932:     VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);
933:   }
934:   if (x) {
935:     PetscObjectReference((PetscObject)x);
936:     used_vec = x;
937:   } else {
938:     PetscObjectReference((PetscObject)pcbddc->temp_solution);
939:     used_vec = pcbddc->temp_solution;
940:     VecSet(used_vec,0.0);
941:   }
942:   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
943:   if (ksp) {
944:     KSPGetInitialGuessNonzero(ksp,&guess_nonzero);
945:     if (!guess_nonzero) {
946:       VecSet(used_vec,0.0);
947:     }
948:   }

950:   /* store the original rhs */
951:   VecCopy(rhs,pcbddc->original_rhs);

953:   /* Take into account zeroed rows -> change rhs and store solution removed */
954:   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
955:   PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);
956:   if (rhs && dirIS) {
957:     PetscInt    dirsize,i,*is_indices;
958:     PetscScalar *array_x,*array_diagonal;

960:     MatGetDiagonal(pc->pmat,pcis->vec1_global);
961:     VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);
962:     VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);
963:     VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);
964:     VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
965:     VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
966:     ISGetLocalSize(dirIS,&dirsize);
967:     VecGetArray(pcis->vec1_N,&array_x);
968:     VecGetArray(pcis->vec2_N,&array_diagonal);
969:     ISGetIndices(dirIS,(const PetscInt**)&is_indices);
970:     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
971:     ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);
972:     VecRestoreArray(pcis->vec2_N,&array_diagonal);
973:     VecRestoreArray(pcis->vec1_N,&array_x);
974:     VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);
975:     VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);

977:     /* remove the computed solution from the rhs */
978:     VecScale(used_vec,-1.0);
979:     MatMultAdd(pc->pmat,used_vec,rhs,rhs);
980:     VecScale(used_vec,-1.0);
981:   }

983:   /* store partially computed solution and set initial guess */
984:   if (x) {
985:     VecCopy(used_vec,pcbddc->temp_solution);
986:     VecSet(used_vec,0.0);
987:     if (pcbddc->use_exact_dirichlet_trick) {
988:       VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
989:       VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
990:       KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
991:       VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);
992:       VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);
993:       if (ksp) {
994:         KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
995:       }
996:     }
997:   }

999:   /* prepare MatMult and rhs for solver */
1000:   if (pcbddc->ChangeOfBasisMatrix) {
1001:     /* swap pointers for local matrices */
1002:     temp_mat = matis->A;
1003:     matis->A = pcbddc->local_mat;
1004:     pcbddc->local_mat = temp_mat;
1005:     if (rhs) {
1006:       /* Get local rhs and apply transformation of basis */
1007:       VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1008:       VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1009:       /* from original basis to modified basis */
1010:       MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);
1011:       /* put back modified values into the global vec using INSERT_VALUES copy mode */
1012:       VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);
1013:       VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);
1014:     }
1015:   }

1017:   /* remove nullspace if present */
1018:   if (ksp && pcbddc->NullSpace) {
1019:     MatNullSpaceRemove(pcbddc->NullSpace,used_vec);
1020:     MatNullSpaceRemove(pcbddc->NullSpace,rhs);
1021:   }
1022:   VecDestroy(&used_vec);
1023:   return(0);
1024: }
1025: /* -------------------------------------------------------------------------- */
1028: /* -------------------------------------------------------------------------- */
1029: /*
1030:    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1031:                      approach has been selected. Also, restores rhs to its original state.

1033:    Input Parameter:
1034: +  pc - the preconditioner contex

1036:    Application Interface Routine: PCPostSolve()

1038:    Notes:
1039:    The interface routine PCPostSolve() is not usually called directly by
1040:    the user, but instead is called by KSPSolve().
1041: */
1042: static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1043: {
1045:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1046:   PC_IS          *pcis   = (PC_IS*)(pc->data);
1047:   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
1048:   Mat            temp_mat;

1051:   if (pcbddc->ChangeOfBasisMatrix) {
1052:     /* swap pointers for local matrices */
1053:     temp_mat = matis->A;
1054:     matis->A = pcbddc->local_mat;
1055:     pcbddc->local_mat = temp_mat;
1056:   }
1057:   if (pcbddc->ChangeOfBasisMatrix && x) {
1058:     /* Get Local boundary and apply transformation of basis to solution vector */
1059:     VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1060:     VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1061:     /* from modified basis to original basis */
1062:     MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);
1063:     /* put back modified values into the global vec using INSERT_VALUES copy mode */
1064:     VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);
1065:     VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);
1066:   }
1067:   /* add solution removed in presolve */
1068:   if (x) {
1069:     VecAXPY(x,1.0,pcbddc->temp_solution);
1070:   }
1071:   /* restore rhs to its original state */
1072:   if (rhs) {
1073:     VecCopy(pcbddc->original_rhs,rhs);
1074:   }
1075:   return(0);
1076: }
1077: /* -------------------------------------------------------------------------- */
1080: /* -------------------------------------------------------------------------- */
1081: /*
1082:    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1083:                   by setting data structures and options.

1085:    Input Parameter:
1086: +  pc - the preconditioner context

1088:    Application Interface Routine: PCSetUp()

1090:    Notes:
1091:    The interface routine PCSetUp() is not usually called directly by
1092:    the user, but instead is called by PCApply() if necessary.
1093: */
1094: PetscErrorCode PCSetUp_BDDC(PC pc)
1095: {
1096:   PetscErrorCode   ierr;
1097:   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
1098:   MatNullSpace     nearnullspace;
1099:   PetscBool        computeis,computetopography,computesolvers;
1100:   PetscBool        new_nearnullspace_provided;

1103:   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1104:   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1105:      Also, BDDC directly build the Dirichlet problem */
1106:   /* split work */
1107:   if (pc->setupcalled) {
1108:     computeis = PETSC_FALSE;
1109:     if (pc->flag == SAME_NONZERO_PATTERN) {
1110:       computetopography = PETSC_FALSE;
1111:       computesolvers = PETSC_TRUE;
1112:     } else { /* DIFFERENT_NONZERO_PATTERN */
1113:       computetopography = PETSC_TRUE;
1114:       computesolvers = PETSC_TRUE;
1115:     }
1116:   } else {
1117:     computeis = PETSC_TRUE;
1118:     computetopography = PETSC_TRUE;
1119:     computesolvers = PETSC_TRUE;
1120:   }
1121:   if (pcbddc->recompute_topography) {
1122:     computetopography = PETSC_TRUE;
1123:   }

1125:   /* Get stdout for dbg */
1126:   if (pcbddc->dbg_flag) {
1127:     if (!pcbddc->dbg_viewer) {
1128:       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1129:       PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);
1130:     }
1131:     PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);
1132:   }

1134:   /* Set up all the "iterative substructuring" common block without computing solvers */
1135:   if (computeis) {
1136:     /* HACK INTO PCIS */
1137:     PC_IS* pcis = (PC_IS*)pc->data;
1138:     pcis->computesolvers = PETSC_FALSE;
1139:     PCISSetUp(pc);
1140:     ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);
1141:   }

1143:   /* check user defined change of basis (if any) */
1144:   if (pcbddc->user_ChangeOfBasisMatrix) {
1145:     PC_IS* pcis= (PC_IS*)pc->data;
1146:     PetscInt n,m;
1147:     MatGetSize(pcbddc->user_ChangeOfBasisMatrix,&n,&m);
1148:     if (n != m) {
1149:       SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Change of basis matrix should be square %d != %d\n",n,m);
1150:     } else if (n != pcis->n_B && n != pcis->n) {
1151:       SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Size of change of basis matrix %d differs from boundary size %d and local size %d\n",n,pcis->n_B,pcis->n);
1152:     }
1153:     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1154:     pcbddc->use_change_of_basis = PETSC_FALSE;
1155:   }
1156:   /* Analyze interface */
1157:   if (computetopography) {
1158:     PCBDDCAnalyzeInterface(pc);
1159:   }

1161:   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1162:   new_nearnullspace_provided = PETSC_FALSE;
1163:   MatGetNearNullSpace(pc->pmat,&nearnullspace);
1164:   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1165:     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1166:       new_nearnullspace_provided = PETSC_TRUE;
1167:     } else {
1168:       /* determine if the two nullspaces are different (should be lightweight) */
1169:       if (nearnullspace != pcbddc->onearnullspace) {
1170:         new_nearnullspace_provided = PETSC_TRUE;
1171:       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1172:         PetscInt         i;
1173:         const Vec        *nearnullvecs;
1174:         PetscObjectState state;
1175:         PetscInt         nnsp_size;
1176:         MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);
1177:         for (i=0;i<nnsp_size;i++) {
1178:           PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);
1179:           if (pcbddc->onearnullvecs_state[i] != state) {
1180:             new_nearnullspace_provided = PETSC_TRUE;
1181:             break;
1182:           }
1183:         }
1184:       }
1185:     }
1186:   } else {
1187:     if (!nearnullspace) { /* both nearnullspaces are null */
1188:       new_nearnullspace_provided = PETSC_FALSE;
1189:     } else { /* nearnullspace attached later */
1190:       new_nearnullspace_provided = PETSC_TRUE;
1191:     }
1192:   }

1194:   /* Setup constraints and related work vectors */
1195:   /* reset primal space flags */
1196:   pcbddc->new_primal_space = PETSC_FALSE;
1197:   pcbddc->new_primal_space_local = PETSC_FALSE;
1198:   if (computetopography || new_nearnullspace_provided) {
1199:     /* It also sets the primal space flags */
1200:     PCBDDCConstraintsSetUp(pc);
1201:     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1202:     PCBDDCSetUpLocalWorkVectors(pc);
1203:   }

1205:   if (computesolvers || pcbddc->new_primal_space) {
1206:     /* reset data */
1207:     PCBDDCScalingDestroy(pc);
1208:     /* Create coarse and local stuffs */
1209:     PCBDDCSetUpSolvers(pc);
1210:     PCBDDCScalingSetUp(pc);
1211:   }

1213:   if (pcbddc->dbg_flag) {
1214:     PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);
1215:   }
1216:   return(0);
1217: }

1219: /* -------------------------------------------------------------------------- */
1220: /*
1221:    PCApply_BDDC - Applies the BDDC operator to a vector.

1223:    Input Parameters:
1224: .  pc - the preconditioner context
1225: .  r - input vector (global)

1227:    Output Parameter:
1228: .  z - output vector (global)

1230:    Application Interface Routine: PCApply()
1231:  */
1234: PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1235: {
1236:   PC_IS             *pcis = (PC_IS*)(pc->data);
1237:   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1238:   PetscErrorCode    ierr;
1239:   const PetscScalar one = 1.0;
1240:   const PetscScalar m_one = -1.0;
1241:   const PetscScalar zero = 0.0;

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

1248:   if (!pcbddc->use_exact_dirichlet_trick) {
1249:     /* First Dirichlet solve */
1250:     VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1251:     VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1252:     KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1253:     /*
1254:       Assembling right hand side for BDDC operator
1255:       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1256:       - pcis->vec1_B the interface part of the global vector z
1257:     */
1258:     VecScale(pcis->vec2_D,m_one);
1259:     MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);
1260:     if (pcbddc->switch_static) { MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D); }
1261:     VecScale(pcis->vec2_D,m_one);
1262:     VecCopy(r,z);
1263:     VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1264:     VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1265:     PCBDDCScalingRestriction(pc,z,pcis->vec1_B);
1266:   } else {
1267:     VecSet(pcis->vec1_D,zero);
1268:     VecSet(pcis->vec2_D,zero);
1269:     PCBDDCScalingRestriction(pc,r,pcis->vec1_B);
1270:   }

1272:   /* Apply interface preconditioner
1273:      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1274:   PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);

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

1279:   /* Second Dirichlet solve and assembling of output */
1280:   VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1281:   VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1282:   MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);
1283:   if (pcbddc->switch_static) { MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D); }
1284:   KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);
1285:   VecScale(pcis->vec4_D,m_one);
1286:   if (pcbddc->switch_static) { VecAXPY(pcis->vec4_D,one,pcis->vec1_D); }
1287:   VecAXPY(pcis->vec2_D,one,pcis->vec4_D);
1288:   VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1289:   VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1290:   return(0);
1291: }

1293: /* -------------------------------------------------------------------------- */
1294: /*
1295:    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.

1297:    Input Parameters:
1298: .  pc - the preconditioner context
1299: .  r - input vector (global)

1301:    Output Parameter:
1302: .  z - output vector (global)

1304:    Application Interface Routine: PCApplyTranspose()
1305:  */
1308: PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
1309: {
1310:   PC_IS             *pcis = (PC_IS*)(pc->data);
1311:   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1312:   PetscErrorCode    ierr;
1313:   const PetscScalar one = 1.0;
1314:   const PetscScalar m_one = -1.0;
1315:   const PetscScalar zero = 0.0;

1318:   if (!pcbddc->use_exact_dirichlet_trick) {
1319:     /* First Dirichlet solve */
1320:     VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1321:     VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1322:     KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1323:     /*
1324:       Assembling right hand side for BDDC operator
1325:       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1326:       - pcis->vec1_B the interface part of the global vector z
1327:     */
1328:     VecScale(pcis->vec2_D,m_one);
1329:     MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);
1330:     if (pcbddc->switch_static) { MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D); }
1331:     VecScale(pcis->vec2_D,m_one);
1332:     VecCopy(r,z);
1333:     VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1334:     VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1335:     PCBDDCScalingRestriction(pc,z,pcis->vec1_B);
1336:   } else {
1337:     VecSet(pcis->vec1_D,zero);
1338:     VecSet(pcis->vec2_D,zero);
1339:     PCBDDCScalingRestriction(pc,r,pcis->vec1_B);
1340:   }

1342:   /* Apply interface preconditioner
1343:      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1344:   PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);

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

1349:   /* Second Dirichlet solve and assembling of output */
1350:   VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1351:   VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1352:   MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);
1353:   if (pcbddc->switch_static) { MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D); }
1354:   KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);
1355:   VecScale(pcis->vec4_D,m_one);
1356:   if (pcbddc->switch_static) { VecAXPY(pcis->vec4_D,one,pcis->vec1_D); }
1357:   VecAXPY(pcis->vec2_D,one,pcis->vec4_D);
1358:   VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1359:   VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1360:   return(0);
1361: }
1362: /* -------------------------------------------------------------------------- */

1366: PetscErrorCode PCDestroy_BDDC(PC pc)
1367: {
1368:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;

1372:   /* free data created by PCIS */
1373:   PCISDestroy(pc);
1374:   /* free BDDC custom data  */
1375:   PCBDDCResetCustomization(pc);
1376:   /* destroy objects related to topography */
1377:   PCBDDCResetTopography(pc);
1378:   /* free allocated graph structure */
1379:   PetscFree(pcbddc->mat_graph);
1380:   /* free data for scaling operator */
1381:   PCBDDCScalingDestroy(pc);
1382:   /* free solvers stuff */
1383:   PCBDDCResetSolvers(pc);
1384:   /* free global vectors needed in presolve */
1385:   VecDestroy(&pcbddc->temp_solution);
1386:   VecDestroy(&pcbddc->original_rhs);
1387:   ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);
1388:   /* remove functions */
1389:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisLocalMat_C",NULL);
1390:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);
1391:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);
1392:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);
1393:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);
1394:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);
1395:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);
1396:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);
1397:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);
1398:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);
1399:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);
1400:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);
1401:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);
1402:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);
1403:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);
1404:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);
1405:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);
1406:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);
1407:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);
1408:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);
1409:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);
1410:   /* Free the private data structure */
1411:   PetscFree(pc->data);
1412:   return(0);
1413: }
1414: /* -------------------------------------------------------------------------- */

1418: static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1419: {
1420:   FETIDPMat_ctx  mat_ctx;
1421:   PC_IS*         pcis;
1422:   PC_BDDC*       pcbddc;

1426:   MatShellGetContext(fetidp_mat,&mat_ctx);
1427:   pcis = (PC_IS*)mat_ctx->pc->data;
1428:   pcbddc = (PC_BDDC*)mat_ctx->pc->data;

1430:   /* change of basis for physical rhs if needed
1431:      It also changes the rhs in case of dirichlet boundaries */
1432:   PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);
1433:   /* store vectors for computation of fetidp final solution */
1434:   VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);
1435:   VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);
1436:   /* scale rhs since it should be unassembled */
1437:   /* TODO use counter scaling? (also below) */
1438:   VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1439:   VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1440:   /* Apply partition of unity */
1441:   VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);
1442:   /* PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B); */
1443:   if (!pcbddc->switch_static) {
1444:     /* compute partially subassembled Schur complement right-hand side */
1445:     KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);
1446:     MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);
1447:     VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);
1448:     VecSet(standard_rhs,0.0);
1449:     VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);
1450:     VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);
1451:     /* PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B); */
1452:     VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1453:     VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1454:     VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);
1455:   }
1456:   /* BDDC rhs */
1457:   VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);
1458:   if (pcbddc->switch_static) {
1459:     VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);
1460:   }
1461:   /* apply BDDC */
1462:   PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
1463:   /* Application of B_delta and assembling of rhs for fetidp fluxes */
1464:   VecSet(fetidp_flux_rhs,0.0);
1465:   MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);
1466:   VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
1467:   VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
1468:   /* restore original rhs */
1469:   VecCopy(pcbddc->original_rhs,standard_rhs);
1470:   return(0);
1471: }

1475: /*@
1476:  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system

1478:    Collective

1480:    Input Parameters:
1481: +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
1482: .  standard_rhs - the right-hand side for your linear system

1484:    Output Parameters:
1485: -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system

1487:    Level: developer

1489:    Notes:

1491: .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1492: @*/
1493: PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1494: {
1495:   FETIDPMat_ctx  mat_ctx;

1499:   MatShellGetContext(fetidp_mat,&mat_ctx);
1500:   PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));
1501:   return(0);
1502: }
1503: /* -------------------------------------------------------------------------- */

1507: static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1508: {
1509:   FETIDPMat_ctx  mat_ctx;
1510:   PC_IS*         pcis;
1511:   PC_BDDC*       pcbddc;

1515:   MatShellGetContext(fetidp_mat,&mat_ctx);
1516:   pcis = (PC_IS*)mat_ctx->pc->data;
1517:   pcbddc = (PC_BDDC*)mat_ctx->pc->data;

1519:   /* apply B_delta^T */
1520:   VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
1521:   VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
1522:   MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);
1523:   /* compute rhs for BDDC application */
1524:   VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);
1525:   if (pcbddc->switch_static) {
1526:     VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);
1527:   }
1528:   /* apply BDDC */
1529:   PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
1530:   /* put values into standard global vector */
1531:   VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1532:   VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1533:   if (!pcbddc->switch_static) {
1534:     /* compute values into the interior if solved for the partially subassembled Schur complement */
1535:     MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);
1536:     VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);
1537:     KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);
1538:   }
1539:   VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1540:   VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1541:   /* final change of basis if needed
1542:      Is also sums the dirichlet part removed during RHS assembling */
1543:   PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);
1544:   return(0);
1545: }

1549: /*@
1550:  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system

1552:    Collective

1554:    Input Parameters:
1555: +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
1556: .  fetidp_flux_sol - the solution of the FETIDP linear system

1558:    Output Parameters:
1559: -  standard_sol      - the solution defined on the physical domain

1561:    Level: developer

1563:    Notes:

1565: .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
1566: @*/
1567: PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1568: {
1569:   FETIDPMat_ctx  mat_ctx;

1573:   MatShellGetContext(fetidp_mat,&mat_ctx);
1574:   PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));
1575:   return(0);
1576: }
1577: /* -------------------------------------------------------------------------- */

1579: extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1580: extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1581: extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1582: extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1583: extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1584: extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);

1588: static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1589: {

1591:   FETIDPMat_ctx  fetidpmat_ctx;
1592:   Mat            newmat;
1593:   FETIDPPC_ctx   fetidppc_ctx;
1594:   PC             newpc;
1595:   MPI_Comm       comm;

1599:   PetscObjectGetComm((PetscObject)pc,&comm);
1600:   /* FETIDP linear matrix */
1601:   PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);
1602:   PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);
1603:   MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);
1604:   MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);
1605:   MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);
1606:   MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);
1607:   MatSetUp(newmat);
1608:   /* FETIDP preconditioner */
1609:   PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);
1610:   PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);
1611:   PCCreate(comm,&newpc);
1612:   PCSetType(newpc,PCSHELL);
1613:   PCShellSetContext(newpc,fetidppc_ctx);
1614:   PCShellSetApply(newpc,FETIDPPCApply);
1615:   PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);
1616:   PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);
1617:   PCSetOperators(newpc,newmat,newmat);
1618:   PCSetUp(newpc);
1619:   /* return pointers for objects created */
1620:   *fetidp_mat=newmat;
1621:   *fetidp_pc=newpc;
1622:   return(0);
1623: }

1627: /*@
1628:  PCBDDCCreateFETIDPOperators - Create operators for FETIDP

1630:    Collective

1632:    Input Parameters:
1633: +  pc - the BDDC preconditioning context already setup

1635:    Output Parameters:
1636: .  fetidp_mat - shell FETIDP matrix object
1637: .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix

1639:    Options Database Keys:
1640: -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers

1642:    Level: developer

1644:    Notes:
1645:      Currently the only operation provided for FETIDP matrix is MatMult

1647: .seealso: PCBDDC
1648: @*/
1649: PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1650: {

1655:   if (pc->setupcalled) {
1656:     PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));
1657:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1658:   return(0);
1659: }
1660: /* -------------------------------------------------------------------------- */
1661: /*MC
1662:    PCBDDC - Balancing Domain Decomposition by Constraints.

1664:    An implementation of the BDDC preconditioner based on

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

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

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

1676:    It also works with unsymmetric and indefinite problems.

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

1680:    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace

1682:    Boundary nodes are split in vertices, edges and faces using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph()

1684:    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace().

1686:    Change of basis is performed similarly to [2] when requested. When more the one constraint is present on a single connected component (i.e. an edge or a face), a robust method based on local QR factorizations is used.

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

1690:    Options Database Keys:

1692: .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
1693: .    -pc_bddc_use_edges <1> - use or not edges in primal space
1694: .    -pc_bddc_use_faces <0> - use or not faces in primal space
1695: .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
1696: .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
1697: .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
1698: .    -pc_bddc_levels <0> - maximum number of levels for multilevel
1699: .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
1700: -    -pc_bddc_check_level <0> - set verbosity level of debugging output

1702:    Options for Dirichlet, Neumann or coarse solver can be set with
1703: .vb
1704:       -pc_bddc_dirichlet_
1705:       -pc_bddc_neumann_
1706:       -pc_bddc_coarse_
1707: .ve
1708:    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg

1710:    When using a multilevel approach, solvers' options at the N-th level can be specified as
1711: .vb
1712:       -pc_bddc_dirichlet_lN_
1713:       -pc_bddc_neumann_lN_
1714:       -pc_bddc_coarse_lN_
1715: .ve
1716:    Note that level number ranges from the finest 0 to the coarsest N.

1718:    Level: intermediate

1720:    Developer notes:

1722:    New deluxe scaling operator will be available soon.

1724:    Contributed by Stefano Zampini

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

1731: PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1732: {
1733:   PetscErrorCode      ierr;
1734:   PC_BDDC             *pcbddc;

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

1741:   /* create PCIS data structure */
1742:   PCISCreate(pc);

1744:   /* BDDC customization */
1745:   pcbddc->use_local_adj       = PETSC_TRUE;
1746:   pcbddc->use_vertices        = PETSC_TRUE;
1747:   pcbddc->use_edges           = PETSC_TRUE;
1748:   pcbddc->use_faces           = PETSC_FALSE;
1749:   pcbddc->use_change_of_basis = PETSC_FALSE;
1750:   pcbddc->use_change_on_faces = PETSC_FALSE;
1751:   pcbddc->switch_static       = PETSC_FALSE;
1752:   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
1753:   pcbddc->dbg_flag            = 0;
1754:   /* private */
1755:   pcbddc->issym                      = PETSC_FALSE;
1756:   pcbddc->BtoNmap                    = 0;
1757:   pcbddc->local_primal_size          = 0;
1758:   pcbddc->n_vertices                 = 0;
1759:   pcbddc->n_actual_vertices          = 0;
1760:   pcbddc->n_constraints              = 0;
1761:   pcbddc->primal_indices_local_idxs  = 0;
1762:   pcbddc->recompute_topography       = PETSC_FALSE;
1763:   pcbddc->coarse_size                = -1;
1764:   pcbddc->new_primal_space           = PETSC_FALSE;
1765:   pcbddc->new_primal_space_local     = PETSC_FALSE;
1766:   pcbddc->global_primal_indices      = 0;
1767:   pcbddc->onearnullspace             = 0;
1768:   pcbddc->onearnullvecs_state        = 0;
1769:   pcbddc->user_primal_vertices       = 0;
1770:   pcbddc->NullSpace                  = 0;
1771:   pcbddc->temp_solution              = 0;
1772:   pcbddc->original_rhs               = 0;
1773:   pcbddc->local_mat                  = 0;
1774:   pcbddc->ChangeOfBasisMatrix        = 0;
1775:   pcbddc->user_ChangeOfBasisMatrix   = 0;
1776:   pcbddc->coarse_vec                 = 0;
1777:   pcbddc->coarse_rhs                 = 0;
1778:   pcbddc->coarse_ksp                 = 0;
1779:   pcbddc->coarse_phi_B               = 0;
1780:   pcbddc->coarse_phi_D               = 0;
1781:   pcbddc->coarse_psi_B               = 0;
1782:   pcbddc->coarse_psi_D               = 0;
1783:   pcbddc->vec1_P                     = 0;
1784:   pcbddc->vec1_R                     = 0;
1785:   pcbddc->vec2_R                     = 0;
1786:   pcbddc->local_auxmat1              = 0;
1787:   pcbddc->local_auxmat2              = 0;
1788:   pcbddc->R_to_B                     = 0;
1789:   pcbddc->R_to_D                     = 0;
1790:   pcbddc->ksp_D                      = 0;
1791:   pcbddc->ksp_R                      = 0;
1792:   pcbddc->NeumannBoundaries          = 0;
1793:   pcbddc->NeumannBoundariesLocal     = 0;
1794:   pcbddc->DirichletBoundaries        = 0;
1795:   pcbddc->DirichletBoundariesLocal   = 0;
1796:   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
1797:   pcbddc->n_ISForDofs                = 0;
1798:   pcbddc->n_ISForDofsLocal           = 0;
1799:   pcbddc->ISForDofs                  = 0;
1800:   pcbddc->ISForDofsLocal             = 0;
1801:   pcbddc->ConstraintMatrix           = 0;
1802:   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
1803:   pcbddc->coarse_loc_to_glob         = 0;
1804:   pcbddc->coarsening_ratio           = 8;
1805:   pcbddc->current_level              = 0;
1806:   pcbddc->max_levels                 = 0;
1807:   pcbddc->use_coarse_estimates       = PETSC_FALSE;
1808:   pcbddc->coarse_subassembling       = 0;
1809:   pcbddc->coarse_subassembling_init  = 0;

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

1814:   /* scaling */
1815:   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1816:   pcbddc->deluxe_ctx                 = 0;
1817:   pcbddc->work_scaling               = 0;

1819:   /* function pointers */
1820:   pc->ops->apply               = PCApply_BDDC;
1821:   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
1822:   pc->ops->setup               = PCSetUp_BDDC;
1823:   pc->ops->destroy             = PCDestroy_BDDC;
1824:   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1825:   pc->ops->view                = 0;
1826:   pc->ops->applyrichardson     = 0;
1827:   pc->ops->applysymmetricleft  = 0;
1828:   pc->ops->applysymmetricright = 0;
1829:   pc->ops->presolve            = PCPreSolve_BDDC;
1830:   pc->ops->postsolve           = PCPostSolve_BDDC;

1832:   /* composing function */
1833:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisLocalMat_C",PCBDDCSetChangeOfBasisLocalMat_BDDC);
1834:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);
1835:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);
1836:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);
1837:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);
1838:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);
1839:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);
1840:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);
1841:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);
1842:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);
1843:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);
1844:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);
1845:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);
1846:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);
1847:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);
1848:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);
1849:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);
1850:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);
1851:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);
1852:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);
1853:   PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);
1854:   return(0);
1855: }