xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 727cdba63cef3c72c66ef27a9d74c90cecc5d6d6)
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    - MAT_IGNORE_ZERO_ENTRIES for Constraints Matrix
7eb97c9d2SStefano Zampini 
8eb97c9d2SStefano Zampini    Solvers
9eb97c9d2SStefano Zampini    - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers)
10eb97c9d2SStefano Zampini    - reuse already allocated coarse matrix if possible
11eb97c9d2SStefano Zampini    - Propagate ksp prefixes for solvers to mat objects?
12eb97c9d2SStefano Zampini    - Propagate nearnullspace info among levels
13eb97c9d2SStefano Zampini 
14eb97c9d2SStefano Zampini    User interface
15eb97c9d2SStefano Zampini    - Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet)
16eb97c9d2SStefano Zampini    - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access)
17eb97c9d2SStefano Zampini    - Provide PCApplyTranpose_BDDC
18eb97c9d2SStefano Zampini    - DofSplitting and DM attached to pc?
19eb97c9d2SStefano Zampini 
20eb97c9d2SStefano Zampini    Debugging output
21eb97c9d2SStefano Zampini    - Better management of verbosity levels of debugging output
22eb97c9d2SStefano Zampini    - Crashes on some architecture -> call SynchronizedAllow before every SynchronizedPrintf
23eb97c9d2SStefano Zampini 
24eb97c9d2SStefano Zampini    Build
25eb97c9d2SStefano Zampini    - make runexe59
26eb97c9d2SStefano Zampini 
27eb97c9d2SStefano Zampini    Extra
28eb97c9d2SStefano Zampini    - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
29eb97c9d2SStefano Zampini    - Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver?
30eb97c9d2SStefano Zampini    - add support for computing h,H and related using coordinates?
31c0b83709SStefano Zampini    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
32eb97c9d2SStefano Zampini    - Better management in PCIS code
33eb97c9d2SStefano Zampini    - BDDC with MG framework?
34eb97c9d2SStefano Zampini 
35eb97c9d2SStefano Zampini    FETIDP
36eb97c9d2SStefano Zampini    - Move FETIDP code to its own classes
37eb97c9d2SStefano Zampini 
38eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
39fb180af4SStefano Zampini    - Add MAT_REUSE in MatConvert_IS_AIJ
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__
3033b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
3043b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
3053b03a366Sstefano_zampini {
3063b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3073b03a366Sstefano_zampini   PetscErrorCode ierr;
3083b03a366Sstefano_zampini 
3093b03a366Sstefano_zampini   PetscFunctionBegin;
3103b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
31136e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
31236e030ebSStefano Zampini   pcbddc->DirichletBoundaries=DirichletBoundaries;
313fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3143b03a366Sstefano_zampini   PetscFunctionReturn(0);
3153b03a366Sstefano_zampini }
3161e6b0712SBarry Smith 
3173b03a366Sstefano_zampini #undef __FUNCT__
3183b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
3193b03a366Sstefano_zampini /*@
32028509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
3213b03a366Sstefano_zampini 
3223b03a366Sstefano_zampini    Not collective
3233b03a366Sstefano_zampini 
3243b03a366Sstefano_zampini    Input Parameters:
3253b03a366Sstefano_zampini +  pc - the preconditioning context
32628509bceSStefano Zampini -  DirichletBoundaries - sequential IS defining the subdomain part of Dirichlet boundaries (in local ordering)
3273b03a366Sstefano_zampini 
3283b03a366Sstefano_zampini    Level: intermediate
3293b03a366Sstefano_zampini 
3303b03a366Sstefano_zampini    Notes:
3313b03a366Sstefano_zampini 
3323b03a366Sstefano_zampini .seealso: PCBDDC
3333b03a366Sstefano_zampini @*/
3343b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3353b03a366Sstefano_zampini {
3363b03a366Sstefano_zampini   PetscErrorCode ierr;
3373b03a366Sstefano_zampini 
3383b03a366Sstefano_zampini   PetscFunctionBegin;
3393b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
340674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
3413b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
3423b03a366Sstefano_zampini   PetscFunctionReturn(0);
3433b03a366Sstefano_zampini }
3443b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
3451e6b0712SBarry Smith 
3463b03a366Sstefano_zampini #undef __FUNCT__
3470c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
34853cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
3490c7d97c5SJed Brown {
3500c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
35153cdbc3dSStefano Zampini   PetscErrorCode ierr;
3520c7d97c5SJed Brown 
3530c7d97c5SJed Brown   PetscFunctionBegin;
35453cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
35536e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
35636e030ebSStefano Zampini   pcbddc->NeumannBoundaries=NeumannBoundaries;
357fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3580c7d97c5SJed Brown   PetscFunctionReturn(0);
3590c7d97c5SJed Brown }
3601e6b0712SBarry Smith 
3610c7d97c5SJed Brown #undef __FUNCT__
3620c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
36357527edcSJed Brown /*@
36428509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
36557527edcSJed Brown 
3669c0446d6SStefano Zampini    Not collective
36757527edcSJed Brown 
36857527edcSJed Brown    Input Parameters:
36957527edcSJed Brown +  pc - the preconditioning context
37028509bceSStefano Zampini -  NeumannBoundaries - sequential IS defining the subdomain part of Neumann boundaries (in local ordering)
37157527edcSJed Brown 
37257527edcSJed Brown    Level: intermediate
37357527edcSJed Brown 
37457527edcSJed Brown    Notes:
37557527edcSJed Brown 
37657527edcSJed Brown .seealso: PCBDDC
37757527edcSJed Brown @*/
37853cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
3790c7d97c5SJed Brown {
3800c7d97c5SJed Brown   PetscErrorCode ierr;
3810c7d97c5SJed Brown 
3820c7d97c5SJed Brown   PetscFunctionBegin;
3830c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
384674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
38553cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
38653cdbc3dSStefano Zampini   PetscFunctionReturn(0);
38753cdbc3dSStefano Zampini }
38853cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
3891e6b0712SBarry Smith 
39053cdbc3dSStefano Zampini #undef __FUNCT__
391da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
392da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
393da1bb401SStefano Zampini {
394da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
395da1bb401SStefano Zampini 
396da1bb401SStefano Zampini   PetscFunctionBegin;
397da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
398da1bb401SStefano Zampini   PetscFunctionReturn(0);
399da1bb401SStefano Zampini }
4001e6b0712SBarry Smith 
401da1bb401SStefano Zampini #undef __FUNCT__
402da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
403da1bb401SStefano Zampini /*@
40428509bceSStefano Zampini  PCBDDCGetDirichletBoundaries - Get IS for local Dirichlet boundaries
405da1bb401SStefano Zampini 
406da1bb401SStefano Zampini    Not collective
407da1bb401SStefano Zampini 
408da1bb401SStefano Zampini    Input Parameters:
40928509bceSStefano Zampini .  pc - the preconditioning context
410da1bb401SStefano Zampini 
411da1bb401SStefano Zampini    Output Parameters:
41228509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
413da1bb401SStefano Zampini 
414da1bb401SStefano Zampini    Level: intermediate
415da1bb401SStefano Zampini 
416da1bb401SStefano Zampini    Notes:
417da1bb401SStefano Zampini 
418da1bb401SStefano Zampini .seealso: PCBDDC
419da1bb401SStefano Zampini @*/
420da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
421da1bb401SStefano Zampini {
422da1bb401SStefano Zampini   PetscErrorCode ierr;
423da1bb401SStefano Zampini 
424da1bb401SStefano Zampini   PetscFunctionBegin;
425da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
426da1bb401SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
427da1bb401SStefano Zampini   PetscFunctionReturn(0);
428da1bb401SStefano Zampini }
429da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
4301e6b0712SBarry Smith 
431da1bb401SStefano Zampini #undef __FUNCT__
43253cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
43353cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
43453cdbc3dSStefano Zampini {
43553cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
43653cdbc3dSStefano Zampini 
43753cdbc3dSStefano Zampini   PetscFunctionBegin;
43853cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
43953cdbc3dSStefano Zampini   PetscFunctionReturn(0);
44053cdbc3dSStefano Zampini }
4411e6b0712SBarry Smith 
44253cdbc3dSStefano Zampini #undef __FUNCT__
44353cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
44453cdbc3dSStefano Zampini /*@
44528509bceSStefano Zampini  PCBDDCGetNeumannBoundaries - Get IS for local Neumann boundaries
44653cdbc3dSStefano Zampini 
4479c0446d6SStefano Zampini    Not collective
44853cdbc3dSStefano Zampini 
44953cdbc3dSStefano Zampini    Input Parameters:
45028509bceSStefano Zampini .  pc - the preconditioning context
45153cdbc3dSStefano Zampini 
45253cdbc3dSStefano Zampini    Output Parameters:
45328509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
45453cdbc3dSStefano Zampini 
45553cdbc3dSStefano Zampini    Level: intermediate
45653cdbc3dSStefano Zampini 
45753cdbc3dSStefano Zampini    Notes:
45853cdbc3dSStefano Zampini 
45953cdbc3dSStefano Zampini .seealso: PCBDDC
46053cdbc3dSStefano Zampini @*/
46153cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
46253cdbc3dSStefano Zampini {
46353cdbc3dSStefano Zampini   PetscErrorCode ierr;
46453cdbc3dSStefano Zampini 
46553cdbc3dSStefano Zampini   PetscFunctionBegin;
46653cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
46753cdbc3dSStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
4680c7d97c5SJed Brown   PetscFunctionReturn(0);
4690c7d97c5SJed Brown }
47036e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
4711e6b0712SBarry Smith 
47236e030ebSStefano Zampini #undef __FUNCT__
473da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
4741a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
47536e030ebSStefano Zampini {
47636e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
477da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
478da1bb401SStefano Zampini   PetscErrorCode ierr;
47936e030ebSStefano Zampini 
48036e030ebSStefano Zampini   PetscFunctionBegin;
481674ae819SStefano Zampini   /* free old CSR */
482674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
483fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
484674ae819SStefano Zampini   /* get CSR into graph structure */
485da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
486674ae819SStefano Zampini     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
487674ae819SStefano Zampini     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
488674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
489674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
490da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
4911a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
4921a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
493674ae819SStefano Zampini   }
494575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
49536e030ebSStefano Zampini   PetscFunctionReturn(0);
49636e030ebSStefano Zampini }
4971e6b0712SBarry Smith 
49836e030ebSStefano Zampini #undef __FUNCT__
499da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
50036e030ebSStefano Zampini /*@
50128509bceSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
50236e030ebSStefano Zampini 
50336e030ebSStefano Zampini    Not collective
50436e030ebSStefano Zampini 
50536e030ebSStefano Zampini    Input Parameters:
50636e030ebSStefano Zampini +  pc - the preconditioning context
50728509bceSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
50828509bceSStefano Zampini .  xadj, adjncy - the CSR graph
50928509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
51036e030ebSStefano Zampini 
51136e030ebSStefano Zampini    Level: intermediate
51236e030ebSStefano Zampini 
51336e030ebSStefano Zampini    Notes:
51436e030ebSStefano Zampini 
51528509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
51636e030ebSStefano Zampini @*/
5171a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
51836e030ebSStefano Zampini {
519575ad6abSStefano Zampini   void (*f)(void) = 0;
52036e030ebSStefano Zampini   PetscErrorCode ierr;
52136e030ebSStefano Zampini 
52236e030ebSStefano Zampini   PetscFunctionBegin;
52336e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
524674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
525674ae819SStefano Zampini   PetscValidIntPointer(xadj,4);
526674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
527674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
528da1bb401SStefano Zampini   }
52936e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
530575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
531575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
532575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
533575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
534575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
53536e030ebSStefano Zampini   }
53636e030ebSStefano Zampini   PetscFunctionReturn(0);
53736e030ebSStefano Zampini }
5389c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
5391e6b0712SBarry Smith 
5409c0446d6SStefano Zampini #undef __FUNCT__
5419c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
5429c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
5439c0446d6SStefano Zampini {
5449c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5459c0446d6SStefano Zampini   PetscInt i;
5469c0446d6SStefano Zampini   PetscErrorCode ierr;
5479c0446d6SStefano Zampini 
5489c0446d6SStefano Zampini   PetscFunctionBegin;
549da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
5509c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
5519c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
5529c0446d6SStefano Zampini   }
553d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
554da1bb401SStefano Zampini   /* allocate space then set */
5559c0446d6SStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
5569c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
557da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
558da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
5599c0446d6SStefano Zampini   }
5609c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
56160d44989SStefano Zampini   pcbddc->user_provided_isfordofs = PETSC_TRUE;
5629c0446d6SStefano Zampini   PetscFunctionReturn(0);
5639c0446d6SStefano Zampini }
5641e6b0712SBarry Smith 
5659c0446d6SStefano Zampini #undef __FUNCT__
5669c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
5679c0446d6SStefano Zampini /*@
56828509bceSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the local Neumann matrix
5699c0446d6SStefano Zampini 
5709c0446d6SStefano Zampini    Not collective
5719c0446d6SStefano Zampini 
5729c0446d6SStefano Zampini    Input Parameters:
5739c0446d6SStefano Zampini +  pc - the preconditioning context
57428509bceSStefano Zampini -  n_is - number of index sets defining the fields
57528509bceSStefano Zampini .  ISForDofs - array of IS describing the fields
5769c0446d6SStefano Zampini 
5779c0446d6SStefano Zampini    Level: intermediate
5789c0446d6SStefano Zampini 
5799c0446d6SStefano Zampini    Notes:
5809c0446d6SStefano Zampini 
5819c0446d6SStefano Zampini .seealso: PCBDDC
5829c0446d6SStefano Zampini @*/
5839c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
5849c0446d6SStefano Zampini {
5852b510759SStefano Zampini   PetscInt       i;
5869c0446d6SStefano Zampini   PetscErrorCode ierr;
5879c0446d6SStefano Zampini 
5889c0446d6SStefano Zampini   PetscFunctionBegin;
5899c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5902b510759SStefano Zampini   for (i=0;i<n_is;i++) {
5912b510759SStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2);
5922b510759SStefano Zampini   }
5939c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
5949c0446d6SStefano Zampini   PetscFunctionReturn(0);
5959c0446d6SStefano Zampini }
596da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
597534831adSStefano Zampini #undef __FUNCT__
598534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
599534831adSStefano Zampini /* -------------------------------------------------------------------------- */
600534831adSStefano Zampini /*
601534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
602534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
6039c0446d6SStefano Zampini 
604534831adSStefano Zampini    Input Parameter:
605534831adSStefano Zampini +  pc - the preconditioner contex
606534831adSStefano Zampini 
607534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
608534831adSStefano Zampini 
609534831adSStefano Zampini    Notes:
610534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
611534831adSStefano Zampini    the user, but instead is called by KSPSolve().
612534831adSStefano Zampini */
613534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
614534831adSStefano Zampini {
615534831adSStefano Zampini   PetscErrorCode ierr;
616534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
617534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
618534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
619534831adSStefano Zampini   Mat            temp_mat;
6203972b0daSStefano Zampini   IS             dirIS;
6213972b0daSStefano Zampini   PetscInt       dirsize,i,*is_indices;
6223972b0daSStefano Zampini   PetscScalar    *array_x,*array_diagonal;
6233972b0daSStefano Zampini   Vec            used_vec;
62492e3dcfbSStefano Zampini   PetscBool      guess_nonzero,flg,bddc_has_dirichlet_boundaries;
625534831adSStefano Zampini 
626534831adSStefano Zampini   PetscFunctionBegin;
62785c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
62885c4d303SStefano Zampini   if (ksp) {
62985c4d303SStefano Zampini     PetscBool iscg;
63085c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
63185c4d303SStefano Zampini     if (!iscg) {
63285c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
63385c4d303SStefano Zampini     }
63485c4d303SStefano Zampini   }
63585c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
63662a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
63762a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
63862a6ff1dSStefano Zampini   }
63962a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
64062a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
64162a6ff1dSStefano Zampini   }
6423972b0daSStefano Zampini   if (x) {
6433972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
6443972b0daSStefano Zampini     used_vec = x;
6453972b0daSStefano Zampini   } else {
6463972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
6473972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
6483972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6493972b0daSStefano Zampini   }
6503972b0daSStefano Zampini   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
6513972b0daSStefano Zampini   if (ksp) {
6523972b0daSStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
6533972b0daSStefano Zampini     if (!guess_nonzero) {
6543972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6553972b0daSStefano Zampini     }
6563972b0daSStefano Zampini   }
6573308cffdSStefano Zampini 
65892e3dcfbSStefano Zampini   /* TODO: remove when Dirichlet boundaries will be shared */
65992e3dcfbSStefano Zampini   ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr);
66092e3dcfbSStefano Zampini   flg = PETSC_FALSE;
66192e3dcfbSStefano Zampini   if (dirIS) flg = PETSC_TRUE;
66292e3dcfbSStefano Zampini   ierr = MPI_Allreduce(&flg,&bddc_has_dirichlet_boundaries,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
66392e3dcfbSStefano Zampini 
6643972b0daSStefano Zampini   /* store the original rhs */
6653972b0daSStefano Zampini   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
6663972b0daSStefano Zampini 
6673972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
66892e3dcfbSStefano Zampini   if (rhs && bddc_has_dirichlet_boundaries) {
6693972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
6703972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
6713972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6723972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6733972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6743972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6753972b0daSStefano Zampini     if (dirIS) {
6763972b0daSStefano Zampini       ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr);
6773972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6783972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6793972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6802fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
6813972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6823972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6833972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6843972b0daSStefano Zampini     }
6853972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6863972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
687b76ba322SStefano Zampini 
6883972b0daSStefano Zampini     /* remove the computed solution from the rhs */
6893972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6903972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
6913972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6923308cffdSStefano Zampini   }
693b76ba322SStefano Zampini 
694b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
6953972b0daSStefano Zampini   if (x) {
6963972b0daSStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
6973972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
69885c4d303SStefano Zampini     if (pcbddc->use_exact_dirichlet_trick) {
699b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
700b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
701b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
702b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
703b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
704b76ba322SStefano Zampini       if (ksp) {
705b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
706b76ba322SStefano Zampini       }
707b76ba322SStefano Zampini     }
7083972b0daSStefano Zampini   }
709b76ba322SStefano Zampini 
71092e3dcfbSStefano Zampini   /* prepare MatMult and rhs for solver */
711674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
712b76ba322SStefano Zampini     /* swap pointers for local matrices */
713b76ba322SStefano Zampini     temp_mat = matis->A;
714b76ba322SStefano Zampini     matis->A = pcbddc->local_mat;
715b76ba322SStefano Zampini     pcbddc->local_mat = temp_mat;
71692e3dcfbSStefano Zampini     if (rhs) {
717b76ba322SStefano Zampini       /* Get local rhs and apply transformation of basis */
718b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
719b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
720b76ba322SStefano Zampini       /* from original basis to modified basis */
721b76ba322SStefano Zampini       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
722b76ba322SStefano Zampini       /* put back modified values into the global vec using INSERT_VALUES copy mode */
723b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
724b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
725674ae819SStefano Zampini     }
72692e3dcfbSStefano Zampini   }
72792e3dcfbSStefano Zampini 
72892e3dcfbSStefano Zampini   /* remove nullspace if present */
7290bdf917eSStefano Zampini   if (ksp && pcbddc->NullSpace) {
730d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
731d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
732b76ba322SStefano Zampini   }
7330bdf917eSStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
734534831adSStefano Zampini   PetscFunctionReturn(0);
735534831adSStefano Zampini }
736534831adSStefano Zampini /* -------------------------------------------------------------------------- */
737534831adSStefano Zampini #undef __FUNCT__
738534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
739534831adSStefano Zampini /* -------------------------------------------------------------------------- */
740534831adSStefano Zampini /*
741534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
742534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
743534831adSStefano Zampini 
744534831adSStefano Zampini    Input Parameter:
745534831adSStefano Zampini +  pc - the preconditioner contex
746534831adSStefano Zampini 
747534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
748534831adSStefano Zampini 
749534831adSStefano Zampini    Notes:
750534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
751534831adSStefano Zampini    the user, but instead is called by KSPSolve().
752534831adSStefano Zampini */
753534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
754534831adSStefano Zampini {
755534831adSStefano Zampini   PetscErrorCode ierr;
756534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
757534831adSStefano Zampini   PC_IS          *pcis   = (PC_IS*)(pc->data);
758534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
759534831adSStefano Zampini   Mat            temp_mat;
760534831adSStefano Zampini 
761534831adSStefano Zampini   PetscFunctionBegin;
762674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
763534831adSStefano Zampini     /* swap pointers for local matrices */
764534831adSStefano Zampini     temp_mat = matis->A;
765534831adSStefano Zampini     matis->A = pcbddc->local_mat;
766534831adSStefano Zampini     pcbddc->local_mat = temp_mat;
7673425bc38SStefano Zampini   }
7683308cffdSStefano Zampini   if (pcbddc->use_change_of_basis && x) {
769534831adSStefano Zampini     /* Get Local boundary and apply transformation of basis to solution vector */
770534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
771534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
772534831adSStefano Zampini     /* from modified basis to original basis */
773534831adSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
774534831adSStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
775534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
776534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
777534831adSStefano Zampini   }
7783972b0daSStefano Zampini   /* add solution removed in presolve */
7793425bc38SStefano Zampini   if (x) {
7803425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
7813425bc38SStefano Zampini   }
782fb223d50SStefano Zampini   /* restore rhs to its original state */
783fb223d50SStefano Zampini   if (rhs) {
784fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
785fb223d50SStefano Zampini   }
786534831adSStefano Zampini   PetscFunctionReturn(0);
787534831adSStefano Zampini }
788534831adSStefano Zampini /* -------------------------------------------------------------------------- */
78953cdbc3dSStefano Zampini #undef __FUNCT__
79053cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
7910c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
7920c7d97c5SJed Brown /*
7930c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
7940c7d97c5SJed Brown                   by setting data structures and options.
7950c7d97c5SJed Brown 
7960c7d97c5SJed Brown    Input Parameter:
79753cdbc3dSStefano Zampini +  pc - the preconditioner context
7980c7d97c5SJed Brown 
7990c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
8000c7d97c5SJed Brown 
8010c7d97c5SJed Brown    Notes:
8020c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
8030c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
8040c7d97c5SJed Brown */
80553cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
8060c7d97c5SJed Brown {
8070c7d97c5SJed Brown   PetscErrorCode   ierr;
8080c7d97c5SJed Brown   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
809f4ddd8eeSStefano Zampini   MatNullSpace     nearnullspace;
810f4ddd8eeSStefano Zampini   const Vec        *nearnullvecs,*onearnullvecs;
811674ae819SStefano Zampini   MatStructure     flag;
812f4ddd8eeSStefano Zampini   PetscObjectState state;
813f4ddd8eeSStefano Zampini   PetscInt         nnsp_size,onnsp_size;
814674ae819SStefano Zampini   PetscBool        computeis,computetopography,computesolvers;
815f4ddd8eeSStefano Zampini   PetscBool        new_nearnullspace_provided,nnsp_has_cnst,onnsp_has_cnst;
8160c7d97c5SJed Brown 
8170c7d97c5SJed Brown   PetscFunctionBegin;
818f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
819f4ddd8eeSStefano Zampini   /* PCIS does not support MatStructure flags different from SAME_PRECONDITIONER */
8203b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
821fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
822f4ddd8eeSStefano Zampini 
823f4ddd8eeSStefano Zampini   /* split work */
824674ae819SStefano Zampini   if (pc->setupcalled) {
825674ae819SStefano Zampini     computeis = PETSC_FALSE;
826674ae819SStefano Zampini     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
827674ae819SStefano Zampini     if (flag == SAME_PRECONDITIONER) {
828f4ddd8eeSStefano Zampini       PetscFunctionReturn(0);
829674ae819SStefano Zampini     } else if (flag == SAME_NONZERO_PATTERN) {
830674ae819SStefano Zampini       computetopography = PETSC_FALSE;
831674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
832674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
833674ae819SStefano Zampini       computetopography = PETSC_TRUE;
834674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
835674ae819SStefano Zampini     }
836674ae819SStefano Zampini   } else {
837674ae819SStefano Zampini     computeis = PETSC_TRUE;
838674ae819SStefano Zampini     computetopography = PETSC_TRUE;
839674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
840674ae819SStefano Zampini   }
841fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
842fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
843fb180af4SStefano Zampini   }
844f4ddd8eeSStefano Zampini 
845f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
846f4ddd8eeSStefano Zampini   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
847f4ddd8eeSStefano Zampini     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
848f4ddd8eeSStefano Zampini     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
849f4ddd8eeSStefano Zampini     if (pcbddc->current_level) {
850f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
851f4ddd8eeSStefano Zampini     }
852f4ddd8eeSStefano Zampini   }
853f4ddd8eeSStefano Zampini 
854fcd409f5SStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
855674ae819SStefano Zampini   if (computeis) {
856fcd409f5SStefano Zampini     /* HACK INTO PCIS */
857fcd409f5SStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
858fcd409f5SStefano Zampini     pcis->computesolvers = PETSC_FALSE;
859674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
860674ae819SStefano Zampini   }
861f4ddd8eeSStefano Zampini 
862f4ddd8eeSStefano Zampini   /* Analyze interface */
863674ae819SStefano Zampini   if (computetopography) {
864674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
865fb8d54d4SStefano Zampini   }
866fb8d54d4SStefano Zampini 
867f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
868fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
869f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
870f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
871f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
872f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
873f4ddd8eeSStefano Zampini     } else {
874f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
875f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
876f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
877f4ddd8eeSStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace */
878f4ddd8eeSStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
879f4ddd8eeSStefano Zampini         ierr = MatNullSpaceGetVecs(pcbddc->onearnullspace,&onnsp_has_cnst,&onnsp_size,&onearnullvecs);CHKERRQ(ierr);
880f4ddd8eeSStefano Zampini         if ( (nnsp_has_cnst != onnsp_has_cnst) || (nnsp_size != onnsp_size) ) {
881f4ddd8eeSStefano Zampini           new_nearnullspace_provided = PETSC_TRUE;
882f4ddd8eeSStefano Zampini         } else { /* nullspaces have the same size, so check vectors or their ObjectStateId */
883f4ddd8eeSStefano Zampini           PetscInt i;
884f4ddd8eeSStefano Zampini           for (i=0;i<nnsp_size;i++) {
885f4ddd8eeSStefano Zampini             ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
886f4ddd8eeSStefano Zampini             if (nearnullvecs[i] != onearnullvecs[i] || pcbddc->onearnullvecs_state[i] != state) {
887f4ddd8eeSStefano Zampini               new_nearnullspace_provided = PETSC_TRUE;
888f4ddd8eeSStefano Zampini               break;
889f4ddd8eeSStefano Zampini             }
890f4ddd8eeSStefano Zampini           }
891f4ddd8eeSStefano Zampini         }
892f4ddd8eeSStefano Zampini       }
893f4ddd8eeSStefano Zampini     }
894f4ddd8eeSStefano Zampini   } else {
895f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
896f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
897f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
898f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
899f4ddd8eeSStefano Zampini     }
900f4ddd8eeSStefano Zampini   }
901f4ddd8eeSStefano Zampini 
902f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
903*727cdba6SStefano Zampini   /* reset primal space flags */
904f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
905*727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
906fb8d54d4SStefano Zampini   if (computetopography || new_nearnullspace_provided) {
907*727cdba6SStefano Zampini     /* It also sets the primal space flags */
908674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
909e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
910f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
911674ae819SStefano Zampini   }
912fb8d54d4SStefano Zampini 
913f4ddd8eeSStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
914674ae819SStefano Zampini     /* reset data */
915674ae819SStefano Zampini     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
91699cc7994SStefano Zampini     /* Create coarse and local stuffs */
91799cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
918674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
9190c7d97c5SJed Brown   }
9202b510759SStefano Zampini   if (pcbddc->dbg_flag && pcbddc->current_level) {
9212b510759SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
9222b510759SStefano Zampini   }
9230c7d97c5SJed Brown   PetscFunctionReturn(0);
9240c7d97c5SJed Brown }
9250c7d97c5SJed Brown 
9260c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
9270c7d97c5SJed Brown /*
9280c7d97c5SJed Brown    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
9290c7d97c5SJed Brown 
9300c7d97c5SJed Brown    Input Parameters:
9310c7d97c5SJed Brown .  pc - the preconditioner context
9320c7d97c5SJed Brown .  r - input vector (global)
9330c7d97c5SJed Brown 
9340c7d97c5SJed Brown    Output Parameter:
9350c7d97c5SJed Brown .  z - output vector (global)
9360c7d97c5SJed Brown 
9370c7d97c5SJed Brown    Application Interface Routine: PCApply()
9380c7d97c5SJed Brown  */
9390c7d97c5SJed Brown #undef __FUNCT__
9400c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
94153cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
9420c7d97c5SJed Brown {
9430c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
9440c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
9450c7d97c5SJed Brown   PetscErrorCode    ierr;
9463b03a366Sstefano_zampini   const PetscScalar one = 1.0;
9473b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
9482617d88aSStefano Zampini   const PetscScalar zero = 0.0;
9490c7d97c5SJed Brown 
9500c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
9510c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
9528eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
9530c7d97c5SJed Brown 
9540c7d97c5SJed Brown   PetscFunctionBegin;
95585c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
9560c7d97c5SJed Brown     /* First Dirichlet solve */
9570c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9580c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
95953cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
9600c7d97c5SJed Brown     /*
9610c7d97c5SJed Brown       Assembling right hand side for BDDC operator
962674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
963674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
9640c7d97c5SJed Brown     */
9650c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
9660c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
9678eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
9680c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
9690c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
9700c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9710c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
972674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
973b76ba322SStefano Zampini   } else {
9740bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
975b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
976674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
977b76ba322SStefano Zampini   }
978b76ba322SStefano Zampini 
9792617d88aSStefano Zampini   /* Apply interface preconditioner
9802617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
9812617d88aSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
9822617d88aSStefano Zampini 
983674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
984674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
9850c7d97c5SJed Brown 
9863b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
9870c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9880c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9890c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
9908eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
991df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
992df187020SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
993df187020SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY (pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
994df187020SStefano Zampini   ierr = VecAXPY (pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
9950c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9960c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9970c7d97c5SJed Brown   PetscFunctionReturn(0);
9980c7d97c5SJed Brown }
999da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1000674ae819SStefano Zampini 
1001da1bb401SStefano Zampini #undef __FUNCT__
1002da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1003da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1004da1bb401SStefano Zampini {
1005da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1006da1bb401SStefano Zampini   PetscErrorCode ierr;
1007da1bb401SStefano Zampini 
1008da1bb401SStefano Zampini   PetscFunctionBegin;
1009da1bb401SStefano Zampini   /* free data created by PCIS */
1010da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1011674ae819SStefano Zampini   /* free BDDC custom data  */
1012674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1013674ae819SStefano Zampini   /* destroy objects related to topography */
1014674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1015674ae819SStefano Zampini   /* free allocated graph structure */
1016da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1017674ae819SStefano Zampini   /* free data for scaling operator */
1018674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1019674ae819SStefano Zampini   /* free solvers stuff */
1020674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
102162a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
102262a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
102362a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
10243425bc38SStefano Zampini   /* remove functions */
1025674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1026bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
10272b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1028b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
10292b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1030bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1031bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1032bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1033bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1034bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1035bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1036bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1037bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1038bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1039bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1040674ae819SStefano Zampini   /* Free the private data structure */
1041674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1042da1bb401SStefano Zampini   PetscFunctionReturn(0);
1043da1bb401SStefano Zampini }
10443425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
10451e6b0712SBarry Smith 
10463425bc38SStefano Zampini #undef __FUNCT__
10473425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
10483425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
10493425bc38SStefano Zampini {
1050674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10513425bc38SStefano Zampini   PC_IS*         pcis;
10523425bc38SStefano Zampini   PC_BDDC*       pcbddc;
10533425bc38SStefano Zampini   PetscErrorCode ierr;
10540c7d97c5SJed Brown 
10553425bc38SStefano Zampini   PetscFunctionBegin;
10563425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10573425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
10583425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
10593425bc38SStefano Zampini 
10603425bc38SStefano Zampini   /* change of basis for physical rhs if needed
10613425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
10623308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
10633425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
10643425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10653425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1066fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1067fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
10683425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10693425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1070674ae819SStefano Zampini   /* Apply partition of unity */
10713425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1072674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
10738eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
10743425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
10753425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10763425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
10773425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
10783425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
10793425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10803425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1081674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
10823425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10833425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10843425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
10853425bc38SStefano Zampini   }
10863425bc38SStefano Zampini   /* BDDC rhs */
10873425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
10888eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
10893425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10903425bc38SStefano Zampini   }
10913425bc38SStefano Zampini   /* apply BDDC */
10923425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
10933425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
10943425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
10953425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
10963425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10973425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10983425bc38SStefano Zampini   /* restore original rhs */
10993425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
11003425bc38SStefano Zampini   PetscFunctionReturn(0);
11013425bc38SStefano Zampini }
11021e6b0712SBarry Smith 
11033425bc38SStefano Zampini #undef __FUNCT__
11043425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
11053425bc38SStefano Zampini /*@
110628509bceSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
11073425bc38SStefano Zampini 
11083425bc38SStefano Zampini    Collective
11093425bc38SStefano Zampini 
11103425bc38SStefano Zampini    Input Parameters:
111128509bceSStefano Zampini +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
111228509bceSStefano Zampini .  standard_rhs - the right-hand side for your linear system
11133425bc38SStefano Zampini 
11143425bc38SStefano Zampini    Output Parameters:
111528509bceSStefano Zampini -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
11163425bc38SStefano Zampini 
11173425bc38SStefano Zampini    Level: developer
11183425bc38SStefano Zampini 
11193425bc38SStefano Zampini    Notes:
11203425bc38SStefano Zampini 
112128509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
11223425bc38SStefano Zampini @*/
11233425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
11243425bc38SStefano Zampini {
1125674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
11263425bc38SStefano Zampini   PetscErrorCode ierr;
11273425bc38SStefano Zampini 
11283425bc38SStefano Zampini   PetscFunctionBegin;
11293425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
11303425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
11313425bc38SStefano Zampini   PetscFunctionReturn(0);
11323425bc38SStefano Zampini }
11333425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
11341e6b0712SBarry Smith 
11353425bc38SStefano Zampini #undef __FUNCT__
11363425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
11373425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
11383425bc38SStefano Zampini {
1139674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
11403425bc38SStefano Zampini   PC_IS*         pcis;
11413425bc38SStefano Zampini   PC_BDDC*       pcbddc;
11423425bc38SStefano Zampini   PetscErrorCode ierr;
11433425bc38SStefano Zampini 
11443425bc38SStefano Zampini   PetscFunctionBegin;
11453425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
11463425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
11473425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
11483425bc38SStefano Zampini 
11493425bc38SStefano Zampini   /* apply B_delta^T */
11503425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11513425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11523425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
11533425bc38SStefano Zampini   /* compute rhs for BDDC application */
11543425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
11558eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
11563425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
11573425bc38SStefano Zampini   }
11583425bc38SStefano Zampini   /* apply BDDC */
11593425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
11603425bc38SStefano Zampini   /* put values into standard global vector */
11613425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11623425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11638eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
11643425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
11653425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
11663425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
11673425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
11683425bc38SStefano Zampini   }
11693425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11703425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11713425bc38SStefano Zampini   /* final change of basis if needed
11723425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
11733308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
11743425bc38SStefano Zampini   PetscFunctionReturn(0);
11753425bc38SStefano Zampini }
11761e6b0712SBarry Smith 
11773425bc38SStefano Zampini #undef __FUNCT__
11783425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
11793425bc38SStefano Zampini /*@
118028509bceSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
11813425bc38SStefano Zampini 
11823425bc38SStefano Zampini    Collective
11833425bc38SStefano Zampini 
11843425bc38SStefano Zampini    Input Parameters:
118528509bceSStefano Zampini +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
118628509bceSStefano Zampini .  fetidp_flux_sol - the solution of the FETIDP linear system
11873425bc38SStefano Zampini 
11883425bc38SStefano Zampini    Output Parameters:
118928509bceSStefano Zampini -  standard_sol      - the solution defined on the physical domain
11903425bc38SStefano Zampini 
11913425bc38SStefano Zampini    Level: developer
11923425bc38SStefano Zampini 
11933425bc38SStefano Zampini    Notes:
11943425bc38SStefano Zampini 
119528509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
11963425bc38SStefano Zampini @*/
11973425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
11983425bc38SStefano Zampini {
1199674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
12003425bc38SStefano Zampini   PetscErrorCode ierr;
12013425bc38SStefano Zampini 
12023425bc38SStefano Zampini   PetscFunctionBegin;
12033425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
12043425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
12053425bc38SStefano Zampini   PetscFunctionReturn(0);
12063425bc38SStefano Zampini }
12073425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
12081e6b0712SBarry Smith 
1209f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1210f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1211f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1212f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1213674ae819SStefano Zampini 
12143425bc38SStefano Zampini #undef __FUNCT__
12153425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
12163425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
12173425bc38SStefano Zampini {
1218674ae819SStefano Zampini 
1219674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
12203425bc38SStefano Zampini   Mat            newmat;
1221674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
12223425bc38SStefano Zampini   PC             newpc;
1223ce94432eSBarry Smith   MPI_Comm       comm;
12243425bc38SStefano Zampini   PetscErrorCode ierr;
12253425bc38SStefano Zampini 
12263425bc38SStefano Zampini   PetscFunctionBegin;
1227ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
12283425bc38SStefano Zampini   /* FETIDP linear matrix */
12293425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
12303425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
12313425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
12323425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
12333425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
12343425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
12353425bc38SStefano Zampini   /* FETIDP preconditioner */
12363425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
12373425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
12383425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
12393425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
12403425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
12413425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
12423425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
12433425bc38SStefano Zampini   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
12443425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
12453425bc38SStefano Zampini   /* return pointers for objects created */
12463425bc38SStefano Zampini   *fetidp_mat=newmat;
12473425bc38SStefano Zampini   *fetidp_pc=newpc;
12483425bc38SStefano Zampini   PetscFunctionReturn(0);
12493425bc38SStefano Zampini }
12501e6b0712SBarry Smith 
12513425bc38SStefano Zampini #undef __FUNCT__
12523425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
12533425bc38SStefano Zampini /*@
125428509bceSStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
12553425bc38SStefano Zampini 
12563425bc38SStefano Zampini    Collective
12573425bc38SStefano Zampini 
12583425bc38SStefano Zampini    Input Parameters:
125928509bceSStefano Zampini +  pc - the BDDC preconditioning context already setup
126028509bceSStefano Zampini 
126128509bceSStefano Zampini    Output Parameters:
126228509bceSStefano Zampini .  fetidp_mat - shell FETIDP matrix object
126328509bceSStefano Zampini .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
126428509bceSStefano Zampini 
126528509bceSStefano Zampini    Options Database Keys:
126628509bceSStefano Zampini -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
12673425bc38SStefano Zampini 
12683425bc38SStefano Zampini    Level: developer
12693425bc38SStefano Zampini 
12703425bc38SStefano Zampini    Notes:
127128509bceSStefano Zampini      Currently the only operation provided for FETIDP matrix is MatMult
12723425bc38SStefano Zampini 
12733425bc38SStefano Zampini .seealso: PCBDDC
12743425bc38SStefano Zampini @*/
12753425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
12763425bc38SStefano Zampini {
12773425bc38SStefano Zampini   PetscErrorCode ierr;
12783425bc38SStefano Zampini 
12793425bc38SStefano Zampini   PetscFunctionBegin;
12803425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12813425bc38SStefano Zampini   if (pc->setupcalled) {
1282516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1283f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
12843425bc38SStefano Zampini   PetscFunctionReturn(0);
12853425bc38SStefano Zampini }
12860c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1287da1bb401SStefano Zampini /*MC
1288da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
12890c7d97c5SJed Brown 
129028509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
129128509bceSStefano Zampini 
129228509bceSStefano Zampini .vb
129328509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
129428509bceSStefano 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
129528509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
129628509bceSStefano Zampini .ve
129728509bceSStefano Zampini 
129828509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
129928509bceSStefano Zampini 
1300b6fdb6dfSStefano Zampini    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
130128509bceSStefano Zampini 
130228509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
130328509bceSStefano Zampini 
1304b6fdb6dfSStefano 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.
1305b6fdb6dfSStefano Zampini 
130628509bceSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
130728509bceSStefano Zampini 
130828509bceSStefano 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
130928509bceSStefano Zampini 
131028509bceSStefano Zampini    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace.
131128509bceSStefano Zampini 
1312b6fdb6dfSStefano 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.
131328509bceSStefano Zampini 
131428509bceSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
131528509bceSStefano Zampini 
1316da1bb401SStefano Zampini    Options Database Keys:
131728509bceSStefano Zampini 
131828509bceSStefano Zampini .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
131928509bceSStefano Zampini .    -pc_bddc_use_edges <1> - use or not edges in primal space
1320b6fdb6dfSStefano Zampini .    -pc_bddc_use_faces <0> - use or not faces in primal space
132128509bceSStefano Zampini .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
132228509bceSStefano Zampini .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
132328509bceSStefano Zampini .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
132428509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
132528509bceSStefano Zampini .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
132628509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
132728509bceSStefano Zampini 
132828509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
132928509bceSStefano Zampini .vb
133028509bceSStefano Zampini       -pc_bddc_dirichlet_
133128509bceSStefano Zampini       -pc_bddc_neumann_
133228509bceSStefano Zampini       -pc_bddc_coarse_
133328509bceSStefano Zampini .ve
133428509bceSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
133528509bceSStefano Zampini 
133628509bceSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level can be specified as
133728509bceSStefano Zampini .vb
133828509bceSStefano Zampini       -pc_bddc_dirichlet_N_
133928509bceSStefano Zampini       -pc_bddc_neumann_N_
134028509bceSStefano Zampini       -pc_bddc_coarse_N_
134128509bceSStefano Zampini .ve
134228509bceSStefano Zampini    Note that level number ranges from the finest 0 to the coarsest N
1343da1bb401SStefano Zampini 
1344da1bb401SStefano Zampini    Level: intermediate
1345da1bb401SStefano Zampini 
1346b6fdb6dfSStefano Zampini    Developer notes:
134728509bceSStefano Zampini      Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose
1348da1bb401SStefano Zampini 
134928509bceSStefano Zampini      New deluxe scaling operator will be available soon.
1350da1bb401SStefano Zampini 
1351da1bb401SStefano Zampini    Contributed by Stefano Zampini
1352da1bb401SStefano Zampini 
1353da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1354da1bb401SStefano Zampini M*/
1355b2573a8aSBarry Smith 
1356da1bb401SStefano Zampini #undef __FUNCT__
1357da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
13588cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1359da1bb401SStefano Zampini {
1360da1bb401SStefano Zampini   PetscErrorCode      ierr;
1361da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1362da1bb401SStefano Zampini 
1363da1bb401SStefano Zampini   PetscFunctionBegin;
1364da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1365da1bb401SStefano Zampini   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1366da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1367da1bb401SStefano Zampini 
1368da1bb401SStefano Zampini   /* create PCIS data structure */
1369da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1370da1bb401SStefano Zampini 
137147d04d0dSStefano Zampini   /* BDDC customization */
137247d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
137347d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
137447d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
137547d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
137647d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
137747d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
137847d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
137947d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
138047d04d0dSStefano Zampini 
1381*727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
1382*727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
1383fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
1384f4ddd8eeSStefano Zampini   pcbddc->coarse_size                = 0;
1385f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
1386*727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
1387f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
1388f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
1389f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
1390674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
13910bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
13923972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1393534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1394534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1395534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1396da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1397da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1398da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1399da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1400da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
140115aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
140215aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1403da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1404da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1405da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1406da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1407da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1408da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1409da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1410da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1411da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1412da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
141360d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
141460d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
1415da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
1416da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
141785c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
141847d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
141947d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
14204fad6a16SStefano Zampini   pcbddc->current_level              = 0;
14212b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1422da1bb401SStefano Zampini 
1423674ae819SStefano Zampini   /* create local graph structure */
1424674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1425674ae819SStefano Zampini 
1426674ae819SStefano Zampini   /* scaling */
1427674ae819SStefano Zampini   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1428674ae819SStefano Zampini   pcbddc->work_scaling               = 0;
1429da1bb401SStefano Zampini 
1430da1bb401SStefano Zampini   /* function pointers */
1431da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
1432da1bb401SStefano Zampini   pc->ops->applytranspose      = 0;
1433da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1434da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1435da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1436da1bb401SStefano Zampini   pc->ops->view                = 0;
1437da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1438da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1439da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1440534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1441534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1442da1bb401SStefano Zampini 
1443da1bb401SStefano Zampini   /* composing function */
1444674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1445bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
14462b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1447b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
14482b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1449bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1450bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1451bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1452bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1453bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1454bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1455bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1456bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1457bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1458bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1459da1bb401SStefano Zampini   PetscFunctionReturn(0);
1460da1bb401SStefano Zampini }
14613425bc38SStefano Zampini 
1462