xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 47d04d0daa25ab4bf9a7bb108bd1a371caa366a3)
153cdbc3dSStefano Zampini /* TODOLIST
2da1bb401SStefano Zampini    DofSplitting and DM attached to pc?
3da1bb401SStefano Zampini    Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet)
4a0ba757dSStefano Zampini    change how to deal with the coarse problem (PCBDDCSetCoarseEnvironment):
5a0ba757dSStefano Zampini      - simplify coarse problem structure -> PCBDDC or PCREDUDANT, nothing else -> same comm for all levels?
6a0ba757dSStefano Zampini      - remove coarse enums and allow use of PCBDDCGetCoarseKSP
7674ae819SStefano Zampini      - remove metis dependency -> use MatPartitioning for multilevel -> Assemble serial adjacency in PCBDDCAnalyzeInterface?
8a0ba757dSStefano Zampini    code refactoring:
9a0ba757dSStefano Zampini      - pick up better names for static functions
10a0ba757dSStefano Zampini    change options structure:
11a0ba757dSStefano Zampini      - insert BDDC into MG framework?
12a0ba757dSStefano Zampini    provide other ops? Ask to developers
13a0ba757dSStefano Zampini    remove all unused printf
14a0ba757dSStefano Zampini    man pages
1553cdbc3dSStefano Zampini */
160c7d97c5SJed Brown 
1753cdbc3dSStefano Zampini /* ----------------------------------------------------------------------------------------------------------------------------------------------
180c7d97c5SJed Brown    Implementation of BDDC preconditioner based on:
190c7d97c5SJed Brown    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
2053cdbc3dSStefano Zampini    ---------------------------------------------------------------------------------------------------------------------------------------------- */
2153cdbc3dSStefano Zampini 
22674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/  /* includes for fortran wrappers */
23674ae819SStefano Zampini #include "bddcprivate.h"
243b03a366Sstefano_zampini #include <petscblaslapack.h>
25674ae819SStefano Zampini 
260c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
270c7d97c5SJed Brown #undef __FUNCT__
280c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
290c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc)
300c7d97c5SJed Brown {
310c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
320c7d97c5SJed Brown   PetscErrorCode ierr;
330c7d97c5SJed Brown 
340c7d97c5SJed Brown   PetscFunctionBegin;
350c7d97c5SJed Brown   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
368eeda7d8SStefano Zampini   /* Verbose debugging */
378eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
388eeda7d8SStefano Zampini   /* Primal space cumstomization */
398eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
408eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
418eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
428eeda7d8SStefano Zampini   /* Change of basis */
438eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr);
448eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr);
45674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
46674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
47674ae819SStefano Zampini   }
488eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
498eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);CHKERRQ(ierr);
508eeda7d8SStefano Zampini   /* Coarse solver context */
518eeda7d8SStefano Zampini   static const char * const avail_coarse_problems[] = {"sequential","replicated","parallel","multilevel","CoarseProblemType","PC_BDDC_",0}; /*order of choiches depends on ENUM defined in bddc.h */
528eeda7d8SStefano Zampini   ierr = PetscOptionsEnum("-pc_bddc_coarse_problem_type","Set coarse problem type","none",avail_coarse_problems,(PetscEnum)pcbddc->coarse_problem_type,(PetscEnum*)&pcbddc->coarse_problem_type,NULL);CHKERRQ(ierr);
530298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
540298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_max_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
55674ae819SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
560c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
570c7d97c5SJed Brown   PetscFunctionReturn(0);
580c7d97c5SJed Brown }
590c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
60674ae819SStefano Zampini #undef __FUNCT__
61674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
62674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
63674ae819SStefano Zampini {
64674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
65674ae819SStefano Zampini   PetscErrorCode ierr;
661e6b0712SBarry Smith 
67674ae819SStefano Zampini   PetscFunctionBegin;
68674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
69674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
70674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
71674ae819SStefano Zampini   PetscFunctionReturn(0);
72674ae819SStefano Zampini }
73674ae819SStefano Zampini #undef __FUNCT__
74674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
75674ae819SStefano Zampini /*@
76674ae819SStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set user defined primal vertices in PCBDDC.
77674ae819SStefano Zampini 
78674ae819SStefano Zampini    Not collective
79674ae819SStefano Zampini 
80674ae819SStefano Zampini    Input Parameters:
81674ae819SStefano Zampini +  pc - the preconditioning context
82674ae819SStefano Zampini -  PrimalVertices - index sets of primal vertices in local numbering
83674ae819SStefano Zampini 
84674ae819SStefano Zampini    Level: intermediate
85674ae819SStefano Zampini 
86674ae819SStefano Zampini    Notes:
87674ae819SStefano Zampini 
88674ae819SStefano Zampini .seealso: PCBDDC
89674ae819SStefano Zampini @*/
90674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
91674ae819SStefano Zampini {
92674ae819SStefano Zampini   PetscErrorCode ierr;
93674ae819SStefano Zampini 
94674ae819SStefano Zampini   PetscFunctionBegin;
95674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
96674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
97674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
98674ae819SStefano Zampini   PetscFunctionReturn(0);
99674ae819SStefano Zampini }
100674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1010c7d97c5SJed Brown #undef __FUNCT__
1020c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetCoarseProblemType_BDDC"
10353cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetCoarseProblemType_BDDC(PC pc, CoarseProblemType CPT)
1040c7d97c5SJed Brown {
1050c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1060c7d97c5SJed Brown 
1070c7d97c5SJed Brown   PetscFunctionBegin;
1080c7d97c5SJed Brown   pcbddc->coarse_problem_type = CPT;
1090c7d97c5SJed Brown   PetscFunctionReturn(0);
1100c7d97c5SJed Brown }
1111e6b0712SBarry Smith 
1120c7d97c5SJed Brown #undef __FUNCT__
1130c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetCoarseProblemType"
11453cdbc3dSStefano Zampini /*@
1159c0446d6SStefano Zampini  PCBDDCSetCoarseProblemType - Set coarse problem type in PCBDDC.
11653cdbc3dSStefano Zampini 
1179c0446d6SStefano Zampini    Not collective
11853cdbc3dSStefano Zampini 
11953cdbc3dSStefano Zampini    Input Parameters:
12053cdbc3dSStefano Zampini +  pc - the preconditioning context
12153cdbc3dSStefano Zampini -  CoarseProblemType - pick a better name and explain what this is
12253cdbc3dSStefano Zampini 
12353cdbc3dSStefano Zampini    Level: intermediate
12453cdbc3dSStefano Zampini 
12553cdbc3dSStefano Zampini    Notes:
126da1bb401SStefano Zampini    Not collective but all procs must call with same arguments.
12753cdbc3dSStefano Zampini 
12853cdbc3dSStefano Zampini .seealso: PCBDDC
12953cdbc3dSStefano Zampini @*/
1300c7d97c5SJed Brown PetscErrorCode PCBDDCSetCoarseProblemType(PC pc, CoarseProblemType CPT)
1310c7d97c5SJed Brown {
1320c7d97c5SJed Brown   PetscErrorCode ierr;
1330c7d97c5SJed Brown 
1340c7d97c5SJed Brown   PetscFunctionBegin;
1350c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1360c7d97c5SJed Brown   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseProblemType_C",(PC,CoarseProblemType),(pc,CPT));CHKERRQ(ierr);
1370c7d97c5SJed Brown   PetscFunctionReturn(0);
1380c7d97c5SJed Brown }
1390c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1400c7d97c5SJed Brown #undef __FUNCT__
1414fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1424fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1434fad6a16SStefano Zampini {
1444fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1454fad6a16SStefano Zampini 
1464fad6a16SStefano Zampini   PetscFunctionBegin;
1474fad6a16SStefano Zampini   pcbddc->coarsening_ratio=k;
1484fad6a16SStefano Zampini   PetscFunctionReturn(0);
1494fad6a16SStefano Zampini }
1501e6b0712SBarry Smith 
1514fad6a16SStefano Zampini #undef __FUNCT__
1524fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1534fad6a16SStefano Zampini /*@
1544fad6a16SStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel coarsening
1554fad6a16SStefano Zampini 
1564fad6a16SStefano Zampini    Logically collective on PC
1574fad6a16SStefano Zampini 
1584fad6a16SStefano Zampini    Input Parameters:
1594fad6a16SStefano Zampini +  pc - the preconditioning context
1604fad6a16SStefano Zampini -  k - coarsening ratio
1614fad6a16SStefano Zampini 
1624fad6a16SStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level.
1634fad6a16SStefano Zampini 
1644fad6a16SStefano Zampini    Level: intermediate
1654fad6a16SStefano Zampini 
1664fad6a16SStefano Zampini    Notes:
1674fad6a16SStefano Zampini 
1684fad6a16SStefano Zampini .seealso: PCBDDC
1694fad6a16SStefano Zampini @*/
1704fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
1714fad6a16SStefano Zampini {
1724fad6a16SStefano Zampini   PetscErrorCode ierr;
1734fad6a16SStefano Zampini 
1744fad6a16SStefano Zampini   PetscFunctionBegin;
1754fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1764fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
1774fad6a16SStefano Zampini   PetscFunctionReturn(0);
1784fad6a16SStefano Zampini }
1794fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
1801e6b0712SBarry Smith 
1814fad6a16SStefano Zampini #undef __FUNCT__
1824fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetMaxLevels_BDDC"
1834fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetMaxLevels_BDDC(PC pc,PetscInt max_levels)
1844fad6a16SStefano Zampini {
1854fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1864fad6a16SStefano Zampini 
1874fad6a16SStefano Zampini   PetscFunctionBegin;
1884fad6a16SStefano Zampini   pcbddc->max_levels=max_levels;
1894fad6a16SStefano Zampini   PetscFunctionReturn(0);
1904fad6a16SStefano Zampini }
1911e6b0712SBarry Smith 
1924fad6a16SStefano Zampini #undef __FUNCT__
1934fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetMaxLevels"
1944fad6a16SStefano Zampini /*@
1954fad6a16SStefano Zampini  PCBDDCSetMaxLevels - Sets the maximum number of levels within the multilevel approach.
1964fad6a16SStefano Zampini 
1974fad6a16SStefano Zampini    Logically collective on PC
1984fad6a16SStefano Zampini 
1994fad6a16SStefano Zampini    Input Parameters:
2004fad6a16SStefano Zampini +  pc - the preconditioning context
2014fad6a16SStefano Zampini -  max_levels - the maximum number of levels
2024fad6a16SStefano Zampini 
2034fad6a16SStefano Zampini    Default value is 1, i.e. coarse problem will be solved inexactly with one application
2044fad6a16SStefano Zampini    of PCBDDC preconditioner if the multilevel approach is requested.
2054fad6a16SStefano Zampini 
2064fad6a16SStefano Zampini    Level: intermediate
2074fad6a16SStefano Zampini 
2084fad6a16SStefano Zampini    Notes:
2094fad6a16SStefano Zampini 
2104fad6a16SStefano Zampini .seealso: PCBDDC
2114fad6a16SStefano Zampini @*/
2124fad6a16SStefano Zampini PetscErrorCode PCBDDCSetMaxLevels(PC pc,PetscInt max_levels)
2134fad6a16SStefano Zampini {
2144fad6a16SStefano Zampini   PetscErrorCode ierr;
2154fad6a16SStefano Zampini 
2164fad6a16SStefano Zampini   PetscFunctionBegin;
2174fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2184fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetMaxLevels_C",(PC,PetscInt),(pc,max_levels));CHKERRQ(ierr);
2194fad6a16SStefano Zampini   PetscFunctionReturn(0);
2204fad6a16SStefano Zampini }
2214fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
2221e6b0712SBarry Smith 
2234fad6a16SStefano Zampini #undef __FUNCT__
2240bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
2250bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
2260bdf917eSStefano Zampini {
2270bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2280bdf917eSStefano Zampini   PetscErrorCode ierr;
2290bdf917eSStefano Zampini 
2300bdf917eSStefano Zampini   PetscFunctionBegin;
2310bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
2320bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
2330bdf917eSStefano Zampini   pcbddc->NullSpace=NullSpace;
2340bdf917eSStefano Zampini   PetscFunctionReturn(0);
2350bdf917eSStefano Zampini }
2361e6b0712SBarry Smith 
2370bdf917eSStefano Zampini #undef __FUNCT__
2380bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
2390bdf917eSStefano Zampini /*@
2400bdf917eSStefano Zampini  PCBDDCSetNullSpace - Set NullSpace of global operator of BDDC preconditioned mat.
2410bdf917eSStefano Zampini 
2420bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
2430bdf917eSStefano Zampini 
2440bdf917eSStefano Zampini    Input Parameters:
2450bdf917eSStefano Zampini +  pc - the preconditioning context
2460bdf917eSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned.
2470bdf917eSStefano Zampini 
2480bdf917eSStefano Zampini    Level: intermediate
2490bdf917eSStefano Zampini 
2500bdf917eSStefano Zampini    Notes:
2510bdf917eSStefano Zampini 
2520bdf917eSStefano Zampini .seealso: PCBDDC
2530bdf917eSStefano Zampini @*/
2540bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
2550bdf917eSStefano Zampini {
2560bdf917eSStefano Zampini   PetscErrorCode ierr;
2570bdf917eSStefano Zampini 
2580bdf917eSStefano Zampini   PetscFunctionBegin;
2590bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
260674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
2610bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
2620bdf917eSStefano Zampini   PetscFunctionReturn(0);
2630bdf917eSStefano Zampini }
2640bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
2651e6b0712SBarry Smith 
2660bdf917eSStefano Zampini #undef __FUNCT__
2673b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
2683b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
2693b03a366Sstefano_zampini {
2703b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2713b03a366Sstefano_zampini   PetscErrorCode ierr;
2723b03a366Sstefano_zampini 
2733b03a366Sstefano_zampini   PetscFunctionBegin;
2743b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
27536e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
27636e030ebSStefano Zampini   pcbddc->DirichletBoundaries=DirichletBoundaries;
2773b03a366Sstefano_zampini   PetscFunctionReturn(0);
2783b03a366Sstefano_zampini }
2791e6b0712SBarry Smith 
2803b03a366Sstefano_zampini #undef __FUNCT__
2813b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
2823b03a366Sstefano_zampini /*@
283da1bb401SStefano Zampini  PCBDDCSetDirichletBoundaries - Set index set defining subdomain part (in local ordering)
284da1bb401SStefano Zampini                               of Dirichlet boundaries for the global problem.
2853b03a366Sstefano_zampini 
2863b03a366Sstefano_zampini    Not collective
2873b03a366Sstefano_zampini 
2883b03a366Sstefano_zampini    Input Parameters:
2893b03a366Sstefano_zampini +  pc - the preconditioning context
2900298fd71SBarry Smith -  DirichletBoundaries - sequential index set defining the subdomain part of Dirichlet boundaries (can be NULL)
2913b03a366Sstefano_zampini 
2923b03a366Sstefano_zampini    Level: intermediate
2933b03a366Sstefano_zampini 
2943b03a366Sstefano_zampini    Notes:
2953b03a366Sstefano_zampini 
2963b03a366Sstefano_zampini .seealso: PCBDDC
2973b03a366Sstefano_zampini @*/
2983b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
2993b03a366Sstefano_zampini {
3003b03a366Sstefano_zampini   PetscErrorCode ierr;
3013b03a366Sstefano_zampini 
3023b03a366Sstefano_zampini   PetscFunctionBegin;
3033b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
304674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
3053b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
3063b03a366Sstefano_zampini   PetscFunctionReturn(0);
3073b03a366Sstefano_zampini }
3083b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
3091e6b0712SBarry Smith 
3103b03a366Sstefano_zampini #undef __FUNCT__
3110c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
31253cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
3130c7d97c5SJed Brown {
3140c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
31553cdbc3dSStefano Zampini   PetscErrorCode ierr;
3160c7d97c5SJed Brown 
3170c7d97c5SJed Brown   PetscFunctionBegin;
31853cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
31936e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
32036e030ebSStefano Zampini   pcbddc->NeumannBoundaries=NeumannBoundaries;
3210c7d97c5SJed Brown   PetscFunctionReturn(0);
3220c7d97c5SJed Brown }
3231e6b0712SBarry Smith 
3240c7d97c5SJed Brown #undef __FUNCT__
3250c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
32657527edcSJed Brown /*@
327da1bb401SStefano Zampini  PCBDDCSetNeumannBoundaries - Set index set defining subdomain part (in local ordering)
328da1bb401SStefano Zampini                               of Neumann boundaries for the global problem.
32957527edcSJed Brown 
3309c0446d6SStefano Zampini    Not collective
33157527edcSJed Brown 
33257527edcSJed Brown    Input Parameters:
33357527edcSJed Brown +  pc - the preconditioning context
3340298fd71SBarry Smith -  NeumannBoundaries - sequential index set defining the subdomain part of Neumann boundaries (can be NULL)
33557527edcSJed Brown 
33657527edcSJed Brown    Level: intermediate
33757527edcSJed Brown 
33857527edcSJed Brown    Notes:
33957527edcSJed Brown 
34057527edcSJed Brown .seealso: PCBDDC
34157527edcSJed Brown @*/
34253cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
3430c7d97c5SJed Brown {
3440c7d97c5SJed Brown   PetscErrorCode ierr;
3450c7d97c5SJed Brown 
3460c7d97c5SJed Brown   PetscFunctionBegin;
3470c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
348674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
34953cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
35053cdbc3dSStefano Zampini   PetscFunctionReturn(0);
35153cdbc3dSStefano Zampini }
35253cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
3531e6b0712SBarry Smith 
35453cdbc3dSStefano Zampini #undef __FUNCT__
355da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
356da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
357da1bb401SStefano Zampini {
358da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
359da1bb401SStefano Zampini 
360da1bb401SStefano Zampini   PetscFunctionBegin;
361da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
362da1bb401SStefano Zampini   PetscFunctionReturn(0);
363da1bb401SStefano Zampini }
3641e6b0712SBarry Smith 
365da1bb401SStefano Zampini #undef __FUNCT__
366da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
367da1bb401SStefano Zampini /*@
368da1bb401SStefano Zampini  PCBDDCGetDirichletBoundaries - Get index set defining subdomain part (in local ordering)
369da1bb401SStefano Zampini                                 of Dirichlet boundaries for the global problem.
370da1bb401SStefano Zampini 
371da1bb401SStefano Zampini    Not collective
372da1bb401SStefano Zampini 
373da1bb401SStefano Zampini    Input Parameters:
374da1bb401SStefano Zampini +  pc - the preconditioning context
375da1bb401SStefano Zampini 
376da1bb401SStefano Zampini    Output Parameters:
377da1bb401SStefano Zampini +  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
378da1bb401SStefano Zampini 
379da1bb401SStefano Zampini    Level: intermediate
380da1bb401SStefano Zampini 
381da1bb401SStefano Zampini    Notes:
382da1bb401SStefano Zampini 
383da1bb401SStefano Zampini .seealso: PCBDDC
384da1bb401SStefano Zampini @*/
385da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
386da1bb401SStefano Zampini {
387da1bb401SStefano Zampini   PetscErrorCode ierr;
388da1bb401SStefano Zampini 
389da1bb401SStefano Zampini   PetscFunctionBegin;
390da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
391da1bb401SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
392da1bb401SStefano Zampini   PetscFunctionReturn(0);
393da1bb401SStefano Zampini }
394da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
3951e6b0712SBarry Smith 
396da1bb401SStefano Zampini #undef __FUNCT__
39753cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
39853cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
39953cdbc3dSStefano Zampini {
40053cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
40153cdbc3dSStefano Zampini 
40253cdbc3dSStefano Zampini   PetscFunctionBegin;
40353cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
40453cdbc3dSStefano Zampini   PetscFunctionReturn(0);
40553cdbc3dSStefano Zampini }
4061e6b0712SBarry Smith 
40753cdbc3dSStefano Zampini #undef __FUNCT__
40853cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
40953cdbc3dSStefano Zampini /*@
410da1bb401SStefano Zampini  PCBDDCGetNeumannBoundaries - Get index set defining subdomain part (in local ordering)
411da1bb401SStefano Zampini                               of Neumann boundaries for the global problem.
41253cdbc3dSStefano Zampini 
4139c0446d6SStefano Zampini    Not collective
41453cdbc3dSStefano Zampini 
41553cdbc3dSStefano Zampini    Input Parameters:
41653cdbc3dSStefano Zampini +  pc - the preconditioning context
41753cdbc3dSStefano Zampini 
41853cdbc3dSStefano Zampini    Output Parameters:
41953cdbc3dSStefano Zampini +  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
42053cdbc3dSStefano Zampini 
42153cdbc3dSStefano Zampini    Level: intermediate
42253cdbc3dSStefano Zampini 
42353cdbc3dSStefano Zampini    Notes:
42453cdbc3dSStefano Zampini 
42553cdbc3dSStefano Zampini .seealso: PCBDDC
42653cdbc3dSStefano Zampini @*/
42753cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
42853cdbc3dSStefano Zampini {
42953cdbc3dSStefano Zampini   PetscErrorCode ierr;
43053cdbc3dSStefano Zampini 
43153cdbc3dSStefano Zampini   PetscFunctionBegin;
43253cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
43353cdbc3dSStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
4340c7d97c5SJed Brown   PetscFunctionReturn(0);
4350c7d97c5SJed Brown }
43636e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
4371e6b0712SBarry Smith 
43836e030ebSStefano Zampini #undef __FUNCT__
439da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
4401a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
44136e030ebSStefano Zampini {
44236e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
443da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
444da1bb401SStefano Zampini   PetscErrorCode ierr;
44536e030ebSStefano Zampini 
44636e030ebSStefano Zampini   PetscFunctionBegin;
447674ae819SStefano Zampini   /* free old CSR */
448674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
449fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
450674ae819SStefano Zampini   /* get CSR into graph structure */
451da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
452674ae819SStefano Zampini     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
453674ae819SStefano Zampini     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
454674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
455674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
456da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
4571a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
4581a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
459674ae819SStefano Zampini   }
460575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
46136e030ebSStefano Zampini   PetscFunctionReturn(0);
46236e030ebSStefano Zampini }
4631e6b0712SBarry Smith 
46436e030ebSStefano Zampini #undef __FUNCT__
465da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
46636e030ebSStefano Zampini /*@
467da1bb401SStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set CSR graph of local matrix for use of PCBDDC.
46836e030ebSStefano Zampini 
46936e030ebSStefano Zampini    Not collective
47036e030ebSStefano Zampini 
47136e030ebSStefano Zampini    Input Parameters:
47236e030ebSStefano Zampini +  pc - the preconditioning context
473da1bb401SStefano Zampini -  nvtxs - number of local vertices of the graph
474da1bb401SStefano Zampini -  xadj, adjncy - the CSR graph
475da1bb401SStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. In the former case the user must free the array passed in;
476da1bb401SStefano Zampini                                                              in the latter case, memory must be obtained with PetscMalloc.
47736e030ebSStefano Zampini 
47836e030ebSStefano Zampini    Level: intermediate
47936e030ebSStefano Zampini 
48036e030ebSStefano Zampini    Notes:
48136e030ebSStefano Zampini 
48236e030ebSStefano Zampini .seealso: PCBDDC
48336e030ebSStefano Zampini @*/
4841a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
48536e030ebSStefano Zampini {
486575ad6abSStefano Zampini   void (*f)(void) = 0;
48736e030ebSStefano Zampini   PetscErrorCode ierr;
48836e030ebSStefano Zampini 
48936e030ebSStefano Zampini   PetscFunctionBegin;
49036e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
491674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
492674ae819SStefano Zampini   PetscValidIntPointer(xadj,4);
493674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
494674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
495da1bb401SStefano Zampini   }
49636e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
497575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
498575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
499575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
500575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
501575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
50236e030ebSStefano Zampini   }
50336e030ebSStefano Zampini   PetscFunctionReturn(0);
50436e030ebSStefano Zampini }
5059c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
5061e6b0712SBarry Smith 
5079c0446d6SStefano Zampini #undef __FUNCT__
5089c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
5099c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
5109c0446d6SStefano Zampini {
5119c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5129c0446d6SStefano Zampini   PetscInt i;
5139c0446d6SStefano Zampini   PetscErrorCode ierr;
5149c0446d6SStefano Zampini 
5159c0446d6SStefano Zampini   PetscFunctionBegin;
516da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
5179c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
5189c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
5199c0446d6SStefano Zampini   }
520d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
521da1bb401SStefano Zampini   /* allocate space then set */
5229c0446d6SStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
5239c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
524da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
525da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
5269c0446d6SStefano Zampini   }
5279c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
5289c0446d6SStefano Zampini   PetscFunctionReturn(0);
5299c0446d6SStefano Zampini }
5301e6b0712SBarry Smith 
5319c0446d6SStefano Zampini #undef __FUNCT__
5329c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
5339c0446d6SStefano Zampini /*@
534da1bb401SStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of local mat.
5359c0446d6SStefano Zampini 
5369c0446d6SStefano Zampini    Not collective
5379c0446d6SStefano Zampini 
5389c0446d6SStefano Zampini    Input Parameters:
5399c0446d6SStefano Zampini +  pc - the preconditioning context
540da1bb401SStefano Zampini -  n - number of index sets defining the fields
541da1bb401SStefano Zampini -  IS[] - array of IS describing the fields
5429c0446d6SStefano Zampini 
5439c0446d6SStefano Zampini    Level: intermediate
5449c0446d6SStefano Zampini 
5459c0446d6SStefano Zampini    Notes:
5469c0446d6SStefano Zampini 
5479c0446d6SStefano Zampini .seealso: PCBDDC
5489c0446d6SStefano Zampini @*/
5499c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
5509c0446d6SStefano Zampini {
5519c0446d6SStefano Zampini   PetscErrorCode ierr;
5529c0446d6SStefano Zampini 
5539c0446d6SStefano Zampini   PetscFunctionBegin;
5549c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5559c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
5569c0446d6SStefano Zampini   PetscFunctionReturn(0);
5579c0446d6SStefano Zampini }
558da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
559534831adSStefano Zampini #undef __FUNCT__
560534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
561534831adSStefano Zampini /* -------------------------------------------------------------------------- */
562534831adSStefano Zampini /*
563534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
564534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
5659c0446d6SStefano Zampini 
566534831adSStefano Zampini    Input Parameter:
567534831adSStefano Zampini +  pc - the preconditioner contex
568534831adSStefano Zampini 
569534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
570534831adSStefano Zampini 
571534831adSStefano Zampini    Notes:
572534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
573534831adSStefano Zampini    the user, but instead is called by KSPSolve().
574534831adSStefano Zampini */
575534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
576534831adSStefano Zampini {
577534831adSStefano Zampini   PetscErrorCode ierr;
578534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
579534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
580534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
581534831adSStefano Zampini   Mat            temp_mat;
5823972b0daSStefano Zampini   IS             dirIS;
5833972b0daSStefano Zampini   PetscInt       dirsize,i,*is_indices;
5843972b0daSStefano Zampini   PetscScalar    *array_x,*array_diagonal;
5853972b0daSStefano Zampini   Vec            used_vec;
5863972b0daSStefano Zampini   PetscBool      guess_nonzero;
587534831adSStefano Zampini 
588534831adSStefano Zampini   PetscFunctionBegin;
58962a6ff1dSStefano Zampini   /* Creates parallel work vectors used in presolve. */
59062a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
59162a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
59262a6ff1dSStefano Zampini   }
59362a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
59462a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
59562a6ff1dSStefano Zampini   }
5963972b0daSStefano Zampini   if (x) {
5973972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
5983972b0daSStefano Zampini     used_vec = x;
5993972b0daSStefano Zampini   } else {
6003972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
6013972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
6023972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6033972b0daSStefano Zampini   }
6043972b0daSStefano Zampini   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
6053972b0daSStefano Zampini   if (ksp) {
6063972b0daSStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
6073972b0daSStefano Zampini     if ( !guess_nonzero ) {
6083972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6093972b0daSStefano Zampini     }
6103972b0daSStefano Zampini   }
6113308cffdSStefano Zampini 
61262a6ff1dSStefano Zampini   if (rhs) { /* TODO: wiser handling of rhs removal, which is only needed in case of zeroed rows */
6133972b0daSStefano Zampini     /* store the original rhs */
6143972b0daSStefano Zampini     ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
6153972b0daSStefano Zampini 
6163972b0daSStefano Zampini     /* Take into account zeroed rows -> change rhs and store solution removed */
6173972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
6183972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
6193972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6203972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6213972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6223972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6233972b0daSStefano Zampini     ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr);
6243972b0daSStefano Zampini     if (dirIS) {
6253972b0daSStefano Zampini       ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr);
6263972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6273972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6283972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6292fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
6303972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6313972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6323972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6333972b0daSStefano Zampini     }
6343972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6353972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
636b76ba322SStefano Zampini 
6373972b0daSStefano Zampini     /* remove the computed solution from the rhs */
6383972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6393972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
6403972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6413308cffdSStefano Zampini   }
642b76ba322SStefano Zampini 
643b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
6443972b0daSStefano Zampini   if (x) {
6453972b0daSStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
6463972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
64715aaf578SStefano Zampini     if (pcbddc->use_exact_dirichlet && !pcbddc->coarse_psi_B) {
648b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
649b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
650b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
651b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
652b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
653b76ba322SStefano Zampini       if (ksp) {
654b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
655b76ba322SStefano Zampini       }
656b76ba322SStefano Zampini     }
6573972b0daSStefano Zampini   }
658b76ba322SStefano Zampini 
659674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
660b76ba322SStefano Zampini     /* swap pointers for local matrices */
661b76ba322SStefano Zampini     temp_mat = matis->A;
662b76ba322SStefano Zampini     matis->A = pcbddc->local_mat;
663b76ba322SStefano Zampini     pcbddc->local_mat = temp_mat;
6643308cffdSStefano Zampini   }
6653308cffdSStefano Zampini   if (pcbddc->use_change_of_basis && rhs) {
666b76ba322SStefano Zampini     /* Get local rhs and apply transformation of basis */
667b76ba322SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
668b76ba322SStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
669b76ba322SStefano Zampini     /* from original basis to modified basis */
670b76ba322SStefano Zampini     ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
671b76ba322SStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
672b76ba322SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
673b76ba322SStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
674674ae819SStefano Zampini   }
6750bdf917eSStefano Zampini   if (ksp && pcbddc->NullSpace) {
676d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
677d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
678b76ba322SStefano Zampini   }
6790bdf917eSStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
680534831adSStefano Zampini   PetscFunctionReturn(0);
681534831adSStefano Zampini }
682534831adSStefano Zampini /* -------------------------------------------------------------------------- */
683534831adSStefano Zampini #undef __FUNCT__
684534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
685534831adSStefano Zampini /* -------------------------------------------------------------------------- */
686534831adSStefano Zampini /*
687534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
688534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
689534831adSStefano Zampini 
690534831adSStefano Zampini    Input Parameter:
691534831adSStefano Zampini +  pc - the preconditioner contex
692534831adSStefano Zampini 
693534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
694534831adSStefano Zampini 
695534831adSStefano Zampini    Notes:
696534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
697534831adSStefano Zampini    the user, but instead is called by KSPSolve().
698534831adSStefano Zampini */
699534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
700534831adSStefano Zampini {
701534831adSStefano Zampini   PetscErrorCode ierr;
702534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
703534831adSStefano Zampini   PC_IS          *pcis   = (PC_IS*)(pc->data);
704534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
705534831adSStefano Zampini   Mat            temp_mat;
706534831adSStefano Zampini 
707534831adSStefano Zampini   PetscFunctionBegin;
708674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
709534831adSStefano Zampini     /* swap pointers for local matrices */
710534831adSStefano Zampini     temp_mat = matis->A;
711534831adSStefano Zampini     matis->A = pcbddc->local_mat;
712534831adSStefano Zampini     pcbddc->local_mat = temp_mat;
7133425bc38SStefano Zampini   }
7143308cffdSStefano Zampini   if (pcbddc->use_change_of_basis && x) {
715534831adSStefano Zampini     /* Get Local boundary and apply transformation of basis to solution vector */
716534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
717534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
718534831adSStefano Zampini     /* from modified basis to original basis */
719534831adSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
720534831adSStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
721534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
722534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
723534831adSStefano Zampini   }
7243972b0daSStefano Zampini   /* add solution removed in presolve */
7253425bc38SStefano Zampini   if (x) {
7263425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
7273425bc38SStefano Zampini   }
728fb223d50SStefano Zampini   /* restore rhs to its original state */
729fb223d50SStefano Zampini   if (rhs) {
730fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
731fb223d50SStefano Zampini   }
732534831adSStefano Zampini   PetscFunctionReturn(0);
733534831adSStefano Zampini }
734534831adSStefano Zampini /* -------------------------------------------------------------------------- */
73553cdbc3dSStefano Zampini #undef __FUNCT__
73653cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
7370c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
7380c7d97c5SJed Brown /*
7390c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
7400c7d97c5SJed Brown                   by setting data structures and options.
7410c7d97c5SJed Brown 
7420c7d97c5SJed Brown    Input Parameter:
74353cdbc3dSStefano Zampini +  pc - the preconditioner context
7440c7d97c5SJed Brown 
7450c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
7460c7d97c5SJed Brown 
7470c7d97c5SJed Brown    Notes:
7480c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
7490c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
7500c7d97c5SJed Brown */
75153cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
7520c7d97c5SJed Brown {
7530c7d97c5SJed Brown   PetscErrorCode ierr;
7540c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
755674ae819SStefano Zampini   MatStructure   flag;
756674ae819SStefano Zampini   PetscBool      computeis,computetopography,computesolvers;
7570c7d97c5SJed Brown 
7580c7d97c5SJed Brown   PetscFunctionBegin;
759674ae819SStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */
7603b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
7619c0446d6SStefano Zampini      So, we set to pcnone the Neumann problem of pcis in order to avoid unneeded computation
7620c7d97c5SJed Brown      Also, we decide to directly build the (same) Dirichlet problem */
7630c7d97c5SJed Brown   ierr = PetscOptionsSetValue("-is_localN_pc_type","none");CHKERRQ(ierr);
7640c7d97c5SJed Brown   ierr = PetscOptionsSetValue("-is_localD_pc_type","none");CHKERRQ(ierr);
7653b03a366Sstefano_zampini   /* Get stdout for dbg */
766674ae819SStefano Zampini   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
767ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
768e269702eSStefano Zampini     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
769e269702eSStefano Zampini   }
770674ae819SStefano Zampini   /* first attempt to split work */
771674ae819SStefano Zampini   if (pc->setupcalled) {
772674ae819SStefano Zampini     computeis = PETSC_FALSE;
773674ae819SStefano Zampini     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
774674ae819SStefano Zampini     if (flag == SAME_PRECONDITIONER) {
775674ae819SStefano Zampini       computetopography = PETSC_FALSE;
776674ae819SStefano Zampini       computesolvers = PETSC_FALSE;
777674ae819SStefano Zampini     } else if (flag == SAME_NONZERO_PATTERN) {
778674ae819SStefano Zampini       computetopography = PETSC_FALSE;
779674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
780674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
781674ae819SStefano Zampini       computetopography = PETSC_TRUE;
782674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
783674ae819SStefano Zampini     }
784674ae819SStefano Zampini   } else {
785674ae819SStefano Zampini     computeis = PETSC_TRUE;
786674ae819SStefano Zampini     computetopography = PETSC_TRUE;
787674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
788674ae819SStefano Zampini   }
789674ae819SStefano Zampini   /* Set up all the "iterative substructuring" common block */
790674ae819SStefano Zampini   if (computeis) {
791674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
792674ae819SStefano Zampini   }
793674ae819SStefano Zampini   /* Analyze interface and set up local constraint and change of basis matrices */
794674ae819SStefano Zampini   if (computetopography) {
795674ae819SStefano Zampini     /* reset data */
796674ae819SStefano Zampini     ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
797674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
798674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
799674ae819SStefano Zampini   }
800674ae819SStefano Zampini   if (computesolvers) {
801674ae819SStefano Zampini     /* reset data */
802674ae819SStefano Zampini     ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
803674ae819SStefano Zampini     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
80499cc7994SStefano Zampini     /* Create coarse and local stuffs */
80599cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
806674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
8070c7d97c5SJed Brown   }
8080c7d97c5SJed Brown   PetscFunctionReturn(0);
8090c7d97c5SJed Brown }
8100c7d97c5SJed Brown 
8110c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
8120c7d97c5SJed Brown /*
8130c7d97c5SJed Brown    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
8140c7d97c5SJed Brown 
8150c7d97c5SJed Brown    Input Parameters:
8160c7d97c5SJed Brown .  pc - the preconditioner context
8170c7d97c5SJed Brown .  r - input vector (global)
8180c7d97c5SJed Brown 
8190c7d97c5SJed Brown    Output Parameter:
8200c7d97c5SJed Brown .  z - output vector (global)
8210c7d97c5SJed Brown 
8220c7d97c5SJed Brown    Application Interface Routine: PCApply()
8230c7d97c5SJed Brown  */
8240c7d97c5SJed Brown #undef __FUNCT__
8250c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
82653cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
8270c7d97c5SJed Brown {
8280c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
8290c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
8300c7d97c5SJed Brown   PetscErrorCode    ierr;
8313b03a366Sstefano_zampini   const PetscScalar one = 1.0;
8323b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
8332617d88aSStefano Zampini   const PetscScalar zero = 0.0;
8340c7d97c5SJed Brown 
8350c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
8360c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
8378eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
8380c7d97c5SJed Brown 
8390c7d97c5SJed Brown   PetscFunctionBegin;
84015aaf578SStefano Zampini   if (!pcbddc->use_exact_dirichlet || pcbddc->coarse_psi_B) {
8410c7d97c5SJed Brown     /* First Dirichlet solve */
8420c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8430c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
84453cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
8450c7d97c5SJed Brown     /*
8460c7d97c5SJed Brown       Assembling right hand side for BDDC operator
847674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
848674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
8490c7d97c5SJed Brown     */
8500c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
8510c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
8528eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
8530c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
8540c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
8550c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8560c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
857674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
858b76ba322SStefano Zampini   } else {
8590bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
860b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
861674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
862b76ba322SStefano Zampini   }
863b76ba322SStefano Zampini 
8642617d88aSStefano Zampini   /* Apply interface preconditioner
8652617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
8662617d88aSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
8672617d88aSStefano Zampini 
868674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
869674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
8700c7d97c5SJed Brown 
8713b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
8720c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8730c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8740c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
8758eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
87653cdbc3dSStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr);
8770c7d97c5SJed Brown   ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr);
8788eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
8790c7d97c5SJed Brown   ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr);
8800c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8810c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8820c7d97c5SJed Brown   PetscFunctionReturn(0);
8830c7d97c5SJed Brown }
884da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
885674ae819SStefano Zampini 
886da1bb401SStefano Zampini #undef __FUNCT__
887da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
888da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
889da1bb401SStefano Zampini {
890da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
891da1bb401SStefano Zampini   PetscErrorCode ierr;
892da1bb401SStefano Zampini 
893da1bb401SStefano Zampini   PetscFunctionBegin;
894da1bb401SStefano Zampini   /* free data created by PCIS */
895da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
896674ae819SStefano Zampini   /* free BDDC custom data  */
897674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
898674ae819SStefano Zampini   /* destroy objects related to topography */
899674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
900674ae819SStefano Zampini   /* free allocated graph structure */
901da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
902674ae819SStefano Zampini   /* free data for scaling operator */
903674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
904674ae819SStefano Zampini   /* free solvers stuff */
905674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
90633bc96a4SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
90733bc96a4SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
908ac78edfcSStefano Zampini   ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
90962a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
91062a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
91162a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
9123425bc38SStefano Zampini   /* remove functions */
913674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
914bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
915bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetMaxLevels_C",NULL);CHKERRQ(ierr);
916bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
917bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
918bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
919bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
920bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
921bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseProblemType_C",NULL);CHKERRQ(ierr);
922bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
923bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
924bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
925bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
926bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
927674ae819SStefano Zampini   /* Free the private data structure */
928674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
929da1bb401SStefano Zampini   PetscFunctionReturn(0);
930da1bb401SStefano Zampini }
9313425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
9321e6b0712SBarry Smith 
9333425bc38SStefano Zampini #undef __FUNCT__
9343425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
9353425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
9363425bc38SStefano Zampini {
937674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
9383425bc38SStefano Zampini   PC_IS*         pcis;
9393425bc38SStefano Zampini   PC_BDDC*       pcbddc;
9403425bc38SStefano Zampini   PetscErrorCode ierr;
9410c7d97c5SJed Brown 
9423425bc38SStefano Zampini   PetscFunctionBegin;
9433425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
9443425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
9453425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
9463425bc38SStefano Zampini 
9473425bc38SStefano Zampini   /* change of basis for physical rhs if needed
9483425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
9493308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
9503425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
9513425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9523425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
953fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
954fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
9553425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9563425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
957674ae819SStefano Zampini   /* Apply partition of unity */
9583425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
959674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
9608eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
9613425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
9623425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
9633425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
9643425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
9653425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
9663425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9673425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
968674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
9693425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9703425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9713425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
9723425bc38SStefano Zampini   }
9733425bc38SStefano Zampini   /* BDDC rhs */
9743425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
9758eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
9763425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
9773425bc38SStefano Zampini   }
9783425bc38SStefano Zampini   /* apply BDDC */
9793425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
9803425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
9813425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
9823425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
9833425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9843425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9853425bc38SStefano Zampini   /* restore original rhs */
9863425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
9873425bc38SStefano Zampini   PetscFunctionReturn(0);
9883425bc38SStefano Zampini }
9891e6b0712SBarry Smith 
9903425bc38SStefano Zampini #undef __FUNCT__
9913425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
9923425bc38SStefano Zampini /*@
9933425bc38SStefano Zampini  PCBDDCMatFETIDPGetRHS - Get rhs for FETIDP linear system.
9943425bc38SStefano Zampini 
9953425bc38SStefano Zampini    Collective
9963425bc38SStefano Zampini 
9973425bc38SStefano Zampini    Input Parameters:
9983425bc38SStefano Zampini +  fetidp_mat   - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators
9993425bc38SStefano Zampini +  standard_rhs - the rhs of your linear system
10003425bc38SStefano Zampini 
10013425bc38SStefano Zampini    Output Parameters:
10023425bc38SStefano Zampini +  fetidp_flux_rhs   - the rhs of the FETIDP linear system
10033425bc38SStefano Zampini 
10043425bc38SStefano Zampini    Level: developer
10053425bc38SStefano Zampini 
10063425bc38SStefano Zampini    Notes:
10073425bc38SStefano Zampini 
10083425bc38SStefano Zampini .seealso: PCBDDC
10093425bc38SStefano Zampini @*/
10103425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
10113425bc38SStefano Zampini {
1012674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10133425bc38SStefano Zampini   PetscErrorCode ierr;
10143425bc38SStefano Zampini 
10153425bc38SStefano Zampini   PetscFunctionBegin;
10163425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10173425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
10183425bc38SStefano Zampini   PetscFunctionReturn(0);
10193425bc38SStefano Zampini }
10203425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
10211e6b0712SBarry Smith 
10223425bc38SStefano Zampini #undef __FUNCT__
10233425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
10243425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
10253425bc38SStefano Zampini {
1026674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10273425bc38SStefano Zampini   PC_IS*         pcis;
10283425bc38SStefano Zampini   PC_BDDC*       pcbddc;
10293425bc38SStefano Zampini   PetscErrorCode ierr;
10303425bc38SStefano Zampini 
10313425bc38SStefano Zampini   PetscFunctionBegin;
10323425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10333425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
10343425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
10353425bc38SStefano Zampini 
10363425bc38SStefano Zampini   /* apply B_delta^T */
10373425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10383425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10393425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
10403425bc38SStefano Zampini   /* compute rhs for BDDC application */
10413425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
10428eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
10433425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10443425bc38SStefano Zampini   }
10453425bc38SStefano Zampini   /* apply BDDC */
10463425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
10473425bc38SStefano Zampini   /* put values into standard global vector */
10483425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10493425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10508eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
10513425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
10523425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
10533425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
10543425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10553425bc38SStefano Zampini   }
10563425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10573425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10583425bc38SStefano Zampini   /* final change of basis if needed
10593425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
10603308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
10613425bc38SStefano Zampini   PetscFunctionReturn(0);
10623425bc38SStefano Zampini }
10631e6b0712SBarry Smith 
10643425bc38SStefano Zampini #undef __FUNCT__
10653425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
10663425bc38SStefano Zampini /*@
10673425bc38SStefano Zampini  PCBDDCMatFETIDPGetSolution - Get Solution for FETIDP linear system.
10683425bc38SStefano Zampini 
10693425bc38SStefano Zampini    Collective
10703425bc38SStefano Zampini 
10713425bc38SStefano Zampini    Input Parameters:
10723425bc38SStefano Zampini +  fetidp_mat        - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators
10733425bc38SStefano Zampini +  fetidp_flux_sol - the solution of the FETIDP linear system
10743425bc38SStefano Zampini 
10753425bc38SStefano Zampini    Output Parameters:
10763425bc38SStefano Zampini +  standard_sol      - the solution on the global domain
10773425bc38SStefano Zampini 
10783425bc38SStefano Zampini    Level: developer
10793425bc38SStefano Zampini 
10803425bc38SStefano Zampini    Notes:
10813425bc38SStefano Zampini 
10823425bc38SStefano Zampini .seealso: PCBDDC
10833425bc38SStefano Zampini @*/
10843425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
10853425bc38SStefano Zampini {
1086674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10873425bc38SStefano Zampini   PetscErrorCode ierr;
10883425bc38SStefano Zampini 
10893425bc38SStefano Zampini   PetscFunctionBegin;
10903425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10913425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
10923425bc38SStefano Zampini   PetscFunctionReturn(0);
10933425bc38SStefano Zampini }
10943425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
10951e6b0712SBarry Smith 
1096f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1097f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1098f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1099f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1100674ae819SStefano Zampini 
11013425bc38SStefano Zampini #undef __FUNCT__
11023425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
11033425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
11043425bc38SStefano Zampini {
1105674ae819SStefano Zampini 
1106674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
11073425bc38SStefano Zampini   Mat            newmat;
1108674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
11093425bc38SStefano Zampini   PC             newpc;
1110ce94432eSBarry Smith   MPI_Comm       comm;
11113425bc38SStefano Zampini   PetscErrorCode ierr;
11123425bc38SStefano Zampini 
11133425bc38SStefano Zampini   PetscFunctionBegin;
1114ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
11153425bc38SStefano Zampini   /* FETIDP linear matrix */
11163425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
11173425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
11183425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
11193425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
11203425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
11213425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
11223425bc38SStefano Zampini   /* FETIDP preconditioner */
11233425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
11243425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
11253425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
11263425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
11273425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
11283425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
11293425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
11303425bc38SStefano Zampini   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
11313425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
11323425bc38SStefano Zampini   /* return pointers for objects created */
11333425bc38SStefano Zampini   *fetidp_mat=newmat;
11343425bc38SStefano Zampini   *fetidp_pc=newpc;
11353425bc38SStefano Zampini   PetscFunctionReturn(0);
11363425bc38SStefano Zampini }
11371e6b0712SBarry Smith 
11383425bc38SStefano Zampini #undef __FUNCT__
11393425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
11403425bc38SStefano Zampini /*@
11413425bc38SStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP.
11423425bc38SStefano Zampini 
11433425bc38SStefano Zampini    Collective
11443425bc38SStefano Zampini 
11453425bc38SStefano Zampini    Input Parameters:
11463425bc38SStefano Zampini +  pc - the BDDC preconditioning context (setup must be already called)
11473425bc38SStefano Zampini 
11483425bc38SStefano Zampini    Level: developer
11493425bc38SStefano Zampini 
11503425bc38SStefano Zampini    Notes:
11513425bc38SStefano Zampini 
11523425bc38SStefano Zampini .seealso: PCBDDC
11533425bc38SStefano Zampini @*/
11543425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
11553425bc38SStefano Zampini {
11563425bc38SStefano Zampini   PetscErrorCode ierr;
11573425bc38SStefano Zampini 
11583425bc38SStefano Zampini   PetscFunctionBegin;
11593425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11603425bc38SStefano Zampini   if (pc->setupcalled) {
1161516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1162f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
11633425bc38SStefano Zampini   PetscFunctionReturn(0);
11643425bc38SStefano Zampini }
11650c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1166da1bb401SStefano Zampini /*MC
1167da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
11680c7d97c5SJed Brown 
1169da1bb401SStefano Zampini    Options Database Keys:
1170da1bb401SStefano Zampini .    -pcbddc ??? -
1171da1bb401SStefano Zampini 
1172da1bb401SStefano Zampini    Level: intermediate
1173da1bb401SStefano Zampini 
1174da1bb401SStefano Zampini    Notes: The matrix used with this preconditioner must be of type MATIS
1175da1bb401SStefano Zampini 
1176da1bb401SStefano Zampini           Unlike more 'conventional' interface preconditioners, this iterates over ALL the
1177da1bb401SStefano Zampini           degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers
1178da1bb401SStefano Zampini           on the subdomains).
1179da1bb401SStefano Zampini 
1180da1bb401SStefano Zampini           Options for the coarse grid preconditioner can be set with -
1181da1bb401SStefano Zampini           Options for the Dirichlet subproblem can be set with -
1182da1bb401SStefano Zampini           Options for the Neumann subproblem can be set with -
1183da1bb401SStefano Zampini 
1184da1bb401SStefano Zampini    Contributed by Stefano Zampini
1185da1bb401SStefano Zampini 
1186da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1187da1bb401SStefano Zampini M*/
1188b2573a8aSBarry Smith 
1189da1bb401SStefano Zampini #undef __FUNCT__
1190da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
11918cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1192da1bb401SStefano Zampini {
1193da1bb401SStefano Zampini   PetscErrorCode      ierr;
1194da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1195da1bb401SStefano Zampini 
1196da1bb401SStefano Zampini   PetscFunctionBegin;
1197da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1198da1bb401SStefano Zampini   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1199da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1200da1bb401SStefano Zampini 
1201da1bb401SStefano Zampini   /* create PCIS data structure */
1202da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1203da1bb401SStefano Zampini 
1204*47d04d0dSStefano Zampini   /* BDDC customization */
1205*47d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
1206*47d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
1207*47d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
1208*47d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
1209*47d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
1210*47d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
1211*47d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
1212*47d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
1213*47d04d0dSStefano Zampini 
1214674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
12150bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
12163972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1217534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1218534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1219534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1220da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1221da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1222da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1223da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1224da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
122515aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
122615aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1227da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1228da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1229da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1230da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1231da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1232da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1233da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1234da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1235da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1236da1bb401SStefano Zampini   pcbddc->local_primal_indices       = 0;
1237da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
1238da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
1239da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
1240da1bb401SStefano Zampini   pcbddc->local_primal_sizes         = 0;
1241da1bb401SStefano Zampini   pcbddc->local_primal_displacements = 0;
1242b76ba322SStefano Zampini   pcbddc->use_exact_dirichlet        = PETSC_TRUE;
1243*47d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
1244*47d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
12454fad6a16SStefano Zampini   pcbddc->current_level              = 0;
12464fad6a16SStefano Zampini   pcbddc->max_levels                 = 1;
1247674ae819SStefano Zampini   pcbddc->replicated_local_primal_indices = 0;
1248674ae819SStefano Zampini   pcbddc->replicated_local_primal_values  = 0;
1249da1bb401SStefano Zampini 
1250674ae819SStefano Zampini   /* create local graph structure */
1251674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1252674ae819SStefano Zampini 
1253674ae819SStefano Zampini   /* scaling */
1254674ae819SStefano Zampini   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1255674ae819SStefano Zampini   pcbddc->work_scaling               = 0;
1256da1bb401SStefano Zampini 
1257da1bb401SStefano Zampini   /* function pointers */
1258da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
1259da1bb401SStefano Zampini   pc->ops->applytranspose      = 0;
1260da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1261da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1262da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1263da1bb401SStefano Zampini   pc->ops->view                = 0;
1264da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1265da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1266da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1267534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1268534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1269da1bb401SStefano Zampini 
1270da1bb401SStefano Zampini   /* composing function */
1271674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1272bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
1273bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetMaxLevels_C",PCBDDCSetMaxLevels_BDDC);CHKERRQ(ierr);
1274bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1275bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1276bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1277bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1278bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1279bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseProblemType_C",PCBDDCSetCoarseProblemType_BDDC);CHKERRQ(ierr);
1280bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1281bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1282bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1283bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1284bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1285da1bb401SStefano Zampini   PetscFunctionReturn(0);
1286da1bb401SStefano Zampini }
12873425bc38SStefano Zampini 
1288