xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 785d12438418abfe1dc590cf0e72c5e5ca69a49b)
153cdbc3dSStefano Zampini /* TODOLIST
2eb97c9d2SStefano Zampini 
3eb97c9d2SStefano Zampini    ConstraintsSetup
4eb97c9d2SStefano Zampini    - assure same constraints between neighbours by sorting vals by global index before SVD!
5eb97c9d2SStefano Zampini    - tolerances for constraints as an option (take care of single precision!)
6eb97c9d2SStefano Zampini    - Allow different constraints customizations among different linear solves (requires also reset/destroy of ksp_R and coarse_ksp)
7eb97c9d2SStefano Zampini    - MAT_IGNORE_ZERO_ENTRIES for Constraints Matrix
8eb97c9d2SStefano Zampini 
9eb97c9d2SStefano Zampini    Solvers
10eb97c9d2SStefano Zampini    - Try to reduce the work when reusing the solvers
11eb97c9d2SStefano Zampini    - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers)
12eb97c9d2SStefano Zampini    - reuse already allocated coarse matrix if possible
13eb97c9d2SStefano Zampini    - Propagate ksp prefixes for solvers to mat objects?
14eb97c9d2SStefano Zampini    - Propagate nearnullspace info among levels
15eb97c9d2SStefano Zampini 
16eb97c9d2SStefano Zampini    User interface
17eb97c9d2SStefano Zampini    - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access)
18eb97c9d2SStefano Zampini    - Provide PCApplyTranpose_BDDC
19eb97c9d2SStefano Zampini    - DofSplitting and DM attached to pc?
20eb97c9d2SStefano Zampini 
21eb97c9d2SStefano Zampini    Debugging output
22eb97c9d2SStefano Zampini    - Better management of verbosity levels of debugging output
23eb97c9d2SStefano Zampini    - Crashes on some architecture -> call SynchronizedAllow before every SynchronizedPrintf
24eb97c9d2SStefano Zampini 
25eb97c9d2SStefano Zampini    Build
26eb97c9d2SStefano Zampini    - make runexe59
27eb97c9d2SStefano Zampini 
28eb97c9d2SStefano Zampini    Extra
29eb97c9d2SStefano Zampini    - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
30eb97c9d2SStefano Zampini    - Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver?
31eb97c9d2SStefano Zampini    - add support for computing h,H and related using coordinates?
32c0b83709SStefano Zampini    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
33eb97c9d2SStefano Zampini    - Better management in PCIS code
34eb97c9d2SStefano Zampini    - BDDC with MG framework?
35eb97c9d2SStefano Zampini 
36eb97c9d2SStefano Zampini    FETIDP
37eb97c9d2SStefano Zampini    - Move FETIDP code to its own classes
38eb97c9d2SStefano Zampini 
39eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
40eb97c9d2SStefano Zampini    - Provide general case for subassembling
41eb97c9d2SStefano Zampini    - Preallocation routines in MatConvert_IS_AIJ
42eb97c9d2SStefano Zampini 
4353cdbc3dSStefano Zampini */
440c7d97c5SJed Brown 
4553cdbc3dSStefano Zampini /* ----------------------------------------------------------------------------------------------------------------------------------------------
460c7d97c5SJed Brown    Implementation of BDDC preconditioner based on:
470c7d97c5SJed Brown    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
4853cdbc3dSStefano Zampini    ---------------------------------------------------------------------------------------------------------------------------------------------- */
4953cdbc3dSStefano Zampini 
50674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/  /* includes for fortran wrappers */
51674ae819SStefano Zampini #include "bddcprivate.h"
523b03a366Sstefano_zampini #include <petscblaslapack.h>
53674ae819SStefano Zampini 
540c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
550c7d97c5SJed Brown #undef __FUNCT__
560c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
570c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc)
580c7d97c5SJed Brown {
590c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
600c7d97c5SJed Brown   PetscErrorCode ierr;
610c7d97c5SJed Brown 
620c7d97c5SJed Brown   PetscFunctionBegin;
630c7d97c5SJed Brown   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
648eeda7d8SStefano Zampini   /* Verbose debugging */
658eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
668eeda7d8SStefano Zampini   /* Primal space cumstomization */
678eeda7d8SStefano 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);
688eeda7d8SStefano 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);
698eeda7d8SStefano 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);
708eeda7d8SStefano Zampini   /* Change of basis */
718eeda7d8SStefano 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);
728eeda7d8SStefano 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);
73674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
74674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
75674ae819SStefano Zampini   }
768eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
778eeda7d8SStefano 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);
780298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
792b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
80674ae819SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
810c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
820c7d97c5SJed Brown   PetscFunctionReturn(0);
830c7d97c5SJed Brown }
840c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
85674ae819SStefano Zampini #undef __FUNCT__
86674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
87674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
88674ae819SStefano Zampini {
89674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
90674ae819SStefano Zampini   PetscErrorCode ierr;
911e6b0712SBarry Smith 
92674ae819SStefano Zampini   PetscFunctionBegin;
93674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
94674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
95674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
96674ae819SStefano Zampini   PetscFunctionReturn(0);
97674ae819SStefano Zampini }
98674ae819SStefano Zampini #undef __FUNCT__
99674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
100674ae819SStefano Zampini /*@
10128509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
102674ae819SStefano Zampini 
103674ae819SStefano Zampini    Not collective
104674ae819SStefano Zampini 
105674ae819SStefano Zampini    Input Parameters:
106674ae819SStefano Zampini +  pc - the preconditioning context
10728509bceSStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering
108674ae819SStefano Zampini 
109674ae819SStefano Zampini    Level: intermediate
110674ae819SStefano Zampini 
111674ae819SStefano Zampini    Notes:
112674ae819SStefano Zampini 
113674ae819SStefano Zampini .seealso: PCBDDC
114674ae819SStefano Zampini @*/
115674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
116674ae819SStefano Zampini {
117674ae819SStefano Zampini   PetscErrorCode ierr;
118674ae819SStefano Zampini 
119674ae819SStefano Zampini   PetscFunctionBegin;
120674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
121674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
122674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
123674ae819SStefano Zampini   PetscFunctionReturn(0);
124674ae819SStefano Zampini }
125674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1260c7d97c5SJed Brown #undef __FUNCT__
1274fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1284fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1294fad6a16SStefano Zampini {
1304fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1314fad6a16SStefano Zampini 
1324fad6a16SStefano Zampini   PetscFunctionBegin;
1334fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1344fad6a16SStefano Zampini   PetscFunctionReturn(0);
1354fad6a16SStefano Zampini }
1361e6b0712SBarry Smith 
1374fad6a16SStefano Zampini #undef __FUNCT__
1384fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1394fad6a16SStefano Zampini /*@
14028509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
1414fad6a16SStefano Zampini 
1424fad6a16SStefano Zampini    Logically collective on PC
1434fad6a16SStefano Zampini 
1444fad6a16SStefano Zampini    Input Parameters:
1454fad6a16SStefano Zampini +  pc - the preconditioning context
14628509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
1474fad6a16SStefano Zampini 
14828509bceSStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
1494fad6a16SStefano Zampini 
1504fad6a16SStefano Zampini    Level: intermediate
1514fad6a16SStefano Zampini 
1524fad6a16SStefano Zampini    Notes:
1534fad6a16SStefano Zampini 
1544fad6a16SStefano Zampini .seealso: PCBDDC
1554fad6a16SStefano Zampini @*/
1564fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
1574fad6a16SStefano Zampini {
1584fad6a16SStefano Zampini   PetscErrorCode ierr;
1594fad6a16SStefano Zampini 
1604fad6a16SStefano Zampini   PetscFunctionBegin;
1614fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1622b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
1634fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
1644fad6a16SStefano Zampini   PetscFunctionReturn(0);
1654fad6a16SStefano Zampini }
1662b510759SStefano Zampini 
167b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
1682b510759SStefano Zampini #undef __FUNCT__
169b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
170b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
171b8ffe317SStefano Zampini {
172b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
173b8ffe317SStefano Zampini 
174b8ffe317SStefano Zampini   PetscFunctionBegin;
17585c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
176b8ffe317SStefano Zampini   PetscFunctionReturn(0);
177b8ffe317SStefano Zampini }
178b8ffe317SStefano Zampini 
179b8ffe317SStefano Zampini #undef __FUNCT__
180b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
181b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
1822b510759SStefano Zampini {
1832b510759SStefano Zampini   PetscErrorCode ierr;
1842b510759SStefano Zampini 
1852b510759SStefano Zampini   PetscFunctionBegin;
1862b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
187b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
188b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1892b510759SStefano Zampini   PetscFunctionReturn(0);
1902b510759SStefano Zampini }
1911e6b0712SBarry Smith 
1924fad6a16SStefano Zampini #undef __FUNCT__
1932b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
1942b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
1954fad6a16SStefano Zampini {
1964fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1974fad6a16SStefano Zampini 
1984fad6a16SStefano Zampini   PetscFunctionBegin;
1992b510759SStefano Zampini   pcbddc->current_level = level;
2004fad6a16SStefano Zampini   PetscFunctionReturn(0);
2014fad6a16SStefano Zampini }
2021e6b0712SBarry Smith 
2034fad6a16SStefano Zampini #undef __FUNCT__
204b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
205b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
206b8ffe317SStefano Zampini {
207b8ffe317SStefano Zampini   PetscErrorCode ierr;
208b8ffe317SStefano Zampini 
209b8ffe317SStefano Zampini   PetscFunctionBegin;
210b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
211b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
212b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
213b8ffe317SStefano Zampini   PetscFunctionReturn(0);
214b8ffe317SStefano Zampini }
215b8ffe317SStefano Zampini 
216b8ffe317SStefano Zampini #undef __FUNCT__
2172b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
2182b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
2192b510759SStefano Zampini {
2202b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2212b510759SStefano Zampini 
2222b510759SStefano Zampini   PetscFunctionBegin;
2232b510759SStefano Zampini   pcbddc->max_levels = levels;
2242b510759SStefano Zampini   PetscFunctionReturn(0);
2252b510759SStefano Zampini }
2262b510759SStefano Zampini 
2272b510759SStefano Zampini #undef __FUNCT__
2282b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2294fad6a16SStefano Zampini /*@
23028509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
2314fad6a16SStefano Zampini 
2324fad6a16SStefano Zampini    Logically collective on PC
2334fad6a16SStefano Zampini 
2344fad6a16SStefano Zampini    Input Parameters:
2354fad6a16SStefano Zampini +  pc - the preconditioning context
23628509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
2374fad6a16SStefano Zampini 
23828509bceSStefano Zampini    Default value is 0, i.e. traditional one-level BDDC
2394fad6a16SStefano Zampini 
2404fad6a16SStefano Zampini    Level: intermediate
2414fad6a16SStefano Zampini 
2424fad6a16SStefano Zampini    Notes:
2434fad6a16SStefano Zampini 
2444fad6a16SStefano Zampini .seealso: PCBDDC
2454fad6a16SStefano Zampini @*/
2462b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
2474fad6a16SStefano Zampini {
2484fad6a16SStefano Zampini   PetscErrorCode ierr;
2494fad6a16SStefano Zampini 
2504fad6a16SStefano Zampini   PetscFunctionBegin;
2514fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2522b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
2532b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
2544fad6a16SStefano Zampini   PetscFunctionReturn(0);
2554fad6a16SStefano Zampini }
2564fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
2571e6b0712SBarry Smith 
2584fad6a16SStefano Zampini #undef __FUNCT__
2590bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
2600bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
2610bdf917eSStefano Zampini {
2620bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2630bdf917eSStefano Zampini   PetscErrorCode ierr;
2640bdf917eSStefano Zampini 
2650bdf917eSStefano Zampini   PetscFunctionBegin;
2660bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
2670bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
2680bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
2690bdf917eSStefano Zampini   PetscFunctionReturn(0);
2700bdf917eSStefano Zampini }
2711e6b0712SBarry Smith 
2720bdf917eSStefano Zampini #undef __FUNCT__
2730bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
2740bdf917eSStefano Zampini /*@
27528509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
2760bdf917eSStefano Zampini 
2770bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
2780bdf917eSStefano Zampini 
2790bdf917eSStefano Zampini    Input Parameters:
2800bdf917eSStefano Zampini +  pc - the preconditioning context
28128509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
2820bdf917eSStefano Zampini 
2830bdf917eSStefano Zampini    Level: intermediate
2840bdf917eSStefano Zampini 
2850bdf917eSStefano Zampini    Notes:
2860bdf917eSStefano Zampini 
2870bdf917eSStefano Zampini .seealso: PCBDDC
2880bdf917eSStefano Zampini @*/
2890bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
2900bdf917eSStefano Zampini {
2910bdf917eSStefano Zampini   PetscErrorCode ierr;
2920bdf917eSStefano Zampini 
2930bdf917eSStefano Zampini   PetscFunctionBegin;
2940bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
295674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
2962b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
2970bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
2980bdf917eSStefano Zampini   PetscFunctionReturn(0);
2990bdf917eSStefano Zampini }
3000bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
3011e6b0712SBarry Smith 
3020bdf917eSStefano Zampini #undef __FUNCT__
303*785d1243SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
304*785d1243SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
305*785d1243SStefano Zampini {
306*785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
307*785d1243SStefano Zampini   PetscErrorCode ierr;
308*785d1243SStefano Zampini 
309*785d1243SStefano Zampini   PetscFunctionBegin;
310*785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
311*785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
312*785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
313*785d1243SStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
314*785d1243SStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
315*785d1243SStefano Zampini   PetscFunctionReturn(0);
316*785d1243SStefano Zampini }
317*785d1243SStefano Zampini 
318*785d1243SStefano Zampini #undef __FUNCT__
319*785d1243SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
320*785d1243SStefano Zampini /*@
321*785d1243SStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
322*785d1243SStefano Zampini 
323*785d1243SStefano Zampini    Collective
324*785d1243SStefano Zampini 
325*785d1243SStefano Zampini    Input Parameters:
326*785d1243SStefano Zampini +  pc - the preconditioning context
327*785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
328*785d1243SStefano Zampini 
329*785d1243SStefano Zampini    Level: intermediate
330*785d1243SStefano Zampini 
331*785d1243SStefano Zampini    Notes: Any process can list any global node
332*785d1243SStefano Zampini 
333*785d1243SStefano Zampini .seealso: PCBDDC
334*785d1243SStefano Zampini @*/
335*785d1243SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
336*785d1243SStefano Zampini {
337*785d1243SStefano Zampini   PetscErrorCode ierr;
338*785d1243SStefano Zampini 
339*785d1243SStefano Zampini   PetscFunctionBegin;
340*785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
341*785d1243SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
342*785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
343*785d1243SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
344*785d1243SStefano Zampini   PetscFunctionReturn(0);
345*785d1243SStefano Zampini }
346*785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
347*785d1243SStefano Zampini 
348*785d1243SStefano Zampini #undef __FUNCT__
34982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
35082ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
3513b03a366Sstefano_zampini {
3523b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3533b03a366Sstefano_zampini   PetscErrorCode ierr;
3543b03a366Sstefano_zampini 
3553b03a366Sstefano_zampini   PetscFunctionBegin;
356*785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
357*785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3583b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
35936e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
360*785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
3613b03a366Sstefano_zampini   PetscFunctionReturn(0);
3623b03a366Sstefano_zampini }
3631e6b0712SBarry Smith 
3643b03a366Sstefano_zampini #undef __FUNCT__
36582ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
3663b03a366Sstefano_zampini /*@
36782ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
3683b03a366Sstefano_zampini 
369*785d1243SStefano Zampini    Collective
3703b03a366Sstefano_zampini 
3713b03a366Sstefano_zampini    Input Parameters:
3723b03a366Sstefano_zampini +  pc - the preconditioning context
37382ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
3743b03a366Sstefano_zampini 
3753b03a366Sstefano_zampini    Level: intermediate
3763b03a366Sstefano_zampini 
3773b03a366Sstefano_zampini    Notes:
3783b03a366Sstefano_zampini 
3793b03a366Sstefano_zampini .seealso: PCBDDC
3803b03a366Sstefano_zampini @*/
38182ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
3823b03a366Sstefano_zampini {
3833b03a366Sstefano_zampini   PetscErrorCode ierr;
3843b03a366Sstefano_zampini 
3853b03a366Sstefano_zampini   PetscFunctionBegin;
3863b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
387674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
38882ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
38982ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
3903b03a366Sstefano_zampini   PetscFunctionReturn(0);
3913b03a366Sstefano_zampini }
3923b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
3931e6b0712SBarry Smith 
3943b03a366Sstefano_zampini #undef __FUNCT__
395*785d1243SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
396*785d1243SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
397*785d1243SStefano Zampini {
398*785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
399*785d1243SStefano Zampini   PetscErrorCode ierr;
400*785d1243SStefano Zampini 
401*785d1243SStefano Zampini   PetscFunctionBegin;
402*785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
403*785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
404*785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
405*785d1243SStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
406*785d1243SStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
407*785d1243SStefano Zampini   PetscFunctionReturn(0);
408*785d1243SStefano Zampini }
409*785d1243SStefano Zampini 
410*785d1243SStefano Zampini #undef __FUNCT__
411*785d1243SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
412*785d1243SStefano Zampini /*@
413*785d1243SStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
414*785d1243SStefano Zampini 
415*785d1243SStefano Zampini    Collective
416*785d1243SStefano Zampini 
417*785d1243SStefano Zampini    Input Parameters:
418*785d1243SStefano Zampini +  pc - the preconditioning context
419*785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
420*785d1243SStefano Zampini 
421*785d1243SStefano Zampini    Level: intermediate
422*785d1243SStefano Zampini 
423*785d1243SStefano Zampini    Notes: Any process can list any global node
424*785d1243SStefano Zampini 
425*785d1243SStefano Zampini .seealso: PCBDDC
426*785d1243SStefano Zampini @*/
427*785d1243SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
428*785d1243SStefano Zampini {
429*785d1243SStefano Zampini   PetscErrorCode ierr;
430*785d1243SStefano Zampini 
431*785d1243SStefano Zampini   PetscFunctionBegin;
432*785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
433*785d1243SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
434*785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
435*785d1243SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
436*785d1243SStefano Zampini   PetscFunctionReturn(0);
437*785d1243SStefano Zampini }
438*785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
439*785d1243SStefano Zampini 
440*785d1243SStefano Zampini #undef __FUNCT__
44182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
44282ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
4430c7d97c5SJed Brown {
4440c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
44553cdbc3dSStefano Zampini   PetscErrorCode ierr;
4460c7d97c5SJed Brown 
4470c7d97c5SJed Brown   PetscFunctionBegin;
448*785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
449*785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
45053cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
45136e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
452*785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
4530c7d97c5SJed Brown   PetscFunctionReturn(0);
4540c7d97c5SJed Brown }
4551e6b0712SBarry Smith 
4560c7d97c5SJed Brown #undef __FUNCT__
45782ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
45857527edcSJed Brown /*@
45982ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
46057527edcSJed Brown 
461*785d1243SStefano Zampini    Collective
46257527edcSJed Brown 
46357527edcSJed Brown    Input Parameters:
46457527edcSJed Brown +  pc - the preconditioning context
46582ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
46657527edcSJed Brown 
46757527edcSJed Brown    Level: intermediate
46857527edcSJed Brown 
46957527edcSJed Brown    Notes:
47057527edcSJed Brown 
47157527edcSJed Brown .seealso: PCBDDC
47257527edcSJed Brown @*/
47382ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
4740c7d97c5SJed Brown {
4750c7d97c5SJed Brown   PetscErrorCode ierr;
4760c7d97c5SJed Brown 
4770c7d97c5SJed Brown   PetscFunctionBegin;
4780c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
479674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
48082ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
48182ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
48253cdbc3dSStefano Zampini   PetscFunctionReturn(0);
48353cdbc3dSStefano Zampini }
48453cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
4851e6b0712SBarry Smith 
48653cdbc3dSStefano Zampini #undef __FUNCT__
487*785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
488*785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
489da1bb401SStefano Zampini {
490da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
491da1bb401SStefano Zampini 
492da1bb401SStefano Zampini   PetscFunctionBegin;
493da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
494da1bb401SStefano Zampini   PetscFunctionReturn(0);
495da1bb401SStefano Zampini }
4961e6b0712SBarry Smith 
497da1bb401SStefano Zampini #undef __FUNCT__
498*785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
499*785d1243SStefano Zampini /*@
500*785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
501*785d1243SStefano Zampini 
502*785d1243SStefano Zampini    Collective
503*785d1243SStefano Zampini 
504*785d1243SStefano Zampini    Input Parameters:
505*785d1243SStefano Zampini .  pc - the preconditioning context
506*785d1243SStefano Zampini 
507*785d1243SStefano Zampini    Output Parameters:
508*785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
509*785d1243SStefano Zampini 
510*785d1243SStefano Zampini    Level: intermediate
511*785d1243SStefano Zampini 
512*785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
513*785d1243SStefano Zampini 
514*785d1243SStefano Zampini .seealso: PCBDDC
515*785d1243SStefano Zampini @*/
516*785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
517*785d1243SStefano Zampini {
518*785d1243SStefano Zampini   PetscErrorCode ierr;
519*785d1243SStefano Zampini 
520*785d1243SStefano Zampini   PetscFunctionBegin;
521*785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
522*785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
523*785d1243SStefano Zampini   PetscFunctionReturn(0);
524*785d1243SStefano Zampini }
525*785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
526*785d1243SStefano Zampini 
527*785d1243SStefano Zampini #undef __FUNCT__
528*785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
529*785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
530*785d1243SStefano Zampini {
531*785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
532*785d1243SStefano Zampini 
533*785d1243SStefano Zampini   PetscFunctionBegin;
534*785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
535*785d1243SStefano Zampini   PetscFunctionReturn(0);
536*785d1243SStefano Zampini }
537*785d1243SStefano Zampini 
538*785d1243SStefano Zampini #undef __FUNCT__
53982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
540da1bb401SStefano Zampini /*@
54182ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
542da1bb401SStefano Zampini 
543*785d1243SStefano Zampini    Collective
544da1bb401SStefano Zampini 
545da1bb401SStefano Zampini    Input Parameters:
54628509bceSStefano Zampini .  pc - the preconditioning context
547da1bb401SStefano Zampini 
548da1bb401SStefano Zampini    Output Parameters:
54928509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
550da1bb401SStefano Zampini 
551da1bb401SStefano Zampini    Level: intermediate
552da1bb401SStefano Zampini 
553da1bb401SStefano Zampini    Notes:
554da1bb401SStefano Zampini 
555da1bb401SStefano Zampini .seealso: PCBDDC
556da1bb401SStefano Zampini @*/
55782ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
558da1bb401SStefano Zampini {
559da1bb401SStefano Zampini   PetscErrorCode ierr;
560da1bb401SStefano Zampini 
561da1bb401SStefano Zampini   PetscFunctionBegin;
562da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
56382ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
564da1bb401SStefano Zampini   PetscFunctionReturn(0);
565da1bb401SStefano Zampini }
566da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
5671e6b0712SBarry Smith 
568da1bb401SStefano Zampini #undef __FUNCT__
569*785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
570*785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
57153cdbc3dSStefano Zampini {
57253cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
57353cdbc3dSStefano Zampini 
57453cdbc3dSStefano Zampini   PetscFunctionBegin;
57553cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
57653cdbc3dSStefano Zampini   PetscFunctionReturn(0);
57753cdbc3dSStefano Zampini }
5781e6b0712SBarry Smith 
57953cdbc3dSStefano Zampini #undef __FUNCT__
580*785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
581*785d1243SStefano Zampini /*@
582*785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
583*785d1243SStefano Zampini 
584*785d1243SStefano Zampini    Collective
585*785d1243SStefano Zampini 
586*785d1243SStefano Zampini    Input Parameters:
587*785d1243SStefano Zampini .  pc - the preconditioning context
588*785d1243SStefano Zampini 
589*785d1243SStefano Zampini    Output Parameters:
590*785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
591*785d1243SStefano Zampini 
592*785d1243SStefano Zampini    Level: intermediate
593*785d1243SStefano Zampini 
594*785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
595*785d1243SStefano Zampini 
596*785d1243SStefano Zampini .seealso: PCBDDC
597*785d1243SStefano Zampini @*/
598*785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
599*785d1243SStefano Zampini {
600*785d1243SStefano Zampini   PetscErrorCode ierr;
601*785d1243SStefano Zampini 
602*785d1243SStefano Zampini   PetscFunctionBegin;
603*785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
604*785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
605*785d1243SStefano Zampini   PetscFunctionReturn(0);
606*785d1243SStefano Zampini }
607*785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
608*785d1243SStefano Zampini 
609*785d1243SStefano Zampini #undef __FUNCT__
610*785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
611*785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
612*785d1243SStefano Zampini {
613*785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
614*785d1243SStefano Zampini 
615*785d1243SStefano Zampini   PetscFunctionBegin;
616*785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
617*785d1243SStefano Zampini   PetscFunctionReturn(0);
618*785d1243SStefano Zampini }
619*785d1243SStefano Zampini 
620*785d1243SStefano Zampini #undef __FUNCT__
62182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
62253cdbc3dSStefano Zampini /*@
62382ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
62453cdbc3dSStefano Zampini 
625*785d1243SStefano Zampini    Collective
62653cdbc3dSStefano Zampini 
62753cdbc3dSStefano Zampini    Input Parameters:
62828509bceSStefano Zampini .  pc - the preconditioning context
62953cdbc3dSStefano Zampini 
63053cdbc3dSStefano Zampini    Output Parameters:
63128509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
63253cdbc3dSStefano Zampini 
63353cdbc3dSStefano Zampini    Level: intermediate
63453cdbc3dSStefano Zampini 
63553cdbc3dSStefano Zampini    Notes:
63653cdbc3dSStefano Zampini 
63753cdbc3dSStefano Zampini .seealso: PCBDDC
63853cdbc3dSStefano Zampini @*/
63982ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
64053cdbc3dSStefano Zampini {
64153cdbc3dSStefano Zampini   PetscErrorCode ierr;
64253cdbc3dSStefano Zampini 
64353cdbc3dSStefano Zampini   PetscFunctionBegin;
64453cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
64582ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
6460c7d97c5SJed Brown   PetscFunctionReturn(0);
6470c7d97c5SJed Brown }
64836e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
6491e6b0712SBarry Smith 
65036e030ebSStefano Zampini #undef __FUNCT__
651da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
6521a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
65336e030ebSStefano Zampini {
65436e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
655da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
656da1bb401SStefano Zampini   PetscErrorCode ierr;
65736e030ebSStefano Zampini 
65836e030ebSStefano Zampini   PetscFunctionBegin;
659674ae819SStefano Zampini   /* free old CSR */
660674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
661fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
662674ae819SStefano Zampini   /* get CSR into graph structure */
663da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
664674ae819SStefano Zampini     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
665674ae819SStefano Zampini     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
666674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
667674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
668da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
6691a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
6701a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
671674ae819SStefano Zampini   }
672575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
67336e030ebSStefano Zampini   PetscFunctionReturn(0);
67436e030ebSStefano Zampini }
6751e6b0712SBarry Smith 
67636e030ebSStefano Zampini #undef __FUNCT__
677da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
67836e030ebSStefano Zampini /*@
67928509bceSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
68036e030ebSStefano Zampini 
68136e030ebSStefano Zampini    Not collective
68236e030ebSStefano Zampini 
68336e030ebSStefano Zampini    Input Parameters:
68436e030ebSStefano Zampini +  pc - the preconditioning context
68528509bceSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
68628509bceSStefano Zampini .  xadj, adjncy - the CSR graph
68728509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
68836e030ebSStefano Zampini 
68936e030ebSStefano Zampini    Level: intermediate
69036e030ebSStefano Zampini 
69136e030ebSStefano Zampini    Notes:
69236e030ebSStefano Zampini 
69328509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
69436e030ebSStefano Zampini @*/
6951a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
69636e030ebSStefano Zampini {
697575ad6abSStefano Zampini   void (*f)(void) = 0;
69836e030ebSStefano Zampini   PetscErrorCode ierr;
69936e030ebSStefano Zampini 
70036e030ebSStefano Zampini   PetscFunctionBegin;
70136e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
702674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
703674ae819SStefano Zampini   PetscValidIntPointer(xadj,4);
704674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
705674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
706da1bb401SStefano Zampini   }
70736e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
708575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
709575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
710575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
711575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
712575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
71336e030ebSStefano Zampini   }
71436e030ebSStefano Zampini   PetscFunctionReturn(0);
71536e030ebSStefano Zampini }
7169c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
7171e6b0712SBarry Smith 
7189c0446d6SStefano Zampini #undef __FUNCT__
7199c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
7209c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
7219c0446d6SStefano Zampini {
7229c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
7239c0446d6SStefano Zampini   PetscInt i;
7249c0446d6SStefano Zampini   PetscErrorCode ierr;
7259c0446d6SStefano Zampini 
7269c0446d6SStefano Zampini   PetscFunctionBegin;
727da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
7289c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
7299c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
7309c0446d6SStefano Zampini   }
731d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
732da1bb401SStefano Zampini   /* allocate space then set */
7339c0446d6SStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
7349c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
735da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
736da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
7379c0446d6SStefano Zampini   }
7389c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
73960d44989SStefano Zampini   pcbddc->user_provided_isfordofs = PETSC_TRUE;
7409c0446d6SStefano Zampini   PetscFunctionReturn(0);
7419c0446d6SStefano Zampini }
7421e6b0712SBarry Smith 
7439c0446d6SStefano Zampini #undef __FUNCT__
7449c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
7459c0446d6SStefano Zampini /*@
74628509bceSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the local Neumann matrix
7479c0446d6SStefano Zampini 
7489c0446d6SStefano Zampini    Not collective
7499c0446d6SStefano Zampini 
7509c0446d6SStefano Zampini    Input Parameters:
7519c0446d6SStefano Zampini +  pc - the preconditioning context
75228509bceSStefano Zampini -  n_is - number of index sets defining the fields
75328509bceSStefano Zampini .  ISForDofs - array of IS describing the fields
7549c0446d6SStefano Zampini 
7559c0446d6SStefano Zampini    Level: intermediate
7569c0446d6SStefano Zampini 
7579c0446d6SStefano Zampini    Notes:
7589c0446d6SStefano Zampini 
7599c0446d6SStefano Zampini .seealso: PCBDDC
7609c0446d6SStefano Zampini @*/
7619c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
7629c0446d6SStefano Zampini {
7632b510759SStefano Zampini   PetscInt       i;
7649c0446d6SStefano Zampini   PetscErrorCode ierr;
7659c0446d6SStefano Zampini 
7669c0446d6SStefano Zampini   PetscFunctionBegin;
7679c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
7682b510759SStefano Zampini   for (i=0;i<n_is;i++) {
7692b510759SStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2);
7702b510759SStefano Zampini   }
7719c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
7729c0446d6SStefano Zampini   PetscFunctionReturn(0);
7739c0446d6SStefano Zampini }
774da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
775534831adSStefano Zampini #undef __FUNCT__
776534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
777534831adSStefano Zampini /* -------------------------------------------------------------------------- */
778534831adSStefano Zampini /*
779534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
780534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
7819c0446d6SStefano Zampini 
782534831adSStefano Zampini    Input Parameter:
783534831adSStefano Zampini +  pc - the preconditioner contex
784534831adSStefano Zampini 
785534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
786534831adSStefano Zampini 
787534831adSStefano Zampini    Notes:
788534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
789534831adSStefano Zampini    the user, but instead is called by KSPSolve().
790534831adSStefano Zampini */
791534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
792534831adSStefano Zampini {
793534831adSStefano Zampini   PetscErrorCode ierr;
794534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
795534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
796534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
797534831adSStefano Zampini   Mat            temp_mat;
7983972b0daSStefano Zampini   IS             dirIS;
7993972b0daSStefano Zampini   Vec            used_vec;
80082ba6b80SStefano Zampini   PetscBool      guess_nonzero;
801534831adSStefano Zampini 
802534831adSStefano Zampini   PetscFunctionBegin;
80385c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
80485c4d303SStefano Zampini   if (ksp) {
80585c4d303SStefano Zampini     PetscBool iscg;
80685c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
80785c4d303SStefano Zampini     if (!iscg) {
80885c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
80985c4d303SStefano Zampini     }
81085c4d303SStefano Zampini   }
81185c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
81262a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
81362a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
81462a6ff1dSStefano Zampini   }
81562a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
81662a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
81762a6ff1dSStefano Zampini   }
8183972b0daSStefano Zampini   if (x) {
8193972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
8203972b0daSStefano Zampini     used_vec = x;
8213972b0daSStefano Zampini   } else {
8223972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
8233972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
8243972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
8253972b0daSStefano Zampini   }
8263972b0daSStefano Zampini   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
8273972b0daSStefano Zampini   if (ksp) {
8283972b0daSStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
8293972b0daSStefano Zampini     if (!guess_nonzero) {
8303972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
8313972b0daSStefano Zampini     }
8323972b0daSStefano Zampini   }
8333308cffdSStefano Zampini 
8343972b0daSStefano Zampini   /* store the original rhs */
8353972b0daSStefano Zampini   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
8363972b0daSStefano Zampini 
8373972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
838*785d1243SStefano Zampini   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
83982ba6b80SStefano Zampini   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
84082ba6b80SStefano Zampini   if (rhs && dirIS) {
841*785d1243SStefano Zampini     PetscInt    dirsize,i,*is_indices;
842*785d1243SStefano Zampini     PetscScalar *array_x,*array_diagonal;
843*785d1243SStefano Zampini 
8443972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
8453972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
8463972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8473972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8483972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8493972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
85082ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
8513972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
8523972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
8533972b0daSStefano Zampini     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
8542fa5cd67SKarl Rupp     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
8553972b0daSStefano Zampini     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
8563972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
8573972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
8583972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8593972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
860b76ba322SStefano Zampini 
8613972b0daSStefano Zampini     /* remove the computed solution from the rhs */
8623972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
8633972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
8643972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
8653308cffdSStefano Zampini   }
866b76ba322SStefano Zampini 
867b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
8683972b0daSStefano Zampini   if (x) {
8693972b0daSStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
8703972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
87185c4d303SStefano Zampini     if (pcbddc->use_exact_dirichlet_trick) {
872b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
873b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
874b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
875b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
876b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
877b76ba322SStefano Zampini       if (ksp) {
878b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
879b76ba322SStefano Zampini       }
880b76ba322SStefano Zampini     }
8813972b0daSStefano Zampini   }
882b76ba322SStefano Zampini 
88392e3dcfbSStefano Zampini   /* prepare MatMult and rhs for solver */
884674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
885b76ba322SStefano Zampini     /* swap pointers for local matrices */
886b76ba322SStefano Zampini     temp_mat = matis->A;
887b76ba322SStefano Zampini     matis->A = pcbddc->local_mat;
888b76ba322SStefano Zampini     pcbddc->local_mat = temp_mat;
88992e3dcfbSStefano Zampini     if (rhs) {
890b76ba322SStefano Zampini       /* Get local rhs and apply transformation of basis */
891b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
892b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
893b76ba322SStefano Zampini       /* from original basis to modified basis */
894b76ba322SStefano Zampini       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
895b76ba322SStefano Zampini       /* put back modified values into the global vec using INSERT_VALUES copy mode */
896b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
897b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
898674ae819SStefano Zampini     }
89992e3dcfbSStefano Zampini   }
90092e3dcfbSStefano Zampini 
90192e3dcfbSStefano Zampini   /* remove nullspace if present */
9020bdf917eSStefano Zampini   if (ksp && pcbddc->NullSpace) {
903d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
904d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
905b76ba322SStefano Zampini   }
9060bdf917eSStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
907534831adSStefano Zampini   PetscFunctionReturn(0);
908534831adSStefano Zampini }
909534831adSStefano Zampini /* -------------------------------------------------------------------------- */
910534831adSStefano Zampini #undef __FUNCT__
911534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
912534831adSStefano Zampini /* -------------------------------------------------------------------------- */
913534831adSStefano Zampini /*
914534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
915534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
916534831adSStefano Zampini 
917534831adSStefano Zampini    Input Parameter:
918534831adSStefano Zampini +  pc - the preconditioner contex
919534831adSStefano Zampini 
920534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
921534831adSStefano Zampini 
922534831adSStefano Zampini    Notes:
923534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
924534831adSStefano Zampini    the user, but instead is called by KSPSolve().
925534831adSStefano Zampini */
926534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
927534831adSStefano Zampini {
928534831adSStefano Zampini   PetscErrorCode ierr;
929534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
930534831adSStefano Zampini   PC_IS          *pcis   = (PC_IS*)(pc->data);
931534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
932534831adSStefano Zampini   Mat            temp_mat;
933534831adSStefano Zampini 
934534831adSStefano Zampini   PetscFunctionBegin;
935674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
936534831adSStefano Zampini     /* swap pointers for local matrices */
937534831adSStefano Zampini     temp_mat = matis->A;
938534831adSStefano Zampini     matis->A = pcbddc->local_mat;
939534831adSStefano Zampini     pcbddc->local_mat = temp_mat;
9403425bc38SStefano Zampini   }
9413308cffdSStefano Zampini   if (pcbddc->use_change_of_basis && x) {
942534831adSStefano Zampini     /* Get Local boundary and apply transformation of basis to solution vector */
943534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
944534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
945534831adSStefano Zampini     /* from modified basis to original basis */
946534831adSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
947534831adSStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
948534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
949534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
950534831adSStefano Zampini   }
9513972b0daSStefano Zampini   /* add solution removed in presolve */
9523425bc38SStefano Zampini   if (x) {
9533425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
9543425bc38SStefano Zampini   }
955fb223d50SStefano Zampini   /* restore rhs to its original state */
956fb223d50SStefano Zampini   if (rhs) {
957fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
958fb223d50SStefano Zampini   }
959534831adSStefano Zampini   PetscFunctionReturn(0);
960534831adSStefano Zampini }
961534831adSStefano Zampini /* -------------------------------------------------------------------------- */
96253cdbc3dSStefano Zampini #undef __FUNCT__
96353cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
9640c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
9650c7d97c5SJed Brown /*
9660c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
9670c7d97c5SJed Brown                   by setting data structures and options.
9680c7d97c5SJed Brown 
9690c7d97c5SJed Brown    Input Parameter:
97053cdbc3dSStefano Zampini +  pc - the preconditioner context
9710c7d97c5SJed Brown 
9720c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
9730c7d97c5SJed Brown 
9740c7d97c5SJed Brown    Notes:
9750c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
9760c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
9770c7d97c5SJed Brown */
97853cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
9790c7d97c5SJed Brown {
9800c7d97c5SJed Brown   PetscErrorCode ierr;
9810c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
982674ae819SStefano Zampini   MatStructure   flag;
983674ae819SStefano Zampini   PetscBool      computeis,computetopography,computesolvers;
9840c7d97c5SJed Brown 
9850c7d97c5SJed Brown   PetscFunctionBegin;
986674ae819SStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */
987fcd409f5SStefano Zampini   /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */
9883b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
989fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
9903b03a366Sstefano_zampini   /* Get stdout for dbg */
991674ae819SStefano Zampini   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
992ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
993e269702eSStefano Zampini     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
9942b510759SStefano Zampini     if (pcbddc->current_level) {
9952b510759SStefano Zampini       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
9962b510759SStefano Zampini     }
997e269702eSStefano Zampini   }
998674ae819SStefano Zampini   /* first attempt to split work */
999674ae819SStefano Zampini   if (pc->setupcalled) {
1000674ae819SStefano Zampini     computeis = PETSC_FALSE;
1001674ae819SStefano Zampini     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
1002674ae819SStefano Zampini     if (flag == SAME_PRECONDITIONER) {
1003674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1004674ae819SStefano Zampini       computesolvers = PETSC_FALSE;
1005674ae819SStefano Zampini     } else if (flag == SAME_NONZERO_PATTERN) {
1006674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1007674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1008674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1009674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1010674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1011674ae819SStefano Zampini     }
1012674ae819SStefano Zampini   } else {
1013674ae819SStefano Zampini     computeis = PETSC_TRUE;
1014674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1015674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1016674ae819SStefano Zampini   }
1017fcd409f5SStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
1018674ae819SStefano Zampini   if (computeis) {
1019fcd409f5SStefano Zampini     /* HACK INTO PCIS */
1020fcd409f5SStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
1021fcd409f5SStefano Zampini     pcis->computesolvers = PETSC_FALSE;
1022674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
1023674ae819SStefano Zampini   }
1024674ae819SStefano Zampini   /* Analyze interface and set up local constraint and change of basis matrices */
1025674ae819SStefano Zampini   if (computetopography) {
1026674ae819SStefano Zampini     /* reset data */
1027674ae819SStefano Zampini     ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1028674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1029674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1030674ae819SStefano Zampini   }
1031674ae819SStefano Zampini   if (computesolvers) {
1032674ae819SStefano Zampini     /* reset data */
1033674ae819SStefano Zampini     ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
1034674ae819SStefano Zampini     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
103599cc7994SStefano Zampini     /* Create coarse and local stuffs */
103699cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1037674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
10380c7d97c5SJed Brown   }
10392b510759SStefano Zampini   if (pcbddc->dbg_flag && pcbddc->current_level) {
10402b510759SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
10412b510759SStefano Zampini   }
10420c7d97c5SJed Brown   PetscFunctionReturn(0);
10430c7d97c5SJed Brown }
10440c7d97c5SJed Brown 
10450c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
10460c7d97c5SJed Brown /*
10470c7d97c5SJed Brown    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
10480c7d97c5SJed Brown 
10490c7d97c5SJed Brown    Input Parameters:
10500c7d97c5SJed Brown .  pc - the preconditioner context
10510c7d97c5SJed Brown .  r - input vector (global)
10520c7d97c5SJed Brown 
10530c7d97c5SJed Brown    Output Parameter:
10540c7d97c5SJed Brown .  z - output vector (global)
10550c7d97c5SJed Brown 
10560c7d97c5SJed Brown    Application Interface Routine: PCApply()
10570c7d97c5SJed Brown  */
10580c7d97c5SJed Brown #undef __FUNCT__
10590c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
106053cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
10610c7d97c5SJed Brown {
10620c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
10630c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
10640c7d97c5SJed Brown   PetscErrorCode    ierr;
10653b03a366Sstefano_zampini   const PetscScalar one = 1.0;
10663b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
10672617d88aSStefano Zampini   const PetscScalar zero = 0.0;
10680c7d97c5SJed Brown 
10690c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
10700c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
10718eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
10720c7d97c5SJed Brown 
10730c7d97c5SJed Brown   PetscFunctionBegin;
107485c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
10750c7d97c5SJed Brown     /* First Dirichlet solve */
10760c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10770c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
107853cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
10790c7d97c5SJed Brown     /*
10800c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1081674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1082674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
10830c7d97c5SJed Brown     */
10840c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
10850c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
10868eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
10870c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
10880c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
10890c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10900c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1091674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1092b76ba322SStefano Zampini   } else {
10930bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1094b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1095674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1096b76ba322SStefano Zampini   }
1097b76ba322SStefano Zampini 
10982617d88aSStefano Zampini   /* Apply interface preconditioner
10992617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
11002617d88aSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
11012617d88aSStefano Zampini 
1102674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1103674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
11040c7d97c5SJed Brown 
11053b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
11060c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11070c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11080c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
11098eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
111053cdbc3dSStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr);
11110c7d97c5SJed Brown   ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr);
11128eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
11130c7d97c5SJed Brown   ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr);
11140c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11150c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11160c7d97c5SJed Brown   PetscFunctionReturn(0);
11170c7d97c5SJed Brown }
1118da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1119674ae819SStefano Zampini 
1120da1bb401SStefano Zampini #undef __FUNCT__
1121da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1122da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1123da1bb401SStefano Zampini {
1124da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1125da1bb401SStefano Zampini   PetscErrorCode ierr;
1126da1bb401SStefano Zampini 
1127da1bb401SStefano Zampini   PetscFunctionBegin;
1128da1bb401SStefano Zampini   /* free data created by PCIS */
1129da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1130674ae819SStefano Zampini   /* free BDDC custom data  */
1131674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1132674ae819SStefano Zampini   /* destroy objects related to topography */
1133674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1134674ae819SStefano Zampini   /* free allocated graph structure */
1135da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1136674ae819SStefano Zampini   /* free data for scaling operator */
1137674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1138674ae819SStefano Zampini   /* free solvers stuff */
1139674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
114033bc96a4SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
114133bc96a4SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
11429881197aSStefano Zampini   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
1143ac78edfcSStefano Zampini   ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
114462a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
114562a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
114662a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
11473425bc38SStefano Zampini   /* remove functions */
1148674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1149bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
11502b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1151b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
11522b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1153bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1154*785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
115582ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1156*785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
115782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1158*785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
115982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1160*785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1161*785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1162bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1163bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1164bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1165bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1166bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1167674ae819SStefano Zampini   /* Free the private data structure */
1168674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1169da1bb401SStefano Zampini   PetscFunctionReturn(0);
1170da1bb401SStefano Zampini }
11713425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
11721e6b0712SBarry Smith 
11733425bc38SStefano Zampini #undef __FUNCT__
11743425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
11753425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
11763425bc38SStefano Zampini {
1177674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
11783425bc38SStefano Zampini   PC_IS*         pcis;
11793425bc38SStefano Zampini   PC_BDDC*       pcbddc;
11803425bc38SStefano Zampini   PetscErrorCode ierr;
11810c7d97c5SJed Brown 
11823425bc38SStefano Zampini   PetscFunctionBegin;
11833425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
11843425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
11853425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
11863425bc38SStefano Zampini 
11873425bc38SStefano Zampini   /* change of basis for physical rhs if needed
11883425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
11893308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
11903425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
11913425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11923425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1193fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1194fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
11953425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11963425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1197674ae819SStefano Zampini   /* Apply partition of unity */
11983425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1199674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
12008eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
12013425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
12023425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
12033425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
12043425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
12053425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
12063425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12073425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1208674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
12093425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12103425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12113425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
12123425bc38SStefano Zampini   }
12133425bc38SStefano Zampini   /* BDDC rhs */
12143425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
12158eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
12163425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
12173425bc38SStefano Zampini   }
12183425bc38SStefano Zampini   /* apply BDDC */
12193425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
12203425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
12213425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
12223425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
12233425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12243425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12253425bc38SStefano Zampini   /* restore original rhs */
12263425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
12273425bc38SStefano Zampini   PetscFunctionReturn(0);
12283425bc38SStefano Zampini }
12291e6b0712SBarry Smith 
12303425bc38SStefano Zampini #undef __FUNCT__
12313425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
12323425bc38SStefano Zampini /*@
123328509bceSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
12343425bc38SStefano Zampini 
12353425bc38SStefano Zampini    Collective
12363425bc38SStefano Zampini 
12373425bc38SStefano Zampini    Input Parameters:
123828509bceSStefano Zampini +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
123928509bceSStefano Zampini .  standard_rhs - the right-hand side for your linear system
12403425bc38SStefano Zampini 
12413425bc38SStefano Zampini    Output Parameters:
124228509bceSStefano Zampini -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
12433425bc38SStefano Zampini 
12443425bc38SStefano Zampini    Level: developer
12453425bc38SStefano Zampini 
12463425bc38SStefano Zampini    Notes:
12473425bc38SStefano Zampini 
124828509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
12493425bc38SStefano Zampini @*/
12503425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
12513425bc38SStefano Zampini {
1252674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
12533425bc38SStefano Zampini   PetscErrorCode ierr;
12543425bc38SStefano Zampini 
12553425bc38SStefano Zampini   PetscFunctionBegin;
12563425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
12573425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
12583425bc38SStefano Zampini   PetscFunctionReturn(0);
12593425bc38SStefano Zampini }
12603425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
12611e6b0712SBarry Smith 
12623425bc38SStefano Zampini #undef __FUNCT__
12633425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
12643425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
12653425bc38SStefano Zampini {
1266674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
12673425bc38SStefano Zampini   PC_IS*         pcis;
12683425bc38SStefano Zampini   PC_BDDC*       pcbddc;
12693425bc38SStefano Zampini   PetscErrorCode ierr;
12703425bc38SStefano Zampini 
12713425bc38SStefano Zampini   PetscFunctionBegin;
12723425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
12733425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
12743425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
12753425bc38SStefano Zampini 
12763425bc38SStefano Zampini   /* apply B_delta^T */
12773425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12783425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12793425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
12803425bc38SStefano Zampini   /* compute rhs for BDDC application */
12813425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
12828eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
12833425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
12843425bc38SStefano Zampini   }
12853425bc38SStefano Zampini   /* apply BDDC */
12863425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
12873425bc38SStefano Zampini   /* put values into standard global vector */
12883425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12893425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12908eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
12913425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
12923425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
12933425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
12943425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
12953425bc38SStefano Zampini   }
12963425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12973425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12983425bc38SStefano Zampini   /* final change of basis if needed
12993425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
13003308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
13013425bc38SStefano Zampini   PetscFunctionReturn(0);
13023425bc38SStefano Zampini }
13031e6b0712SBarry Smith 
13043425bc38SStefano Zampini #undef __FUNCT__
13053425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
13063425bc38SStefano Zampini /*@
130728509bceSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
13083425bc38SStefano Zampini 
13093425bc38SStefano Zampini    Collective
13103425bc38SStefano Zampini 
13113425bc38SStefano Zampini    Input Parameters:
131228509bceSStefano Zampini +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
131328509bceSStefano Zampini .  fetidp_flux_sol - the solution of the FETIDP linear system
13143425bc38SStefano Zampini 
13153425bc38SStefano Zampini    Output Parameters:
131628509bceSStefano Zampini -  standard_sol      - the solution defined on the physical domain
13173425bc38SStefano Zampini 
13183425bc38SStefano Zampini    Level: developer
13193425bc38SStefano Zampini 
13203425bc38SStefano Zampini    Notes:
13213425bc38SStefano Zampini 
132228509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
13233425bc38SStefano Zampini @*/
13243425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
13253425bc38SStefano Zampini {
1326674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
13273425bc38SStefano Zampini   PetscErrorCode ierr;
13283425bc38SStefano Zampini 
13293425bc38SStefano Zampini   PetscFunctionBegin;
13303425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
13313425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
13323425bc38SStefano Zampini   PetscFunctionReturn(0);
13333425bc38SStefano Zampini }
13343425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
13351e6b0712SBarry Smith 
1336f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1337f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1338f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1339f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1340674ae819SStefano Zampini 
13413425bc38SStefano Zampini #undef __FUNCT__
13423425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
13433425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
13443425bc38SStefano Zampini {
1345674ae819SStefano Zampini 
1346674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
13473425bc38SStefano Zampini   Mat            newmat;
1348674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
13493425bc38SStefano Zampini   PC             newpc;
1350ce94432eSBarry Smith   MPI_Comm       comm;
13513425bc38SStefano Zampini   PetscErrorCode ierr;
13523425bc38SStefano Zampini 
13533425bc38SStefano Zampini   PetscFunctionBegin;
1354ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
13553425bc38SStefano Zampini   /* FETIDP linear matrix */
13563425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
13573425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
13583425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
13593425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
13603425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
13613425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
13623425bc38SStefano Zampini   /* FETIDP preconditioner */
13633425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
13643425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
13653425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
13663425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
13673425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
13683425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
13693425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
13703425bc38SStefano Zampini   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
13713425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
13723425bc38SStefano Zampini   /* return pointers for objects created */
13733425bc38SStefano Zampini   *fetidp_mat=newmat;
13743425bc38SStefano Zampini   *fetidp_pc=newpc;
13753425bc38SStefano Zampini   PetscFunctionReturn(0);
13763425bc38SStefano Zampini }
13771e6b0712SBarry Smith 
13783425bc38SStefano Zampini #undef __FUNCT__
13793425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
13803425bc38SStefano Zampini /*@
138128509bceSStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
13823425bc38SStefano Zampini 
13833425bc38SStefano Zampini    Collective
13843425bc38SStefano Zampini 
13853425bc38SStefano Zampini    Input Parameters:
138628509bceSStefano Zampini +  pc - the BDDC preconditioning context already setup
138728509bceSStefano Zampini 
138828509bceSStefano Zampini    Output Parameters:
138928509bceSStefano Zampini .  fetidp_mat - shell FETIDP matrix object
139028509bceSStefano Zampini .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
139128509bceSStefano Zampini 
139228509bceSStefano Zampini    Options Database Keys:
139328509bceSStefano Zampini -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
13943425bc38SStefano Zampini 
13953425bc38SStefano Zampini    Level: developer
13963425bc38SStefano Zampini 
13973425bc38SStefano Zampini    Notes:
139828509bceSStefano Zampini      Currently the only operation provided for FETIDP matrix is MatMult
13993425bc38SStefano Zampini 
14003425bc38SStefano Zampini .seealso: PCBDDC
14013425bc38SStefano Zampini @*/
14023425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
14033425bc38SStefano Zampini {
14043425bc38SStefano Zampini   PetscErrorCode ierr;
14053425bc38SStefano Zampini 
14063425bc38SStefano Zampini   PetscFunctionBegin;
14073425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14083425bc38SStefano Zampini   if (pc->setupcalled) {
1409516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1410f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
14113425bc38SStefano Zampini   PetscFunctionReturn(0);
14123425bc38SStefano Zampini }
14130c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1414da1bb401SStefano Zampini /*MC
1415da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
14160c7d97c5SJed Brown 
141728509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
141828509bceSStefano Zampini 
141928509bceSStefano Zampini .vb
142028509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
142128509bceSStefano Zampini    [2] A. Klawonn and O. B. Widlund. "Dual-Primal FETI Methods for Linear Elasticity", http://cs.nyu.edu/csweb/Research/TechReports/TR2004-855/TR2004-855.pdf
142228509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
142328509bceSStefano Zampini .ve
142428509bceSStefano Zampini 
142528509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
142628509bceSStefano Zampini 
1427b6fdb6dfSStefano Zampini    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
142828509bceSStefano Zampini 
142928509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
143028509bceSStefano Zampini 
1431b6fdb6dfSStefano Zampini    Unlike 'conventional' interface preconditioners, PCBDDC iterates over all degrees of freedom, not just those on the interface. This allows the use of approximate solvers on the subdomains.
1432b6fdb6dfSStefano Zampini 
143328509bceSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
143428509bceSStefano Zampini 
143528509bceSStefano Zampini    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
143628509bceSStefano Zampini 
143728509bceSStefano Zampini    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace.
143828509bceSStefano Zampini 
1439b6fdb6dfSStefano Zampini    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.
144028509bceSStefano Zampini 
144128509bceSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
144228509bceSStefano Zampini 
1443da1bb401SStefano Zampini    Options Database Keys:
144428509bceSStefano Zampini 
144528509bceSStefano Zampini .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
144628509bceSStefano Zampini .    -pc_bddc_use_edges <1> - use or not edges in primal space
1447b6fdb6dfSStefano Zampini .    -pc_bddc_use_faces <0> - use or not faces in primal space
144828509bceSStefano Zampini .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
144928509bceSStefano Zampini .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
145028509bceSStefano Zampini .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
145128509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
145228509bceSStefano Zampini .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
145328509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
145428509bceSStefano Zampini 
145528509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
145628509bceSStefano Zampini .vb
145728509bceSStefano Zampini       -pc_bddc_dirichlet_
145828509bceSStefano Zampini       -pc_bddc_neumann_
145928509bceSStefano Zampini       -pc_bddc_coarse_
146028509bceSStefano Zampini .ve
146128509bceSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
146228509bceSStefano Zampini 
146328509bceSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level can be specified as
146428509bceSStefano Zampini .vb
146528509bceSStefano Zampini       -pc_bddc_dirichlet_N_
146628509bceSStefano Zampini       -pc_bddc_neumann_N_
146728509bceSStefano Zampini       -pc_bddc_coarse_N_
146828509bceSStefano Zampini .ve
146928509bceSStefano Zampini    Note that level number ranges from the finest 0 to the coarsest N
1470da1bb401SStefano Zampini 
1471da1bb401SStefano Zampini    Level: intermediate
1472da1bb401SStefano Zampini 
1473b6fdb6dfSStefano Zampini    Developer notes:
147428509bceSStefano Zampini      Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose
1475da1bb401SStefano Zampini 
147628509bceSStefano Zampini      New deluxe scaling operator will be available soon.
1477da1bb401SStefano Zampini 
1478da1bb401SStefano Zampini    Contributed by Stefano Zampini
1479da1bb401SStefano Zampini 
1480da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1481da1bb401SStefano Zampini M*/
1482b2573a8aSBarry Smith 
1483da1bb401SStefano Zampini #undef __FUNCT__
1484da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
14858cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1486da1bb401SStefano Zampini {
1487da1bb401SStefano Zampini   PetscErrorCode      ierr;
1488da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1489da1bb401SStefano Zampini 
1490da1bb401SStefano Zampini   PetscFunctionBegin;
1491da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1492da1bb401SStefano Zampini   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1493da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1494da1bb401SStefano Zampini 
1495da1bb401SStefano Zampini   /* create PCIS data structure */
1496da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1497da1bb401SStefano Zampini 
149847d04d0dSStefano Zampini   /* BDDC customization */
149947d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
150047d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
150147d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
150247d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
150347d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
150447d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
150547d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
150647d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
150747d04d0dSStefano Zampini 
1508674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
15090bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
15103972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1511534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1512534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1513534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1514da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1515da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1516da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1517da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1518da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
151915aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
152015aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1521da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1522da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1523da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1524da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1525da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1526da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1527da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1528da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1529da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1530da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
1531*785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
1532*785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
1533*785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
153460d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
153560d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
1536da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
1537da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
153885c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
153947d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
154047d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
15414fad6a16SStefano Zampini   pcbddc->current_level              = 0;
15422b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1543da1bb401SStefano Zampini 
1544674ae819SStefano Zampini   /* create local graph structure */
1545674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1546674ae819SStefano Zampini 
1547674ae819SStefano Zampini   /* scaling */
1548674ae819SStefano Zampini   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1549674ae819SStefano Zampini   pcbddc->work_scaling               = 0;
1550da1bb401SStefano Zampini 
1551da1bb401SStefano Zampini   /* function pointers */
1552da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
1553da1bb401SStefano Zampini   pc->ops->applytranspose      = 0;
1554da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1555da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1556da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1557da1bb401SStefano Zampini   pc->ops->view                = 0;
1558da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1559da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1560da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1561534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1562534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1563da1bb401SStefano Zampini 
1564da1bb401SStefano Zampini   /* composing function */
1565674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1566bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
15672b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1568b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
15692b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1570bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1571*785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
157282ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1573*785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
157482ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1575*785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
157682ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1577*785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
157882ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1579bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1580bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1581bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1582bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1583bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1584da1bb401SStefano Zampini   PetscFunctionReturn(0);
1585da1bb401SStefano Zampini }
15863425bc38SStefano Zampini 
1587