xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision c1c8e73620870d0684eaac3c19d9c866a35aa888)
153cdbc3dSStefano Zampini /* TODOLIST
2eb97c9d2SStefano Zampini 
3eb97c9d2SStefano Zampini    ConstraintsSetup
4eb97c9d2SStefano Zampini    - tolerances for constraints as an option (take care of single precision!)
5*c1c8e736SStefano Zampini    - Can MAT_IGNORE_ZERO_ENTRIES be used for Constraints Matrix?
6eb97c9d2SStefano Zampini 
7eb97c9d2SStefano Zampini    Solvers
8eb97c9d2SStefano Zampini    - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers)
9eb97c9d2SStefano Zampini    - Propagate ksp prefixes for solvers to mat objects?
10eb97c9d2SStefano Zampini    - Propagate nearnullspace info among levels
11eb97c9d2SStefano Zampini 
12eb97c9d2SStefano Zampini    User interface
13eb97c9d2SStefano Zampini    - Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet)
14*c1c8e736SStefano Zampini    - Negative indices in dirichlet and Neumann ISs should be skipped (now they cause out-of-bounds access)
15eb97c9d2SStefano Zampini    - Provide PCApplyTranpose_BDDC
16eb97c9d2SStefano Zampini    - DofSplitting and DM attached to pc?
17eb97c9d2SStefano Zampini 
18eb97c9d2SStefano Zampini    Debugging output
19eb97c9d2SStefano Zampini    - Better management of verbosity levels of debugging output
20eb97c9d2SStefano Zampini 
21eb97c9d2SStefano Zampini    Build
22eb97c9d2SStefano Zampini    - make runexe59
23eb97c9d2SStefano Zampini 
24eb97c9d2SStefano Zampini    Extra
25eb97c9d2SStefano Zampini    - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
26eb97c9d2SStefano Zampini    - Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver?
27eb97c9d2SStefano Zampini    - add support for computing h,H and related using coordinates?
28c0b83709SStefano Zampini    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
29eb97c9d2SStefano Zampini    - Better management in PCIS code
30eb97c9d2SStefano Zampini    - BDDC with MG framework?
31eb97c9d2SStefano Zampini 
32eb97c9d2SStefano Zampini    FETIDP
33eb97c9d2SStefano Zampini    - Move FETIDP code to its own classes
34eb97c9d2SStefano Zampini 
35eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
36eb97c9d2SStefano Zampini    - Provide general case for subassembling
3739e2fb2aSStefano Zampini    - Preallocation routines in MatISGetMPIAXAIJ
38eb97c9d2SStefano Zampini 
3953cdbc3dSStefano Zampini */
400c7d97c5SJed Brown 
4153cdbc3dSStefano Zampini /* ----------------------------------------------------------------------------------------------------------------------------------------------
420c7d97c5SJed Brown    Implementation of BDDC preconditioner based on:
430c7d97c5SJed Brown    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
4453cdbc3dSStefano Zampini    ---------------------------------------------------------------------------------------------------------------------------------------------- */
4553cdbc3dSStefano Zampini 
46674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/  /* includes for fortran wrappers */
47674ae819SStefano Zampini #include "bddcprivate.h"
483b03a366Sstefano_zampini #include <petscblaslapack.h>
49674ae819SStefano Zampini 
500c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
510c7d97c5SJed Brown #undef __FUNCT__
520c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
530c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc)
540c7d97c5SJed Brown {
550c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
560c7d97c5SJed Brown   PetscErrorCode ierr;
570c7d97c5SJed Brown 
580c7d97c5SJed Brown   PetscFunctionBegin;
590c7d97c5SJed Brown   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
608eeda7d8SStefano Zampini   /* Verbose debugging */
618eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
628eeda7d8SStefano Zampini   /* Primal space cumstomization */
638eeda7d8SStefano 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);
648eeda7d8SStefano 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);
658eeda7d8SStefano 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);
668eeda7d8SStefano Zampini   /* Change of basis */
678eeda7d8SStefano 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);
688eeda7d8SStefano 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);
69674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
70674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
71674ae819SStefano Zampini   }
728eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
738eeda7d8SStefano 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);
740298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
752b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
76674ae819SStefano 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);
770c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
780c7d97c5SJed Brown   PetscFunctionReturn(0);
790c7d97c5SJed Brown }
800c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
81674ae819SStefano Zampini #undef __FUNCT__
82674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
83674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
84674ae819SStefano Zampini {
85674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
86674ae819SStefano Zampini   PetscErrorCode ierr;
871e6b0712SBarry Smith 
88674ae819SStefano Zampini   PetscFunctionBegin;
89674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
90674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
91674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
92674ae819SStefano Zampini   PetscFunctionReturn(0);
93674ae819SStefano Zampini }
94674ae819SStefano Zampini #undef __FUNCT__
95674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
96674ae819SStefano Zampini /*@
9728509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
98674ae819SStefano Zampini 
99674ae819SStefano Zampini    Not collective
100674ae819SStefano Zampini 
101674ae819SStefano Zampini    Input Parameters:
102674ae819SStefano Zampini +  pc - the preconditioning context
10328509bceSStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering
104674ae819SStefano Zampini 
105674ae819SStefano Zampini    Level: intermediate
106674ae819SStefano Zampini 
107674ae819SStefano Zampini    Notes:
108674ae819SStefano Zampini 
109674ae819SStefano Zampini .seealso: PCBDDC
110674ae819SStefano Zampini @*/
111674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
112674ae819SStefano Zampini {
113674ae819SStefano Zampini   PetscErrorCode ierr;
114674ae819SStefano Zampini 
115674ae819SStefano Zampini   PetscFunctionBegin;
116674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
117674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
118674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
119674ae819SStefano Zampini   PetscFunctionReturn(0);
120674ae819SStefano Zampini }
121674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1220c7d97c5SJed Brown #undef __FUNCT__
1234fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1244fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1254fad6a16SStefano Zampini {
1264fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1274fad6a16SStefano Zampini 
1284fad6a16SStefano Zampini   PetscFunctionBegin;
1294fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1304fad6a16SStefano Zampini   PetscFunctionReturn(0);
1314fad6a16SStefano Zampini }
1321e6b0712SBarry Smith 
1334fad6a16SStefano Zampini #undef __FUNCT__
1344fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1354fad6a16SStefano Zampini /*@
13628509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
1374fad6a16SStefano Zampini 
1384fad6a16SStefano Zampini    Logically collective on PC
1394fad6a16SStefano Zampini 
1404fad6a16SStefano Zampini    Input Parameters:
1414fad6a16SStefano Zampini +  pc - the preconditioning context
14228509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
1434fad6a16SStefano Zampini 
14428509bceSStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
1454fad6a16SStefano Zampini 
1464fad6a16SStefano Zampini    Level: intermediate
1474fad6a16SStefano Zampini 
1484fad6a16SStefano Zampini    Notes:
1494fad6a16SStefano Zampini 
1504fad6a16SStefano Zampini .seealso: PCBDDC
1514fad6a16SStefano Zampini @*/
1524fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
1534fad6a16SStefano Zampini {
1544fad6a16SStefano Zampini   PetscErrorCode ierr;
1554fad6a16SStefano Zampini 
1564fad6a16SStefano Zampini   PetscFunctionBegin;
1574fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1582b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
1594fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
1604fad6a16SStefano Zampini   PetscFunctionReturn(0);
1614fad6a16SStefano Zampini }
1622b510759SStefano Zampini 
163b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
1642b510759SStefano Zampini #undef __FUNCT__
165b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
166b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
167b8ffe317SStefano Zampini {
168b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
169b8ffe317SStefano Zampini 
170b8ffe317SStefano Zampini   PetscFunctionBegin;
17185c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
172b8ffe317SStefano Zampini   PetscFunctionReturn(0);
173b8ffe317SStefano Zampini }
174b8ffe317SStefano Zampini 
175b8ffe317SStefano Zampini #undef __FUNCT__
176b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
177b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
1782b510759SStefano Zampini {
1792b510759SStefano Zampini   PetscErrorCode ierr;
1802b510759SStefano Zampini 
1812b510759SStefano Zampini   PetscFunctionBegin;
1822b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
183b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
184b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1852b510759SStefano Zampini   PetscFunctionReturn(0);
1862b510759SStefano Zampini }
1871e6b0712SBarry Smith 
1884fad6a16SStefano Zampini #undef __FUNCT__
1892b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
1902b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
1914fad6a16SStefano Zampini {
1924fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1934fad6a16SStefano Zampini 
1944fad6a16SStefano Zampini   PetscFunctionBegin;
1952b510759SStefano Zampini   pcbddc->current_level = level;
1964fad6a16SStefano Zampini   PetscFunctionReturn(0);
1974fad6a16SStefano Zampini }
1981e6b0712SBarry Smith 
1994fad6a16SStefano Zampini #undef __FUNCT__
200b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
201b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
202b8ffe317SStefano Zampini {
203b8ffe317SStefano Zampini   PetscErrorCode ierr;
204b8ffe317SStefano Zampini 
205b8ffe317SStefano Zampini   PetscFunctionBegin;
206b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
207b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
208b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
209b8ffe317SStefano Zampini   PetscFunctionReturn(0);
210b8ffe317SStefano Zampini }
211b8ffe317SStefano Zampini 
212b8ffe317SStefano Zampini #undef __FUNCT__
2132b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
2142b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
2152b510759SStefano Zampini {
2162b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2172b510759SStefano Zampini 
2182b510759SStefano Zampini   PetscFunctionBegin;
2192b510759SStefano Zampini   pcbddc->max_levels = levels;
2202b510759SStefano Zampini   PetscFunctionReturn(0);
2212b510759SStefano Zampini }
2222b510759SStefano Zampini 
2232b510759SStefano Zampini #undef __FUNCT__
2242b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2254fad6a16SStefano Zampini /*@
22628509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
2274fad6a16SStefano Zampini 
2284fad6a16SStefano Zampini    Logically collective on PC
2294fad6a16SStefano Zampini 
2304fad6a16SStefano Zampini    Input Parameters:
2314fad6a16SStefano Zampini +  pc - the preconditioning context
23228509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
2334fad6a16SStefano Zampini 
23428509bceSStefano Zampini    Default value is 0, i.e. traditional one-level BDDC
2354fad6a16SStefano Zampini 
2364fad6a16SStefano Zampini    Level: intermediate
2374fad6a16SStefano Zampini 
2384fad6a16SStefano Zampini    Notes:
2394fad6a16SStefano Zampini 
2404fad6a16SStefano Zampini .seealso: PCBDDC
2414fad6a16SStefano Zampini @*/
2422b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
2434fad6a16SStefano Zampini {
2444fad6a16SStefano Zampini   PetscErrorCode ierr;
2454fad6a16SStefano Zampini 
2464fad6a16SStefano Zampini   PetscFunctionBegin;
2474fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2482b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
2492b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
2504fad6a16SStefano Zampini   PetscFunctionReturn(0);
2514fad6a16SStefano Zampini }
2524fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
2531e6b0712SBarry Smith 
2544fad6a16SStefano Zampini #undef __FUNCT__
2550bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
2560bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
2570bdf917eSStefano Zampini {
2580bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2590bdf917eSStefano Zampini   PetscErrorCode ierr;
2600bdf917eSStefano Zampini 
2610bdf917eSStefano Zampini   PetscFunctionBegin;
2620bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
2630bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
2640bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
2650bdf917eSStefano Zampini   PetscFunctionReturn(0);
2660bdf917eSStefano Zampini }
2671e6b0712SBarry Smith 
2680bdf917eSStefano Zampini #undef __FUNCT__
2690bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
2700bdf917eSStefano Zampini /*@
27128509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
2720bdf917eSStefano Zampini 
2730bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
2740bdf917eSStefano Zampini 
2750bdf917eSStefano Zampini    Input Parameters:
2760bdf917eSStefano Zampini +  pc - the preconditioning context
27728509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
2780bdf917eSStefano Zampini 
2790bdf917eSStefano Zampini    Level: intermediate
2800bdf917eSStefano Zampini 
2810bdf917eSStefano Zampini    Notes:
2820bdf917eSStefano Zampini 
2830bdf917eSStefano Zampini .seealso: PCBDDC
2840bdf917eSStefano Zampini @*/
2850bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
2860bdf917eSStefano Zampini {
2870bdf917eSStefano Zampini   PetscErrorCode ierr;
2880bdf917eSStefano Zampini 
2890bdf917eSStefano Zampini   PetscFunctionBegin;
2900bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
291674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
2922b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
2930bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
2940bdf917eSStefano Zampini   PetscFunctionReturn(0);
2950bdf917eSStefano Zampini }
2960bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
2971e6b0712SBarry Smith 
2980bdf917eSStefano Zampini #undef __FUNCT__
2993b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
3003b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
3013b03a366Sstefano_zampini {
3023b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3033b03a366Sstefano_zampini   PetscErrorCode ierr;
3043b03a366Sstefano_zampini 
3053b03a366Sstefano_zampini   PetscFunctionBegin;
3063b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
30736e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
30836e030ebSStefano Zampini   pcbddc->DirichletBoundaries=DirichletBoundaries;
309fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3103b03a366Sstefano_zampini   PetscFunctionReturn(0);
3113b03a366Sstefano_zampini }
3121e6b0712SBarry Smith 
3133b03a366Sstefano_zampini #undef __FUNCT__
3143b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
3153b03a366Sstefano_zampini /*@
31628509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
3173b03a366Sstefano_zampini 
3183b03a366Sstefano_zampini    Not collective
3193b03a366Sstefano_zampini 
3203b03a366Sstefano_zampini    Input Parameters:
3213b03a366Sstefano_zampini +  pc - the preconditioning context
32228509bceSStefano Zampini -  DirichletBoundaries - sequential IS defining the subdomain part of Dirichlet boundaries (in local ordering)
3233b03a366Sstefano_zampini 
3243b03a366Sstefano_zampini    Level: intermediate
3253b03a366Sstefano_zampini 
3263b03a366Sstefano_zampini    Notes:
3273b03a366Sstefano_zampini 
3283b03a366Sstefano_zampini .seealso: PCBDDC
3293b03a366Sstefano_zampini @*/
3303b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3313b03a366Sstefano_zampini {
3323b03a366Sstefano_zampini   PetscErrorCode ierr;
3333b03a366Sstefano_zampini 
3343b03a366Sstefano_zampini   PetscFunctionBegin;
3353b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
336674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
3373b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
3383b03a366Sstefano_zampini   PetscFunctionReturn(0);
3393b03a366Sstefano_zampini }
3403b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
3411e6b0712SBarry Smith 
3423b03a366Sstefano_zampini #undef __FUNCT__
3430c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
34453cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
3450c7d97c5SJed Brown {
3460c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
34753cdbc3dSStefano Zampini   PetscErrorCode ierr;
3480c7d97c5SJed Brown 
3490c7d97c5SJed Brown   PetscFunctionBegin;
35053cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
35136e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
35236e030ebSStefano Zampini   pcbddc->NeumannBoundaries=NeumannBoundaries;
353fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3540c7d97c5SJed Brown   PetscFunctionReturn(0);
3550c7d97c5SJed Brown }
3561e6b0712SBarry Smith 
3570c7d97c5SJed Brown #undef __FUNCT__
3580c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
35957527edcSJed Brown /*@
36028509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
36157527edcSJed Brown 
3629c0446d6SStefano Zampini    Not collective
36357527edcSJed Brown 
36457527edcSJed Brown    Input Parameters:
36557527edcSJed Brown +  pc - the preconditioning context
36628509bceSStefano Zampini -  NeumannBoundaries - sequential IS defining the subdomain part of Neumann boundaries (in local ordering)
36757527edcSJed Brown 
36857527edcSJed Brown    Level: intermediate
36957527edcSJed Brown 
37057527edcSJed Brown    Notes:
37157527edcSJed Brown 
37257527edcSJed Brown .seealso: PCBDDC
37357527edcSJed Brown @*/
37453cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
3750c7d97c5SJed Brown {
3760c7d97c5SJed Brown   PetscErrorCode ierr;
3770c7d97c5SJed Brown 
3780c7d97c5SJed Brown   PetscFunctionBegin;
3790c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
380674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
38153cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
38253cdbc3dSStefano Zampini   PetscFunctionReturn(0);
38353cdbc3dSStefano Zampini }
38453cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
3851e6b0712SBarry Smith 
38653cdbc3dSStefano Zampini #undef __FUNCT__
387da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
388da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
389da1bb401SStefano Zampini {
390da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
391da1bb401SStefano Zampini 
392da1bb401SStefano Zampini   PetscFunctionBegin;
393da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
394da1bb401SStefano Zampini   PetscFunctionReturn(0);
395da1bb401SStefano Zampini }
3961e6b0712SBarry Smith 
397da1bb401SStefano Zampini #undef __FUNCT__
398da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
399da1bb401SStefano Zampini /*@
40028509bceSStefano Zampini  PCBDDCGetDirichletBoundaries - Get IS for local Dirichlet boundaries
401da1bb401SStefano Zampini 
402da1bb401SStefano Zampini    Not collective
403da1bb401SStefano Zampini 
404da1bb401SStefano Zampini    Input Parameters:
40528509bceSStefano Zampini .  pc - the preconditioning context
406da1bb401SStefano Zampini 
407da1bb401SStefano Zampini    Output Parameters:
40828509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
409da1bb401SStefano Zampini 
410da1bb401SStefano Zampini    Level: intermediate
411da1bb401SStefano Zampini 
412da1bb401SStefano Zampini    Notes:
413da1bb401SStefano Zampini 
414da1bb401SStefano Zampini .seealso: PCBDDC
415da1bb401SStefano Zampini @*/
416da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
417da1bb401SStefano Zampini {
418da1bb401SStefano Zampini   PetscErrorCode ierr;
419da1bb401SStefano Zampini 
420da1bb401SStefano Zampini   PetscFunctionBegin;
421da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
422da1bb401SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
423da1bb401SStefano Zampini   PetscFunctionReturn(0);
424da1bb401SStefano Zampini }
425da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
4261e6b0712SBarry Smith 
427da1bb401SStefano Zampini #undef __FUNCT__
42853cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
42953cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
43053cdbc3dSStefano Zampini {
43153cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
43253cdbc3dSStefano Zampini 
43353cdbc3dSStefano Zampini   PetscFunctionBegin;
43453cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
43553cdbc3dSStefano Zampini   PetscFunctionReturn(0);
43653cdbc3dSStefano Zampini }
4371e6b0712SBarry Smith 
43853cdbc3dSStefano Zampini #undef __FUNCT__
43953cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
44053cdbc3dSStefano Zampini /*@
44128509bceSStefano Zampini  PCBDDCGetNeumannBoundaries - Get IS for local Neumann boundaries
44253cdbc3dSStefano Zampini 
4439c0446d6SStefano Zampini    Not collective
44453cdbc3dSStefano Zampini 
44553cdbc3dSStefano Zampini    Input Parameters:
44628509bceSStefano Zampini .  pc - the preconditioning context
44753cdbc3dSStefano Zampini 
44853cdbc3dSStefano Zampini    Output Parameters:
44928509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
45053cdbc3dSStefano Zampini 
45153cdbc3dSStefano Zampini    Level: intermediate
45253cdbc3dSStefano Zampini 
45353cdbc3dSStefano Zampini    Notes:
45453cdbc3dSStefano Zampini 
45553cdbc3dSStefano Zampini .seealso: PCBDDC
45653cdbc3dSStefano Zampini @*/
45753cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
45853cdbc3dSStefano Zampini {
45953cdbc3dSStefano Zampini   PetscErrorCode ierr;
46053cdbc3dSStefano Zampini 
46153cdbc3dSStefano Zampini   PetscFunctionBegin;
46253cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
46353cdbc3dSStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
4640c7d97c5SJed Brown   PetscFunctionReturn(0);
4650c7d97c5SJed Brown }
46636e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
4671e6b0712SBarry Smith 
46836e030ebSStefano Zampini #undef __FUNCT__
469da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
4701a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
47136e030ebSStefano Zampini {
47236e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
473da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
474da1bb401SStefano Zampini   PetscErrorCode ierr;
47536e030ebSStefano Zampini 
47636e030ebSStefano Zampini   PetscFunctionBegin;
477674ae819SStefano Zampini   /* free old CSR */
478674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
479fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
480674ae819SStefano Zampini   /* get CSR into graph structure */
481da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
482674ae819SStefano Zampini     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
483674ae819SStefano Zampini     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
484674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
485674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
486da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
4871a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
4881a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
489674ae819SStefano Zampini   }
490575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
49136e030ebSStefano Zampini   PetscFunctionReturn(0);
49236e030ebSStefano Zampini }
4931e6b0712SBarry Smith 
49436e030ebSStefano Zampini #undef __FUNCT__
495da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
49636e030ebSStefano Zampini /*@
49728509bceSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
49836e030ebSStefano Zampini 
49936e030ebSStefano Zampini    Not collective
50036e030ebSStefano Zampini 
50136e030ebSStefano Zampini    Input Parameters:
50236e030ebSStefano Zampini +  pc - the preconditioning context
50328509bceSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
50428509bceSStefano Zampini .  xadj, adjncy - the CSR graph
50528509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
50636e030ebSStefano Zampini 
50736e030ebSStefano Zampini    Level: intermediate
50836e030ebSStefano Zampini 
50936e030ebSStefano Zampini    Notes:
51036e030ebSStefano Zampini 
51128509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
51236e030ebSStefano Zampini @*/
5131a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
51436e030ebSStefano Zampini {
515575ad6abSStefano Zampini   void (*f)(void) = 0;
51636e030ebSStefano Zampini   PetscErrorCode ierr;
51736e030ebSStefano Zampini 
51836e030ebSStefano Zampini   PetscFunctionBegin;
51936e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
520674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
521674ae819SStefano Zampini   PetscValidIntPointer(xadj,4);
522674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
523674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
524da1bb401SStefano Zampini   }
52536e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
526575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
527575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
528575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
529575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
530575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
53136e030ebSStefano Zampini   }
53236e030ebSStefano Zampini   PetscFunctionReturn(0);
53336e030ebSStefano Zampini }
5349c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
5351e6b0712SBarry Smith 
5369c0446d6SStefano Zampini #undef __FUNCT__
5379c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
5389c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
5399c0446d6SStefano Zampini {
5409c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5419c0446d6SStefano Zampini   PetscInt i;
5429c0446d6SStefano Zampini   PetscErrorCode ierr;
5439c0446d6SStefano Zampini 
5449c0446d6SStefano Zampini   PetscFunctionBegin;
545da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
5469c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
5479c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
5489c0446d6SStefano Zampini   }
549d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
550da1bb401SStefano Zampini   /* allocate space then set */
5519c0446d6SStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
5529c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
553da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
554da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
5559c0446d6SStefano Zampini   }
5569c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
55760d44989SStefano Zampini   pcbddc->user_provided_isfordofs = PETSC_TRUE;
5589c0446d6SStefano Zampini   PetscFunctionReturn(0);
5599c0446d6SStefano Zampini }
5601e6b0712SBarry Smith 
5619c0446d6SStefano Zampini #undef __FUNCT__
5629c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
5639c0446d6SStefano Zampini /*@
56428509bceSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the local Neumann matrix
5659c0446d6SStefano Zampini 
5669c0446d6SStefano Zampini    Not collective
5679c0446d6SStefano Zampini 
5689c0446d6SStefano Zampini    Input Parameters:
5699c0446d6SStefano Zampini +  pc - the preconditioning context
57028509bceSStefano Zampini -  n_is - number of index sets defining the fields
57128509bceSStefano Zampini .  ISForDofs - array of IS describing the fields
5729c0446d6SStefano Zampini 
5739c0446d6SStefano Zampini    Level: intermediate
5749c0446d6SStefano Zampini 
5759c0446d6SStefano Zampini    Notes:
5769c0446d6SStefano Zampini 
5779c0446d6SStefano Zampini .seealso: PCBDDC
5789c0446d6SStefano Zampini @*/
5799c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
5809c0446d6SStefano Zampini {
5812b510759SStefano Zampini   PetscInt       i;
5829c0446d6SStefano Zampini   PetscErrorCode ierr;
5839c0446d6SStefano Zampini 
5849c0446d6SStefano Zampini   PetscFunctionBegin;
5859c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5862b510759SStefano Zampini   for (i=0;i<n_is;i++) {
5872b510759SStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2);
5882b510759SStefano Zampini   }
5899c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
5909c0446d6SStefano Zampini   PetscFunctionReturn(0);
5919c0446d6SStefano Zampini }
592da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
593534831adSStefano Zampini #undef __FUNCT__
594534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
595534831adSStefano Zampini /* -------------------------------------------------------------------------- */
596534831adSStefano Zampini /*
597534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
598534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
5999c0446d6SStefano Zampini 
600534831adSStefano Zampini    Input Parameter:
601534831adSStefano Zampini +  pc - the preconditioner contex
602534831adSStefano Zampini 
603534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
604534831adSStefano Zampini 
605534831adSStefano Zampini    Notes:
606534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
607534831adSStefano Zampini    the user, but instead is called by KSPSolve().
608534831adSStefano Zampini */
609534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
610534831adSStefano Zampini {
611534831adSStefano Zampini   PetscErrorCode ierr;
612534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
613534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
614534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
615534831adSStefano Zampini   Mat            temp_mat;
6163972b0daSStefano Zampini   IS             dirIS;
6173972b0daSStefano Zampini   PetscInt       dirsize,i,*is_indices;
6183972b0daSStefano Zampini   PetscScalar    *array_x,*array_diagonal;
6193972b0daSStefano Zampini   Vec            used_vec;
62092e3dcfbSStefano Zampini   PetscBool      guess_nonzero,flg,bddc_has_dirichlet_boundaries;
621534831adSStefano Zampini 
622534831adSStefano Zampini   PetscFunctionBegin;
62385c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
62485c4d303SStefano Zampini   if (ksp) {
62585c4d303SStefano Zampini     PetscBool iscg;
62685c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
62785c4d303SStefano Zampini     if (!iscg) {
62885c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
62985c4d303SStefano Zampini     }
63085c4d303SStefano Zampini   }
63185c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
63262a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
63362a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
63462a6ff1dSStefano Zampini   }
63562a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
63662a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
63762a6ff1dSStefano Zampini   }
6383972b0daSStefano Zampini   if (x) {
6393972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
6403972b0daSStefano Zampini     used_vec = x;
6413972b0daSStefano Zampini   } else {
6423972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
6433972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
6443972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6453972b0daSStefano Zampini   }
6463972b0daSStefano Zampini   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
6473972b0daSStefano Zampini   if (ksp) {
6483972b0daSStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
6493972b0daSStefano Zampini     if (!guess_nonzero) {
6503972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6513972b0daSStefano Zampini     }
6523972b0daSStefano Zampini   }
6533308cffdSStefano Zampini 
65492e3dcfbSStefano Zampini   /* TODO: remove when Dirichlet boundaries will be shared */
65592e3dcfbSStefano Zampini   ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr);
65692e3dcfbSStefano Zampini   flg = PETSC_FALSE;
65792e3dcfbSStefano Zampini   if (dirIS) flg = PETSC_TRUE;
65892e3dcfbSStefano Zampini   ierr = MPI_Allreduce(&flg,&bddc_has_dirichlet_boundaries,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
65992e3dcfbSStefano Zampini 
6603972b0daSStefano Zampini   /* store the original rhs */
6613972b0daSStefano Zampini   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
6623972b0daSStefano Zampini 
6633972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
66492e3dcfbSStefano Zampini   if (rhs && bddc_has_dirichlet_boundaries) {
6653972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
6663972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
6673972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6683972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6693972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6703972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6713972b0daSStefano Zampini     if (dirIS) {
6723972b0daSStefano Zampini       ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr);
6733972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6743972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6753972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6762fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
6773972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6783972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6793972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6803972b0daSStefano Zampini     }
6813972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6823972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
683b76ba322SStefano Zampini 
6843972b0daSStefano Zampini     /* remove the computed solution from the rhs */
6853972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6863972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
6873972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6883308cffdSStefano Zampini   }
689b76ba322SStefano Zampini 
690b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
6913972b0daSStefano Zampini   if (x) {
6923972b0daSStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
6933972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
69485c4d303SStefano Zampini     if (pcbddc->use_exact_dirichlet_trick) {
695b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
696b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
697b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
698b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
699b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
700b76ba322SStefano Zampini       if (ksp) {
701b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
702b76ba322SStefano Zampini       }
703b76ba322SStefano Zampini     }
7043972b0daSStefano Zampini   }
705b76ba322SStefano Zampini 
70692e3dcfbSStefano Zampini   /* prepare MatMult and rhs for solver */
707674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
708b76ba322SStefano Zampini     /* swap pointers for local matrices */
709b76ba322SStefano Zampini     temp_mat = matis->A;
710b76ba322SStefano Zampini     matis->A = pcbddc->local_mat;
711b76ba322SStefano Zampini     pcbddc->local_mat = temp_mat;
71292e3dcfbSStefano Zampini     if (rhs) {
713b76ba322SStefano Zampini       /* Get local rhs and apply transformation of basis */
714b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
715b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
716b76ba322SStefano Zampini       /* from original basis to modified basis */
717b76ba322SStefano Zampini       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
718b76ba322SStefano Zampini       /* put back modified values into the global vec using INSERT_VALUES copy mode */
719b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
720b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
721674ae819SStefano Zampini     }
72292e3dcfbSStefano Zampini   }
72392e3dcfbSStefano Zampini 
72492e3dcfbSStefano Zampini   /* remove nullspace if present */
7250bdf917eSStefano Zampini   if (ksp && pcbddc->NullSpace) {
726d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
727d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
728b76ba322SStefano Zampini   }
7290bdf917eSStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
730534831adSStefano Zampini   PetscFunctionReturn(0);
731534831adSStefano Zampini }
732534831adSStefano Zampini /* -------------------------------------------------------------------------- */
733534831adSStefano Zampini #undef __FUNCT__
734534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
735534831adSStefano Zampini /* -------------------------------------------------------------------------- */
736534831adSStefano Zampini /*
737534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
738534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
739534831adSStefano Zampini 
740534831adSStefano Zampini    Input Parameter:
741534831adSStefano Zampini +  pc - the preconditioner contex
742534831adSStefano Zampini 
743534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
744534831adSStefano Zampini 
745534831adSStefano Zampini    Notes:
746534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
747534831adSStefano Zampini    the user, but instead is called by KSPSolve().
748534831adSStefano Zampini */
749534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
750534831adSStefano Zampini {
751534831adSStefano Zampini   PetscErrorCode ierr;
752534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
753534831adSStefano Zampini   PC_IS          *pcis   = (PC_IS*)(pc->data);
754534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
755534831adSStefano Zampini   Mat            temp_mat;
756534831adSStefano Zampini 
757534831adSStefano Zampini   PetscFunctionBegin;
758674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
759534831adSStefano Zampini     /* swap pointers for local matrices */
760534831adSStefano Zampini     temp_mat = matis->A;
761534831adSStefano Zampini     matis->A = pcbddc->local_mat;
762534831adSStefano Zampini     pcbddc->local_mat = temp_mat;
7633425bc38SStefano Zampini   }
7643308cffdSStefano Zampini   if (pcbddc->use_change_of_basis && x) {
765534831adSStefano Zampini     /* Get Local boundary and apply transformation of basis to solution vector */
766534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
767534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
768534831adSStefano Zampini     /* from modified basis to original basis */
769534831adSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
770534831adSStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
771534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
772534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
773534831adSStefano Zampini   }
7743972b0daSStefano Zampini   /* add solution removed in presolve */
7753425bc38SStefano Zampini   if (x) {
7763425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
7773425bc38SStefano Zampini   }
778fb223d50SStefano Zampini   /* restore rhs to its original state */
779fb223d50SStefano Zampini   if (rhs) {
780fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
781fb223d50SStefano Zampini   }
782534831adSStefano Zampini   PetscFunctionReturn(0);
783534831adSStefano Zampini }
784534831adSStefano Zampini /* -------------------------------------------------------------------------- */
78553cdbc3dSStefano Zampini #undef __FUNCT__
78653cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
7870c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
7880c7d97c5SJed Brown /*
7890c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
7900c7d97c5SJed Brown                   by setting data structures and options.
7910c7d97c5SJed Brown 
7920c7d97c5SJed Brown    Input Parameter:
79353cdbc3dSStefano Zampini +  pc - the preconditioner context
7940c7d97c5SJed Brown 
7950c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
7960c7d97c5SJed Brown 
7970c7d97c5SJed Brown    Notes:
7980c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
7990c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
8000c7d97c5SJed Brown */
80153cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
8020c7d97c5SJed Brown {
8030c7d97c5SJed Brown   PetscErrorCode   ierr;
8040c7d97c5SJed Brown   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
805f4ddd8eeSStefano Zampini   MatNullSpace     nearnullspace;
806f4ddd8eeSStefano Zampini   const Vec        *nearnullvecs,*onearnullvecs;
807674ae819SStefano Zampini   MatStructure     flag;
808f4ddd8eeSStefano Zampini   PetscObjectState state;
809f4ddd8eeSStefano Zampini   PetscInt         nnsp_size,onnsp_size;
810674ae819SStefano Zampini   PetscBool        computeis,computetopography,computesolvers;
811f4ddd8eeSStefano Zampini   PetscBool        new_nearnullspace_provided,nnsp_has_cnst,onnsp_has_cnst;
8120c7d97c5SJed Brown 
8130c7d97c5SJed Brown   PetscFunctionBegin;
814f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
815f4ddd8eeSStefano Zampini   /* PCIS does not support MatStructure flags different from SAME_PRECONDITIONER */
8163b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
817fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
818f4ddd8eeSStefano Zampini 
819f4ddd8eeSStefano Zampini   /* split work */
820674ae819SStefano Zampini   if (pc->setupcalled) {
821674ae819SStefano Zampini     computeis = PETSC_FALSE;
822674ae819SStefano Zampini     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
823674ae819SStefano Zampini     if (flag == SAME_PRECONDITIONER) {
824f4ddd8eeSStefano Zampini       PetscFunctionReturn(0);
825674ae819SStefano Zampini     } else if (flag == SAME_NONZERO_PATTERN) {
826674ae819SStefano Zampini       computetopography = PETSC_FALSE;
827674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
828674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
829674ae819SStefano Zampini       computetopography = PETSC_TRUE;
830674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
831674ae819SStefano Zampini     }
832674ae819SStefano Zampini   } else {
833674ae819SStefano Zampini     computeis = PETSC_TRUE;
834674ae819SStefano Zampini     computetopography = PETSC_TRUE;
835674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
836674ae819SStefano Zampini   }
837fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
838fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
839fb180af4SStefano Zampini   }
840f4ddd8eeSStefano Zampini 
841f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
842f4ddd8eeSStefano Zampini   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
843f4ddd8eeSStefano Zampini     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
844f4ddd8eeSStefano Zampini     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
845f4ddd8eeSStefano Zampini     if (pcbddc->current_level) {
846f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
847f4ddd8eeSStefano Zampini     }
848f4ddd8eeSStefano Zampini   }
849f4ddd8eeSStefano Zampini 
850fcd409f5SStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
851674ae819SStefano Zampini   if (computeis) {
852fcd409f5SStefano Zampini     /* HACK INTO PCIS */
853fcd409f5SStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
854fcd409f5SStefano Zampini     pcis->computesolvers = PETSC_FALSE;
855674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
85639e2fb2aSStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr);
857674ae819SStefano Zampini   }
858f4ddd8eeSStefano Zampini 
859f4ddd8eeSStefano Zampini   /* Analyze interface */
860674ae819SStefano Zampini   if (computetopography) {
861674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
862fb8d54d4SStefano Zampini   }
863fb8d54d4SStefano Zampini 
864f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
865fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
866f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
867f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
868f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
869f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
870f4ddd8eeSStefano Zampini     } else {
871f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
872f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
873f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
874f4ddd8eeSStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace */
875f4ddd8eeSStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
876f4ddd8eeSStefano Zampini         ierr = MatNullSpaceGetVecs(pcbddc->onearnullspace,&onnsp_has_cnst,&onnsp_size,&onearnullvecs);CHKERRQ(ierr);
877f4ddd8eeSStefano Zampini         if ( (nnsp_has_cnst != onnsp_has_cnst) || (nnsp_size != onnsp_size) ) {
878f4ddd8eeSStefano Zampini           new_nearnullspace_provided = PETSC_TRUE;
879f4ddd8eeSStefano Zampini         } else { /* nullspaces have the same size, so check vectors or their ObjectStateId */
880f4ddd8eeSStefano Zampini           PetscInt i;
881f4ddd8eeSStefano Zampini           for (i=0;i<nnsp_size;i++) {
882f4ddd8eeSStefano Zampini             ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
883f4ddd8eeSStefano Zampini             if (nearnullvecs[i] != onearnullvecs[i] || pcbddc->onearnullvecs_state[i] != state) {
884f4ddd8eeSStefano Zampini               new_nearnullspace_provided = PETSC_TRUE;
885f4ddd8eeSStefano Zampini               break;
886f4ddd8eeSStefano Zampini             }
887f4ddd8eeSStefano Zampini           }
888f4ddd8eeSStefano Zampini         }
889f4ddd8eeSStefano Zampini       }
890f4ddd8eeSStefano Zampini     }
891f4ddd8eeSStefano Zampini   } else {
892f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
893f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
894f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
895f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
896f4ddd8eeSStefano Zampini     }
897f4ddd8eeSStefano Zampini   }
898f4ddd8eeSStefano Zampini 
899f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
900727cdba6SStefano Zampini   /* reset primal space flags */
901f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
902727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
903fb8d54d4SStefano Zampini   if (computetopography || new_nearnullspace_provided) {
904727cdba6SStefano Zampini     /* It also sets the primal space flags */
905674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
906e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
907f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
908674ae819SStefano Zampini   }
909fb8d54d4SStefano Zampini 
910f4ddd8eeSStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
911674ae819SStefano Zampini     /* reset data */
912674ae819SStefano Zampini     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
91399cc7994SStefano Zampini     /* Create coarse and local stuffs */
91499cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
915674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
9160c7d97c5SJed Brown   }
9172b510759SStefano Zampini   if (pcbddc->dbg_flag && pcbddc->current_level) {
9182b510759SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
9192b510759SStefano Zampini   }
9200c7d97c5SJed Brown   PetscFunctionReturn(0);
9210c7d97c5SJed Brown }
9220c7d97c5SJed Brown 
9230c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
9240c7d97c5SJed Brown /*
9250c7d97c5SJed Brown    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
9260c7d97c5SJed Brown 
9270c7d97c5SJed Brown    Input Parameters:
9280c7d97c5SJed Brown .  pc - the preconditioner context
9290c7d97c5SJed Brown .  r - input vector (global)
9300c7d97c5SJed Brown 
9310c7d97c5SJed Brown    Output Parameter:
9320c7d97c5SJed Brown .  z - output vector (global)
9330c7d97c5SJed Brown 
9340c7d97c5SJed Brown    Application Interface Routine: PCApply()
9350c7d97c5SJed Brown  */
9360c7d97c5SJed Brown #undef __FUNCT__
9370c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
93853cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
9390c7d97c5SJed Brown {
9400c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
9410c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
9420c7d97c5SJed Brown   PetscErrorCode    ierr;
9433b03a366Sstefano_zampini   const PetscScalar one = 1.0;
9443b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
9452617d88aSStefano Zampini   const PetscScalar zero = 0.0;
9460c7d97c5SJed Brown 
9470c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
9480c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
9498eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
9500c7d97c5SJed Brown 
9510c7d97c5SJed Brown   PetscFunctionBegin;
95285c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
9530c7d97c5SJed Brown     /* First Dirichlet solve */
9540c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9550c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
95653cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
9570c7d97c5SJed Brown     /*
9580c7d97c5SJed Brown       Assembling right hand side for BDDC operator
959674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
960674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
9610c7d97c5SJed Brown     */
9620c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
9630c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
9648eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
9650c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
9660c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
9670c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9680c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
969674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
970b76ba322SStefano Zampini   } else {
9710bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
972b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
973674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
974b76ba322SStefano Zampini   }
975b76ba322SStefano Zampini 
9762617d88aSStefano Zampini   /* Apply interface preconditioner
9772617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
9782617d88aSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
9792617d88aSStefano Zampini 
980674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
981674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
9820c7d97c5SJed Brown 
9833b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
9840c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9850c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9860c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
9878eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
988df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
989df187020SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
990df187020SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY (pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
991df187020SStefano Zampini   ierr = VecAXPY (pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
9920c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9930c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9940c7d97c5SJed Brown   PetscFunctionReturn(0);
9950c7d97c5SJed Brown }
996da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
997674ae819SStefano Zampini 
998da1bb401SStefano Zampini #undef __FUNCT__
999da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1000da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1001da1bb401SStefano Zampini {
1002da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1003da1bb401SStefano Zampini   PetscErrorCode ierr;
1004da1bb401SStefano Zampini 
1005da1bb401SStefano Zampini   PetscFunctionBegin;
1006da1bb401SStefano Zampini   /* free data created by PCIS */
1007da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1008674ae819SStefano Zampini   /* free BDDC custom data  */
1009674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1010674ae819SStefano Zampini   /* destroy objects related to topography */
1011674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1012674ae819SStefano Zampini   /* free allocated graph structure */
1013da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1014674ae819SStefano Zampini   /* free data for scaling operator */
1015674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1016674ae819SStefano Zampini   /* free solvers stuff */
1017674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
101862a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
101962a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
102062a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
102139e2fb2aSStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr);
10223425bc38SStefano Zampini   /* remove functions */
1023674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1024bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
10252b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1026b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
10272b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1028bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1029bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1030bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1031bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
1032bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1033bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
1034bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1035bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1036bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1037bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1038674ae819SStefano Zampini   /* Free the private data structure */
1039674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1040da1bb401SStefano Zampini   PetscFunctionReturn(0);
1041da1bb401SStefano Zampini }
10423425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
10431e6b0712SBarry Smith 
10443425bc38SStefano Zampini #undef __FUNCT__
10453425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
10463425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
10473425bc38SStefano Zampini {
1048674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10493425bc38SStefano Zampini   PC_IS*         pcis;
10503425bc38SStefano Zampini   PC_BDDC*       pcbddc;
10513425bc38SStefano Zampini   PetscErrorCode ierr;
10520c7d97c5SJed Brown 
10533425bc38SStefano Zampini   PetscFunctionBegin;
10543425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10553425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
10563425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
10573425bc38SStefano Zampini 
10583425bc38SStefano Zampini   /* change of basis for physical rhs if needed
10593425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
10603308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
10613425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
10623425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10633425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1064fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1065fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
10663425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10673425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1068674ae819SStefano Zampini   /* Apply partition of unity */
10693425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1070674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
10718eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
10723425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
10733425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10743425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
10753425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
10763425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
10773425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10783425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1079674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
10803425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10813425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10823425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
10833425bc38SStefano Zampini   }
10843425bc38SStefano Zampini   /* BDDC rhs */
10853425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
10868eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
10873425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10883425bc38SStefano Zampini   }
10893425bc38SStefano Zampini   /* apply BDDC */
10903425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
10913425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
10923425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
10933425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
10943425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10953425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10963425bc38SStefano Zampini   /* restore original rhs */
10973425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
10983425bc38SStefano Zampini   PetscFunctionReturn(0);
10993425bc38SStefano Zampini }
11001e6b0712SBarry Smith 
11013425bc38SStefano Zampini #undef __FUNCT__
11023425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
11033425bc38SStefano Zampini /*@
110428509bceSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
11053425bc38SStefano Zampini 
11063425bc38SStefano Zampini    Collective
11073425bc38SStefano Zampini 
11083425bc38SStefano Zampini    Input Parameters:
110928509bceSStefano Zampini +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
111028509bceSStefano Zampini .  standard_rhs - the right-hand side for your linear system
11113425bc38SStefano Zampini 
11123425bc38SStefano Zampini    Output Parameters:
111328509bceSStefano Zampini -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
11143425bc38SStefano Zampini 
11153425bc38SStefano Zampini    Level: developer
11163425bc38SStefano Zampini 
11173425bc38SStefano Zampini    Notes:
11183425bc38SStefano Zampini 
111928509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
11203425bc38SStefano Zampini @*/
11213425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
11223425bc38SStefano Zampini {
1123674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
11243425bc38SStefano Zampini   PetscErrorCode ierr;
11253425bc38SStefano Zampini 
11263425bc38SStefano Zampini   PetscFunctionBegin;
11273425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
11283425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
11293425bc38SStefano Zampini   PetscFunctionReturn(0);
11303425bc38SStefano Zampini }
11313425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
11321e6b0712SBarry Smith 
11333425bc38SStefano Zampini #undef __FUNCT__
11343425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
11353425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
11363425bc38SStefano Zampini {
1137674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
11383425bc38SStefano Zampini   PC_IS*         pcis;
11393425bc38SStefano Zampini   PC_BDDC*       pcbddc;
11403425bc38SStefano Zampini   PetscErrorCode ierr;
11413425bc38SStefano Zampini 
11423425bc38SStefano Zampini   PetscFunctionBegin;
11433425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
11443425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
11453425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
11463425bc38SStefano Zampini 
11473425bc38SStefano Zampini   /* apply B_delta^T */
11483425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11493425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11503425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
11513425bc38SStefano Zampini   /* compute rhs for BDDC application */
11523425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
11538eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
11543425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
11553425bc38SStefano Zampini   }
11563425bc38SStefano Zampini   /* apply BDDC */
11573425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
11583425bc38SStefano Zampini   /* put values into standard global vector */
11593425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11603425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11618eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
11623425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
11633425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
11643425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
11653425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
11663425bc38SStefano Zampini   }
11673425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11683425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11693425bc38SStefano Zampini   /* final change of basis if needed
11703425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
11713308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
11723425bc38SStefano Zampini   PetscFunctionReturn(0);
11733425bc38SStefano Zampini }
11741e6b0712SBarry Smith 
11753425bc38SStefano Zampini #undef __FUNCT__
11763425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
11773425bc38SStefano Zampini /*@
117828509bceSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
11793425bc38SStefano Zampini 
11803425bc38SStefano Zampini    Collective
11813425bc38SStefano Zampini 
11823425bc38SStefano Zampini    Input Parameters:
118328509bceSStefano Zampini +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
118428509bceSStefano Zampini .  fetidp_flux_sol - the solution of the FETIDP linear system
11853425bc38SStefano Zampini 
11863425bc38SStefano Zampini    Output Parameters:
118728509bceSStefano Zampini -  standard_sol      - the solution defined on the physical domain
11883425bc38SStefano Zampini 
11893425bc38SStefano Zampini    Level: developer
11903425bc38SStefano Zampini 
11913425bc38SStefano Zampini    Notes:
11923425bc38SStefano Zampini 
119328509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
11943425bc38SStefano Zampini @*/
11953425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
11963425bc38SStefano Zampini {
1197674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
11983425bc38SStefano Zampini   PetscErrorCode ierr;
11993425bc38SStefano Zampini 
12003425bc38SStefano Zampini   PetscFunctionBegin;
12013425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
12023425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
12033425bc38SStefano Zampini   PetscFunctionReturn(0);
12043425bc38SStefano Zampini }
12053425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
12061e6b0712SBarry Smith 
1207f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1208f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1209f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1210f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1211674ae819SStefano Zampini 
12123425bc38SStefano Zampini #undef __FUNCT__
12133425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
12143425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
12153425bc38SStefano Zampini {
1216674ae819SStefano Zampini 
1217674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
12183425bc38SStefano Zampini   Mat            newmat;
1219674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
12203425bc38SStefano Zampini   PC             newpc;
1221ce94432eSBarry Smith   MPI_Comm       comm;
12223425bc38SStefano Zampini   PetscErrorCode ierr;
12233425bc38SStefano Zampini 
12243425bc38SStefano Zampini   PetscFunctionBegin;
1225ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
12263425bc38SStefano Zampini   /* FETIDP linear matrix */
12273425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
12283425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
12293425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
12303425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
12313425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
12323425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
12333425bc38SStefano Zampini   /* FETIDP preconditioner */
12343425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
12353425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
12363425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
12373425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
12383425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
12393425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
12403425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
12413425bc38SStefano Zampini   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
12423425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
12433425bc38SStefano Zampini   /* return pointers for objects created */
12443425bc38SStefano Zampini   *fetidp_mat=newmat;
12453425bc38SStefano Zampini   *fetidp_pc=newpc;
12463425bc38SStefano Zampini   PetscFunctionReturn(0);
12473425bc38SStefano Zampini }
12481e6b0712SBarry Smith 
12493425bc38SStefano Zampini #undef __FUNCT__
12503425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
12513425bc38SStefano Zampini /*@
125228509bceSStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
12533425bc38SStefano Zampini 
12543425bc38SStefano Zampini    Collective
12553425bc38SStefano Zampini 
12563425bc38SStefano Zampini    Input Parameters:
125728509bceSStefano Zampini +  pc - the BDDC preconditioning context already setup
125828509bceSStefano Zampini 
125928509bceSStefano Zampini    Output Parameters:
126028509bceSStefano Zampini .  fetidp_mat - shell FETIDP matrix object
126128509bceSStefano Zampini .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
126228509bceSStefano Zampini 
126328509bceSStefano Zampini    Options Database Keys:
126428509bceSStefano Zampini -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
12653425bc38SStefano Zampini 
12663425bc38SStefano Zampini    Level: developer
12673425bc38SStefano Zampini 
12683425bc38SStefano Zampini    Notes:
126928509bceSStefano Zampini      Currently the only operation provided for FETIDP matrix is MatMult
12703425bc38SStefano Zampini 
12713425bc38SStefano Zampini .seealso: PCBDDC
12723425bc38SStefano Zampini @*/
12733425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
12743425bc38SStefano Zampini {
12753425bc38SStefano Zampini   PetscErrorCode ierr;
12763425bc38SStefano Zampini 
12773425bc38SStefano Zampini   PetscFunctionBegin;
12783425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12793425bc38SStefano Zampini   if (pc->setupcalled) {
1280516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1281f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
12823425bc38SStefano Zampini   PetscFunctionReturn(0);
12833425bc38SStefano Zampini }
12840c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1285da1bb401SStefano Zampini /*MC
1286da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
12870c7d97c5SJed Brown 
128828509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
128928509bceSStefano Zampini 
129028509bceSStefano Zampini .vb
129128509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
129228509bceSStefano 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
129328509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
129428509bceSStefano Zampini .ve
129528509bceSStefano Zampini 
129628509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
129728509bceSStefano Zampini 
1298b6fdb6dfSStefano Zampini    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
129928509bceSStefano Zampini 
130028509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
130128509bceSStefano Zampini 
1302b6fdb6dfSStefano 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.
1303b6fdb6dfSStefano Zampini 
130428509bceSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
130528509bceSStefano Zampini 
130628509bceSStefano 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
130728509bceSStefano Zampini 
130828509bceSStefano Zampini    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace.
130928509bceSStefano Zampini 
1310b6fdb6dfSStefano 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.
131128509bceSStefano Zampini 
131228509bceSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
131328509bceSStefano Zampini 
1314da1bb401SStefano Zampini    Options Database Keys:
131528509bceSStefano Zampini 
131628509bceSStefano Zampini .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
131728509bceSStefano Zampini .    -pc_bddc_use_edges <1> - use or not edges in primal space
1318b6fdb6dfSStefano Zampini .    -pc_bddc_use_faces <0> - use or not faces in primal space
131928509bceSStefano Zampini .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
132028509bceSStefano Zampini .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
132128509bceSStefano Zampini .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
132228509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
132328509bceSStefano Zampini .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
132428509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
132528509bceSStefano Zampini 
132628509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
132728509bceSStefano Zampini .vb
132828509bceSStefano Zampini       -pc_bddc_dirichlet_
132928509bceSStefano Zampini       -pc_bddc_neumann_
133028509bceSStefano Zampini       -pc_bddc_coarse_
133128509bceSStefano Zampini .ve
133228509bceSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
133328509bceSStefano Zampini 
133428509bceSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level can be specified as
133528509bceSStefano Zampini .vb
133628509bceSStefano Zampini       -pc_bddc_dirichlet_N_
133728509bceSStefano Zampini       -pc_bddc_neumann_N_
133828509bceSStefano Zampini       -pc_bddc_coarse_N_
133928509bceSStefano Zampini .ve
134028509bceSStefano Zampini    Note that level number ranges from the finest 0 to the coarsest N
1341da1bb401SStefano Zampini 
1342da1bb401SStefano Zampini    Level: intermediate
1343da1bb401SStefano Zampini 
1344b6fdb6dfSStefano Zampini    Developer notes:
134528509bceSStefano Zampini      Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose
1346da1bb401SStefano Zampini 
134728509bceSStefano Zampini      New deluxe scaling operator will be available soon.
1348da1bb401SStefano Zampini 
1349da1bb401SStefano Zampini    Contributed by Stefano Zampini
1350da1bb401SStefano Zampini 
1351da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1352da1bb401SStefano Zampini M*/
1353b2573a8aSBarry Smith 
1354da1bb401SStefano Zampini #undef __FUNCT__
1355da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
13568cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1357da1bb401SStefano Zampini {
1358da1bb401SStefano Zampini   PetscErrorCode      ierr;
1359da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1360da1bb401SStefano Zampini 
1361da1bb401SStefano Zampini   PetscFunctionBegin;
1362da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1363da1bb401SStefano Zampini   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1364da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1365da1bb401SStefano Zampini 
1366da1bb401SStefano Zampini   /* create PCIS data structure */
1367da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1368da1bb401SStefano Zampini 
136947d04d0dSStefano Zampini   /* BDDC customization */
137047d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
137147d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
137247d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
137347d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
137447d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
137547d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
137647d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
137747d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
137847d04d0dSStefano Zampini 
137939e2fb2aSStefano Zampini   pcbddc->BtoNmap                    = 0;
1380727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
1381e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
1382e9189074SStefano Zampini   pcbddc->n_actual_vertices          = 0;
1383e9189074SStefano Zampini   pcbddc->n_constraints              = 0;
1384727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
1385fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
1386f4ddd8eeSStefano Zampini   pcbddc->coarse_size                = 0;
1387f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
1388727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
1389f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
1390f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
1391f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
1392674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
13930bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
13943972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1395534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1396534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1397534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1398da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1399da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1400da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1401da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1402da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
140315aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
140415aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1405da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1406da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1407da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1408da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1409da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1410da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1411da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1412da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1413da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1414da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
141560d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
141660d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
1417da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
1418da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
141985c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
142047d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
142147d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
14224fad6a16SStefano Zampini   pcbddc->current_level              = 0;
14232b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1424da1bb401SStefano Zampini 
1425674ae819SStefano Zampini   /* create local graph structure */
1426674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1427674ae819SStefano Zampini 
1428674ae819SStefano Zampini   /* scaling */
1429674ae819SStefano Zampini   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1430674ae819SStefano Zampini   pcbddc->work_scaling               = 0;
1431da1bb401SStefano Zampini 
1432da1bb401SStefano Zampini   /* function pointers */
1433da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
1434da1bb401SStefano Zampini   pc->ops->applytranspose      = 0;
1435da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1436da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1437da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1438da1bb401SStefano Zampini   pc->ops->view                = 0;
1439da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1440da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1441da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1442534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1443534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1444da1bb401SStefano Zampini 
1445da1bb401SStefano Zampini   /* composing function */
1446674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1447bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
14482b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1449b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
14502b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1451bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1452bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1453bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1454bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1455bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1456bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1457bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1458bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1459bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1460bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1461da1bb401SStefano Zampini   PetscFunctionReturn(0);
1462da1bb401SStefano Zampini }
14633425bc38SStefano Zampini 
1464