xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 70cf5478945a9d09b9da66ef3d9730962f546e98)
153cdbc3dSStefano Zampini /* TODOLIST
2eb97c9d2SStefano Zampini 
3eb97c9d2SStefano Zampini    ConstraintsSetup
4eb97c9d2SStefano Zampini    - tolerances for constraints as an option (take care of single precision!)
5c1c8e736SStefano 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    - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access)
14eb97c9d2SStefano Zampini    - Provide PCApplyTranpose_BDDC
15eb97c9d2SStefano Zampini    - DofSplitting and DM attached to pc?
16eb97c9d2SStefano Zampini 
17eb97c9d2SStefano Zampini    Debugging output
18eb97c9d2SStefano Zampini    - Better management of verbosity levels of debugging output
19eb97c9d2SStefano Zampini 
20eb97c9d2SStefano Zampini    Build
21eb97c9d2SStefano Zampini    - make runexe59
22eb97c9d2SStefano Zampini 
23eb97c9d2SStefano Zampini    Extra
24eb97c9d2SStefano Zampini    - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
25eb97c9d2SStefano Zampini    - Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver?
26eb97c9d2SStefano Zampini    - add support for computing h,H and related using coordinates?
27c0b83709SStefano Zampini    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
28eb97c9d2SStefano Zampini    - Better management in PCIS code
29eb97c9d2SStefano Zampini    - BDDC with MG framework?
30eb97c9d2SStefano Zampini 
31eb97c9d2SStefano Zampini    FETIDP
32eb97c9d2SStefano Zampini    - Move FETIDP code to its own classes
33eb97c9d2SStefano Zampini 
34eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
35eb97c9d2SStefano Zampini    - Provide general case for subassembling
3639e2fb2aSStefano Zampini    - Preallocation routines in MatISGetMPIAXAIJ
37eb97c9d2SStefano Zampini 
3853cdbc3dSStefano Zampini */
390c7d97c5SJed Brown 
4053cdbc3dSStefano Zampini /* ----------------------------------------------------------------------------------------------------------------------------------------------
410c7d97c5SJed Brown    Implementation of BDDC preconditioner based on:
420c7d97c5SJed Brown    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
4353cdbc3dSStefano Zampini    ---------------------------------------------------------------------------------------------------------------------------------------------- */
4453cdbc3dSStefano Zampini 
45674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/  /* includes for fortran wrappers */
46674ae819SStefano Zampini #include "bddcprivate.h"
473b03a366Sstefano_zampini #include <petscblaslapack.h>
48674ae819SStefano Zampini 
490c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
500c7d97c5SJed Brown #undef __FUNCT__
510c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
520c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc)
530c7d97c5SJed Brown {
540c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
550c7d97c5SJed Brown   PetscErrorCode ierr;
560c7d97c5SJed Brown 
570c7d97c5SJed Brown   PetscFunctionBegin;
580c7d97c5SJed Brown   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
598eeda7d8SStefano Zampini   /* Verbose debugging */
608eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
618eeda7d8SStefano Zampini   /* Primal space cumstomization */
628eeda7d8SStefano 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);
638eeda7d8SStefano 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);
648eeda7d8SStefano 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);
658eeda7d8SStefano Zampini   /* Change of basis */
668eeda7d8SStefano 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);
678eeda7d8SStefano 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);
68674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
69674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
70674ae819SStefano Zampini   }
718eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
728eeda7d8SStefano 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);
730298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
742b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
75674ae819SStefano 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);
760c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
770c7d97c5SJed Brown   PetscFunctionReturn(0);
780c7d97c5SJed Brown }
790c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
80674ae819SStefano Zampini #undef __FUNCT__
81674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
82674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
83674ae819SStefano Zampini {
84674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
85674ae819SStefano Zampini   PetscErrorCode ierr;
861e6b0712SBarry Smith 
87674ae819SStefano Zampini   PetscFunctionBegin;
88674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
89674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
90674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
91674ae819SStefano Zampini   PetscFunctionReturn(0);
92674ae819SStefano Zampini }
93674ae819SStefano Zampini #undef __FUNCT__
94674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
95674ae819SStefano Zampini /*@
9628509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
97674ae819SStefano Zampini 
98674ae819SStefano Zampini    Not collective
99674ae819SStefano Zampini 
100674ae819SStefano Zampini    Input Parameters:
101674ae819SStefano Zampini +  pc - the preconditioning context
10228509bceSStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering
103674ae819SStefano Zampini 
104674ae819SStefano Zampini    Level: intermediate
105674ae819SStefano Zampini 
106674ae819SStefano Zampini    Notes:
107674ae819SStefano Zampini 
108674ae819SStefano Zampini .seealso: PCBDDC
109674ae819SStefano Zampini @*/
110674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
111674ae819SStefano Zampini {
112674ae819SStefano Zampini   PetscErrorCode ierr;
113674ae819SStefano Zampini 
114674ae819SStefano Zampini   PetscFunctionBegin;
115674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
116674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
117674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
118674ae819SStefano Zampini   PetscFunctionReturn(0);
119674ae819SStefano Zampini }
120674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1210c7d97c5SJed Brown #undef __FUNCT__
1224fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1234fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1244fad6a16SStefano Zampini {
1254fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1264fad6a16SStefano Zampini 
1274fad6a16SStefano Zampini   PetscFunctionBegin;
1284fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1294fad6a16SStefano Zampini   PetscFunctionReturn(0);
1304fad6a16SStefano Zampini }
1311e6b0712SBarry Smith 
1324fad6a16SStefano Zampini #undef __FUNCT__
1334fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1344fad6a16SStefano Zampini /*@
13528509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
1364fad6a16SStefano Zampini 
1374fad6a16SStefano Zampini    Logically collective on PC
1384fad6a16SStefano Zampini 
1394fad6a16SStefano Zampini    Input Parameters:
1404fad6a16SStefano Zampini +  pc - the preconditioning context
14128509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
1424fad6a16SStefano Zampini 
14328509bceSStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
1444fad6a16SStefano Zampini 
1454fad6a16SStefano Zampini    Level: intermediate
1464fad6a16SStefano Zampini 
1474fad6a16SStefano Zampini    Notes:
1484fad6a16SStefano Zampini 
1494fad6a16SStefano Zampini .seealso: PCBDDC
1504fad6a16SStefano Zampini @*/
1514fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
1524fad6a16SStefano Zampini {
1534fad6a16SStefano Zampini   PetscErrorCode ierr;
1544fad6a16SStefano Zampini 
1554fad6a16SStefano Zampini   PetscFunctionBegin;
1564fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1572b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
1584fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
1594fad6a16SStefano Zampini   PetscFunctionReturn(0);
1604fad6a16SStefano Zampini }
1612b510759SStefano Zampini 
162b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
1632b510759SStefano Zampini #undef __FUNCT__
164b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
165b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
166b8ffe317SStefano Zampini {
167b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
168b8ffe317SStefano Zampini 
169b8ffe317SStefano Zampini   PetscFunctionBegin;
17085c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
171b8ffe317SStefano Zampini   PetscFunctionReturn(0);
172b8ffe317SStefano Zampini }
173b8ffe317SStefano Zampini 
174b8ffe317SStefano Zampini #undef __FUNCT__
175b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
176b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
1772b510759SStefano Zampini {
1782b510759SStefano Zampini   PetscErrorCode ierr;
1792b510759SStefano Zampini 
1802b510759SStefano Zampini   PetscFunctionBegin;
1812b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
182b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
183b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1842b510759SStefano Zampini   PetscFunctionReturn(0);
1852b510759SStefano Zampini }
1861e6b0712SBarry Smith 
1874fad6a16SStefano Zampini #undef __FUNCT__
1882b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
1892b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
1904fad6a16SStefano Zampini {
1914fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1924fad6a16SStefano Zampini 
1934fad6a16SStefano Zampini   PetscFunctionBegin;
1942b510759SStefano Zampini   pcbddc->current_level = level;
1954fad6a16SStefano Zampini   PetscFunctionReturn(0);
1964fad6a16SStefano Zampini }
1971e6b0712SBarry Smith 
1984fad6a16SStefano Zampini #undef __FUNCT__
199b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
200b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
201b8ffe317SStefano Zampini {
202b8ffe317SStefano Zampini   PetscErrorCode ierr;
203b8ffe317SStefano Zampini 
204b8ffe317SStefano Zampini   PetscFunctionBegin;
205b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
206b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
207b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
208b8ffe317SStefano Zampini   PetscFunctionReturn(0);
209b8ffe317SStefano Zampini }
210b8ffe317SStefano Zampini 
211b8ffe317SStefano Zampini #undef __FUNCT__
2122b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
2132b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
2142b510759SStefano Zampini {
2152b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2162b510759SStefano Zampini 
2172b510759SStefano Zampini   PetscFunctionBegin;
2182b510759SStefano Zampini   pcbddc->max_levels = levels;
2192b510759SStefano Zampini   PetscFunctionReturn(0);
2202b510759SStefano Zampini }
2212b510759SStefano Zampini 
2222b510759SStefano Zampini #undef __FUNCT__
2232b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2244fad6a16SStefano Zampini /*@
22528509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
2264fad6a16SStefano Zampini 
2274fad6a16SStefano Zampini    Logically collective on PC
2284fad6a16SStefano Zampini 
2294fad6a16SStefano Zampini    Input Parameters:
2304fad6a16SStefano Zampini +  pc - the preconditioning context
23128509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
2324fad6a16SStefano Zampini 
23328509bceSStefano Zampini    Default value is 0, i.e. traditional one-level BDDC
2344fad6a16SStefano Zampini 
2354fad6a16SStefano Zampini    Level: intermediate
2364fad6a16SStefano Zampini 
2374fad6a16SStefano Zampini    Notes:
2384fad6a16SStefano Zampini 
2394fad6a16SStefano Zampini .seealso: PCBDDC
2404fad6a16SStefano Zampini @*/
2412b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
2424fad6a16SStefano Zampini {
2434fad6a16SStefano Zampini   PetscErrorCode ierr;
2444fad6a16SStefano Zampini 
2454fad6a16SStefano Zampini   PetscFunctionBegin;
2464fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2472b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
2482b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
2494fad6a16SStefano Zampini   PetscFunctionReturn(0);
2504fad6a16SStefano Zampini }
2514fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
2521e6b0712SBarry Smith 
2534fad6a16SStefano Zampini #undef __FUNCT__
2540bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
2550bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
2560bdf917eSStefano Zampini {
2570bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2580bdf917eSStefano Zampini   PetscErrorCode ierr;
2590bdf917eSStefano Zampini 
2600bdf917eSStefano Zampini   PetscFunctionBegin;
2610bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
2620bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
2630bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
2640bdf917eSStefano Zampini   PetscFunctionReturn(0);
2650bdf917eSStefano Zampini }
2661e6b0712SBarry Smith 
2670bdf917eSStefano Zampini #undef __FUNCT__
2680bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
2690bdf917eSStefano Zampini /*@
27028509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
2710bdf917eSStefano Zampini 
2720bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
2730bdf917eSStefano Zampini 
2740bdf917eSStefano Zampini    Input Parameters:
2750bdf917eSStefano Zampini +  pc - the preconditioning context
27628509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
2770bdf917eSStefano Zampini 
2780bdf917eSStefano Zampini    Level: intermediate
2790bdf917eSStefano Zampini 
2800bdf917eSStefano Zampini    Notes:
2810bdf917eSStefano Zampini 
2820bdf917eSStefano Zampini .seealso: PCBDDC
2830bdf917eSStefano Zampini @*/
2840bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
2850bdf917eSStefano Zampini {
2860bdf917eSStefano Zampini   PetscErrorCode ierr;
2870bdf917eSStefano Zampini 
2880bdf917eSStefano Zampini   PetscFunctionBegin;
2890bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
290674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
2912b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
2920bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
2930bdf917eSStefano Zampini   PetscFunctionReturn(0);
2940bdf917eSStefano Zampini }
2950bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
2961e6b0712SBarry Smith 
2970bdf917eSStefano Zampini #undef __FUNCT__
2983b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
2993b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
3003b03a366Sstefano_zampini {
3013b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3023b03a366Sstefano_zampini   PetscErrorCode ierr;
3033b03a366Sstefano_zampini 
3043b03a366Sstefano_zampini   PetscFunctionBegin;
305785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
306785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3073b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
30836e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
30936e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
310fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3113b03a366Sstefano_zampini   PetscFunctionReturn(0);
3123b03a366Sstefano_zampini }
3131e6b0712SBarry Smith 
3143b03a366Sstefano_zampini #undef __FUNCT__
3153b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
3163b03a366Sstefano_zampini /*@
31728509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
3183b03a366Sstefano_zampini 
319785d1243SStefano Zampini    Collective
3203b03a366Sstefano_zampini 
3213b03a366Sstefano_zampini    Input Parameters:
3223b03a366Sstefano_zampini +  pc - the preconditioning context
323785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
3243b03a366Sstefano_zampini 
3253b03a366Sstefano_zampini    Level: intermediate
3263b03a366Sstefano_zampini 
327785d1243SStefano Zampini    Notes: Any process can list any global node
3283b03a366Sstefano_zampini 
3293b03a366Sstefano_zampini .seealso: PCBDDC
3303b03a366Sstefano_zampini @*/
3313b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3323b03a366Sstefano_zampini {
3333b03a366Sstefano_zampini   PetscErrorCode ierr;
3343b03a366Sstefano_zampini 
3353b03a366Sstefano_zampini   PetscFunctionBegin;
3363b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
337674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
338785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
3393b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
3403b03a366Sstefano_zampini   PetscFunctionReturn(0);
3413b03a366Sstefano_zampini }
3423b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
3431e6b0712SBarry Smith 
3443b03a366Sstefano_zampini #undef __FUNCT__
34582ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
34682ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
3470c7d97c5SJed Brown {
3480c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3490c7d97c5SJed Brown   PetscErrorCode ierr;
3500c7d97c5SJed Brown 
3510c7d97c5SJed Brown   PetscFunctionBegin;
352785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
353785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3540c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
3550c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
356785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
3577d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3580c7d97c5SJed Brown   PetscFunctionReturn(0);
3590c7d97c5SJed Brown }
3600c7d97c5SJed Brown 
3610c7d97c5SJed Brown #undef __FUNCT__
36282ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
3639c0446d6SStefano Zampini /*@
36482ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
3659c0446d6SStefano Zampini 
366785d1243SStefano Zampini    Collective
36753cdbc3dSStefano Zampini 
36853cdbc3dSStefano Zampini    Input Parameters:
36953cdbc3dSStefano Zampini +  pc - the preconditioning context
37082ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
37153cdbc3dSStefano Zampini 
37253cdbc3dSStefano Zampini    Level: intermediate
37353cdbc3dSStefano Zampini 
3749c0446d6SStefano Zampini    Notes:
37553cdbc3dSStefano Zampini 
37653cdbc3dSStefano Zampini .seealso: PCBDDC
37753cdbc3dSStefano Zampini @*/
37882ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
3790c7d97c5SJed Brown {
3800c7d97c5SJed Brown   PetscErrorCode ierr;
3810c7d97c5SJed Brown 
3820c7d97c5SJed Brown   PetscFunctionBegin;
3830c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3840c7d97c5SJed Brown   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
38582ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
38682ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
3870c7d97c5SJed Brown   PetscFunctionReturn(0);
3880c7d97c5SJed Brown }
3890c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
3900c7d97c5SJed Brown 
3910c7d97c5SJed Brown #undef __FUNCT__
3920c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
39353cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
3940c7d97c5SJed Brown {
3950c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
39653cdbc3dSStefano Zampini   PetscErrorCode ierr;
3970c7d97c5SJed Brown 
3980c7d97c5SJed Brown   PetscFunctionBegin;
399785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
400785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
40153cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
40236e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
40336e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
404fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4050c7d97c5SJed Brown   PetscFunctionReturn(0);
4060c7d97c5SJed Brown }
4071e6b0712SBarry Smith 
4080c7d97c5SJed Brown #undef __FUNCT__
4090c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
41057527edcSJed Brown /*@
41128509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
41257527edcSJed Brown 
413785d1243SStefano Zampini    Collective
41457527edcSJed Brown 
41557527edcSJed Brown    Input Parameters:
41657527edcSJed Brown +  pc - the preconditioning context
417785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
41857527edcSJed Brown 
41957527edcSJed Brown    Level: intermediate
42057527edcSJed Brown 
421785d1243SStefano Zampini    Notes: Any process can list any global node
42257527edcSJed Brown 
42357527edcSJed Brown .seealso: PCBDDC
42457527edcSJed Brown @*/
42553cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
4260c7d97c5SJed Brown {
4270c7d97c5SJed Brown   PetscErrorCode ierr;
4280c7d97c5SJed Brown 
4290c7d97c5SJed Brown   PetscFunctionBegin;
4300c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
431674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
432785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
43353cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
43453cdbc3dSStefano Zampini   PetscFunctionReturn(0);
43553cdbc3dSStefano Zampini }
43653cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
4371e6b0712SBarry Smith 
43853cdbc3dSStefano Zampini #undef __FUNCT__
43982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
44082ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
4410c7d97c5SJed Brown {
4420c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4430c7d97c5SJed Brown   PetscErrorCode ierr;
4440c7d97c5SJed Brown 
4450c7d97c5SJed Brown   PetscFunctionBegin;
446785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
447785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
4480c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
4490c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
450785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
4517d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4520c7d97c5SJed Brown   PetscFunctionReturn(0);
4530c7d97c5SJed Brown }
4540c7d97c5SJed Brown 
4550c7d97c5SJed Brown #undef __FUNCT__
45682ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
4570c7d97c5SJed Brown /*@
45882ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
4590c7d97c5SJed Brown 
460785d1243SStefano Zampini    Collective
4610c7d97c5SJed Brown 
4620c7d97c5SJed Brown    Input Parameters:
4630c7d97c5SJed Brown +  pc - the preconditioning context
46482ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
4650c7d97c5SJed Brown 
4660c7d97c5SJed Brown    Level: intermediate
4670c7d97c5SJed Brown 
4680c7d97c5SJed Brown    Notes:
4690c7d97c5SJed Brown 
4700c7d97c5SJed Brown .seealso: PCBDDC
4710c7d97c5SJed Brown @*/
47282ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
4730c7d97c5SJed Brown {
4740c7d97c5SJed Brown   PetscErrorCode ierr;
4750c7d97c5SJed Brown 
4760c7d97c5SJed Brown   PetscFunctionBegin;
4770c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4780c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
47982ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
48082ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
48153cdbc3dSStefano Zampini   PetscFunctionReturn(0);
48253cdbc3dSStefano Zampini }
48353cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
48453cdbc3dSStefano Zampini 
48553cdbc3dSStefano Zampini #undef __FUNCT__
486da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
487da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
488da1bb401SStefano Zampini {
489da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
490da1bb401SStefano Zampini 
491da1bb401SStefano Zampini   PetscFunctionBegin;
492da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
493da1bb401SStefano Zampini   PetscFunctionReturn(0);
494da1bb401SStefano Zampini }
4951e6b0712SBarry Smith 
496da1bb401SStefano Zampini #undef __FUNCT__
497da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
498da1bb401SStefano Zampini /*@
499785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
500da1bb401SStefano Zampini 
501785d1243SStefano Zampini    Collective
502785d1243SStefano Zampini 
503785d1243SStefano Zampini    Input Parameters:
504785d1243SStefano Zampini .  pc - the preconditioning context
505785d1243SStefano Zampini 
506785d1243SStefano Zampini    Output Parameters:
507785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
508785d1243SStefano Zampini 
509785d1243SStefano Zampini    Level: intermediate
510785d1243SStefano Zampini 
511785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
512785d1243SStefano Zampini 
513785d1243SStefano Zampini .seealso: PCBDDC
514785d1243SStefano Zampini @*/
515785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
516785d1243SStefano Zampini {
517785d1243SStefano Zampini   PetscErrorCode ierr;
518785d1243SStefano Zampini 
519785d1243SStefano Zampini   PetscFunctionBegin;
520785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
521785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
522785d1243SStefano Zampini   PetscFunctionReturn(0);
523785d1243SStefano Zampini }
524785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
525785d1243SStefano Zampini 
526785d1243SStefano Zampini #undef __FUNCT__
527785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
528785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
529785d1243SStefano Zampini {
530785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
531785d1243SStefano Zampini 
532785d1243SStefano Zampini   PetscFunctionBegin;
533785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
534785d1243SStefano Zampini   PetscFunctionReturn(0);
535785d1243SStefano Zampini }
536785d1243SStefano Zampini 
537785d1243SStefano Zampini #undef __FUNCT__
53882ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
539da1bb401SStefano Zampini /*@
54082ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
541da1bb401SStefano Zampini 
542785d1243SStefano Zampini    Collective
543da1bb401SStefano Zampini 
544da1bb401SStefano Zampini    Input Parameters:
54528509bceSStefano Zampini .  pc - the preconditioning context
546da1bb401SStefano Zampini 
547da1bb401SStefano Zampini    Output Parameters:
54828509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
549da1bb401SStefano Zampini 
550da1bb401SStefano Zampini    Level: intermediate
551da1bb401SStefano Zampini 
552da1bb401SStefano Zampini    Notes:
553da1bb401SStefano Zampini 
554da1bb401SStefano Zampini .seealso: PCBDDC
555da1bb401SStefano Zampini @*/
55682ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
557da1bb401SStefano Zampini {
558da1bb401SStefano Zampini   PetscErrorCode ierr;
559da1bb401SStefano Zampini 
560da1bb401SStefano Zampini   PetscFunctionBegin;
561da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
56282ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
563da1bb401SStefano Zampini   PetscFunctionReturn(0);
564da1bb401SStefano Zampini }
565da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
5661e6b0712SBarry Smith 
567da1bb401SStefano Zampini #undef __FUNCT__
56853cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
56953cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
57053cdbc3dSStefano Zampini {
57153cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
57253cdbc3dSStefano Zampini 
57353cdbc3dSStefano Zampini   PetscFunctionBegin;
57453cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
57553cdbc3dSStefano Zampini   PetscFunctionReturn(0);
57653cdbc3dSStefano Zampini }
5771e6b0712SBarry Smith 
57853cdbc3dSStefano Zampini #undef __FUNCT__
57953cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
58053cdbc3dSStefano Zampini /*@
581785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
58253cdbc3dSStefano Zampini 
583785d1243SStefano Zampini    Collective
584785d1243SStefano Zampini 
585785d1243SStefano Zampini    Input Parameters:
586785d1243SStefano Zampini .  pc - the preconditioning context
587785d1243SStefano Zampini 
588785d1243SStefano Zampini    Output Parameters:
589785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
590785d1243SStefano Zampini 
591785d1243SStefano Zampini    Level: intermediate
592785d1243SStefano Zampini 
593785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
594785d1243SStefano Zampini 
595785d1243SStefano Zampini .seealso: PCBDDC
596785d1243SStefano Zampini @*/
597785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
598785d1243SStefano Zampini {
599785d1243SStefano Zampini   PetscErrorCode ierr;
600785d1243SStefano Zampini 
601785d1243SStefano Zampini   PetscFunctionBegin;
602785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
603785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
604785d1243SStefano Zampini   PetscFunctionReturn(0);
605785d1243SStefano Zampini }
606785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
607785d1243SStefano Zampini 
608785d1243SStefano Zampini #undef __FUNCT__
609785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
610785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
611785d1243SStefano Zampini {
612785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
613785d1243SStefano Zampini 
614785d1243SStefano Zampini   PetscFunctionBegin;
615785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
616785d1243SStefano Zampini   PetscFunctionReturn(0);
617785d1243SStefano Zampini }
618785d1243SStefano Zampini 
619785d1243SStefano Zampini #undef __FUNCT__
62082ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
62153cdbc3dSStefano Zampini /*@
62282ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
62353cdbc3dSStefano Zampini 
624785d1243SStefano Zampini    Collective
62553cdbc3dSStefano Zampini 
62653cdbc3dSStefano Zampini    Input Parameters:
62728509bceSStefano Zampini .  pc - the preconditioning context
62853cdbc3dSStefano Zampini 
62953cdbc3dSStefano Zampini    Output Parameters:
63028509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
63153cdbc3dSStefano Zampini 
63253cdbc3dSStefano Zampini    Level: intermediate
63353cdbc3dSStefano Zampini 
63453cdbc3dSStefano Zampini    Notes:
63553cdbc3dSStefano Zampini 
63653cdbc3dSStefano Zampini .seealso: PCBDDC
63753cdbc3dSStefano Zampini @*/
63882ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
63953cdbc3dSStefano Zampini {
64053cdbc3dSStefano Zampini   PetscErrorCode ierr;
64153cdbc3dSStefano Zampini 
64253cdbc3dSStefano Zampini   PetscFunctionBegin;
64353cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
64482ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
6450c7d97c5SJed Brown   PetscFunctionReturn(0);
6460c7d97c5SJed Brown }
64736e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
6481e6b0712SBarry Smith 
64936e030ebSStefano Zampini #undef __FUNCT__
650da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
6511a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
65236e030ebSStefano Zampini {
65336e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
654da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
655da1bb401SStefano Zampini   PetscErrorCode ierr;
65636e030ebSStefano Zampini 
65736e030ebSStefano Zampini   PetscFunctionBegin;
658674ae819SStefano Zampini   /* free old CSR */
659674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
660fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
661674ae819SStefano Zampini   /* get CSR into graph structure */
662da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
663674ae819SStefano Zampini     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
664674ae819SStefano Zampini     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
665674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
666674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
667da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
6681a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
6691a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
670674ae819SStefano Zampini   }
671575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
67236e030ebSStefano Zampini   PetscFunctionReturn(0);
67336e030ebSStefano Zampini }
6741e6b0712SBarry Smith 
67536e030ebSStefano Zampini #undef __FUNCT__
676da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
67736e030ebSStefano Zampini /*@
67828509bceSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
67936e030ebSStefano Zampini 
68036e030ebSStefano Zampini    Not collective
68136e030ebSStefano Zampini 
68236e030ebSStefano Zampini    Input Parameters:
68336e030ebSStefano Zampini +  pc - the preconditioning context
68428509bceSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
68528509bceSStefano Zampini .  xadj, adjncy - the CSR graph
68628509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
68736e030ebSStefano Zampini 
68836e030ebSStefano Zampini    Level: intermediate
68936e030ebSStefano Zampini 
69036e030ebSStefano Zampini    Notes:
69136e030ebSStefano Zampini 
69228509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
69336e030ebSStefano Zampini @*/
6941a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
69536e030ebSStefano Zampini {
696575ad6abSStefano Zampini   void (*f)(void) = 0;
69736e030ebSStefano Zampini   PetscErrorCode ierr;
69836e030ebSStefano Zampini 
69936e030ebSStefano Zampini   PetscFunctionBegin;
70036e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
701674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
702674ae819SStefano Zampini   PetscValidIntPointer(xadj,4);
703674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
704674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
705da1bb401SStefano Zampini   }
70636e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
707575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
708575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
709575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
710575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
711575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
71236e030ebSStefano Zampini   }
71336e030ebSStefano Zampini   PetscFunctionReturn(0);
71436e030ebSStefano Zampini }
7159c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
7161e6b0712SBarry Smith 
7179c0446d6SStefano Zampini #undef __FUNCT__
71863602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
71963602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
72063602bcaSStefano Zampini {
72163602bcaSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
72263602bcaSStefano Zampini   PetscInt i;
72363602bcaSStefano Zampini   PetscErrorCode ierr;
72463602bcaSStefano Zampini 
72563602bcaSStefano Zampini   PetscFunctionBegin;
72663602bcaSStefano Zampini   /* Destroy ISes if they were already set */
72763602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
72863602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
72963602bcaSStefano Zampini   }
73063602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
73163602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
73263602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
73363602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
73463602bcaSStefano Zampini   }
73563602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
73663602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
73763602bcaSStefano Zampini   /* allocate space then set */
73863602bcaSStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
73963602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
74063602bcaSStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
74163602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
74263602bcaSStefano Zampini   }
74363602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal=n_is;
74463602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
7451a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
74663602bcaSStefano Zampini   PetscFunctionReturn(0);
74763602bcaSStefano Zampini }
74863602bcaSStefano Zampini 
74963602bcaSStefano Zampini #undef __FUNCT__
75063602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
75163602bcaSStefano Zampini /*@
75263602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
75363602bcaSStefano Zampini 
75463602bcaSStefano Zampini    Collective
75563602bcaSStefano Zampini 
75663602bcaSStefano Zampini    Input Parameters:
75763602bcaSStefano Zampini +  pc - the preconditioning context
75863602bcaSStefano Zampini -  n_is - number of index sets defining the fields
75963602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in local ordering
76063602bcaSStefano Zampini 
76163602bcaSStefano Zampini    Level: intermediate
76263602bcaSStefano Zampini 
76363602bcaSStefano Zampini    Notes: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to a different field.
76463602bcaSStefano Zampini 
76563602bcaSStefano Zampini .seealso: PCBDDC
76663602bcaSStefano Zampini @*/
76763602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
76863602bcaSStefano Zampini {
76963602bcaSStefano Zampini   PetscInt       i;
77063602bcaSStefano Zampini   PetscErrorCode ierr;
77163602bcaSStefano Zampini 
77263602bcaSStefano Zampini   PetscFunctionBegin;
77363602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
77463602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
77563602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
77663602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
77763602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
77863602bcaSStefano Zampini   }
77963602bcaSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
78063602bcaSStefano Zampini   PetscFunctionReturn(0);
78163602bcaSStefano Zampini }
78263602bcaSStefano Zampini /* -------------------------------------------------------------------------- */
78363602bcaSStefano Zampini 
78463602bcaSStefano Zampini #undef __FUNCT__
7859c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
7869c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
7879c0446d6SStefano Zampini {
7889c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
7899c0446d6SStefano Zampini   PetscInt i;
7909c0446d6SStefano Zampini   PetscErrorCode ierr;
7919c0446d6SStefano Zampini 
7929c0446d6SStefano Zampini   PetscFunctionBegin;
793da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
7949c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
7959c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
7969c0446d6SStefano Zampini   }
797d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
79863602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
79963602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
80063602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
80163602bcaSStefano Zampini   }
80263602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
80363602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
804da1bb401SStefano Zampini   /* allocate space then set */
8059c0446d6SStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
8069c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
807da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
808da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
8099c0446d6SStefano Zampini   }
8109c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
81163602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8121a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
8139c0446d6SStefano Zampini   PetscFunctionReturn(0);
8149c0446d6SStefano Zampini }
8151e6b0712SBarry Smith 
8169c0446d6SStefano Zampini #undef __FUNCT__
8179c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
8189c0446d6SStefano Zampini /*@
81963602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
8209c0446d6SStefano Zampini 
82163602bcaSStefano Zampini    Collective
8229c0446d6SStefano Zampini 
8239c0446d6SStefano Zampini    Input Parameters:
8249c0446d6SStefano Zampini +  pc - the preconditioning context
82528509bceSStefano Zampini -  n_is - number of index sets defining the fields
82663602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in global ordering
8279c0446d6SStefano Zampini 
8289c0446d6SStefano Zampini    Level: intermediate
8299c0446d6SStefano Zampini 
83063602bcaSStefano Zampini    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.
8319c0446d6SStefano Zampini 
8329c0446d6SStefano Zampini .seealso: PCBDDC
8339c0446d6SStefano Zampini @*/
8349c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
8359c0446d6SStefano Zampini {
8362b510759SStefano Zampini   PetscInt       i;
8379c0446d6SStefano Zampini   PetscErrorCode ierr;
8389c0446d6SStefano Zampini 
8399c0446d6SStefano Zampini   PetscFunctionBegin;
8409c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
84163602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
8422b510759SStefano Zampini   for (i=0;i<n_is;i++) {
84363602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
84463602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
8452b510759SStefano Zampini   }
8469c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
8479c0446d6SStefano Zampini   PetscFunctionReturn(0);
8489c0446d6SStefano Zampini }
849da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
850534831adSStefano Zampini #undef __FUNCT__
851534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
852534831adSStefano Zampini /* -------------------------------------------------------------------------- */
853534831adSStefano Zampini /*
854534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
855534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
8569c0446d6SStefano Zampini 
857534831adSStefano Zampini    Input Parameter:
858534831adSStefano Zampini +  pc - the preconditioner contex
859534831adSStefano Zampini 
860534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
861534831adSStefano Zampini 
862534831adSStefano Zampini    Notes:
863534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
864534831adSStefano Zampini    the user, but instead is called by KSPSolve().
865534831adSStefano Zampini */
866534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
867534831adSStefano Zampini {
868534831adSStefano Zampini   PetscErrorCode ierr;
869534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
870534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
871534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
872534831adSStefano Zampini   Mat            temp_mat;
8733972b0daSStefano Zampini   IS             dirIS;
8743972b0daSStefano Zampini   Vec            used_vec;
87582ba6b80SStefano Zampini   PetscBool      guess_nonzero;
876534831adSStefano Zampini 
877534831adSStefano Zampini   PetscFunctionBegin;
87885c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
87985c4d303SStefano Zampini   if (ksp) {
88085c4d303SStefano Zampini     PetscBool iscg;
88185c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
88285c4d303SStefano Zampini     if (!iscg) {
88385c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
88485c4d303SStefano Zampini     }
88585c4d303SStefano Zampini   }
88685c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
88762a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
88862a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
88962a6ff1dSStefano Zampini   }
89062a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
89162a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
89262a6ff1dSStefano Zampini   }
8933972b0daSStefano Zampini   if (x) {
8943972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
8953972b0daSStefano Zampini     used_vec = x;
8963972b0daSStefano Zampini   } else {
8973972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
8983972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
8993972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9003972b0daSStefano Zampini   }
9013972b0daSStefano Zampini   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
9023972b0daSStefano Zampini   if (ksp) {
9033972b0daSStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
9043972b0daSStefano Zampini     if (!guess_nonzero) {
9053972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9063972b0daSStefano Zampini     }
9073972b0daSStefano Zampini   }
9083308cffdSStefano Zampini 
9093972b0daSStefano Zampini   /* store the original rhs */
9103972b0daSStefano Zampini   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
9113972b0daSStefano Zampini 
9123972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
913785d1243SStefano Zampini   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
91482ba6b80SStefano Zampini   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
91582ba6b80SStefano Zampini   if (rhs && dirIS) {
916785d1243SStefano Zampini     PetscInt    dirsize,i,*is_indices;
917785d1243SStefano Zampini     PetscScalar *array_x,*array_diagonal;
918785d1243SStefano Zampini 
9193972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
9203972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
9213972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9223972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9233972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9243972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
92582ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
9263972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
9273972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
9283972b0daSStefano Zampini     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9292fa5cd67SKarl Rupp     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
9303972b0daSStefano Zampini     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9313972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
9323972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
9333972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9343972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
935b76ba322SStefano Zampini 
9363972b0daSStefano Zampini     /* remove the computed solution from the rhs */
9373972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
9383972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
9393972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
9403308cffdSStefano Zampini   }
941b76ba322SStefano Zampini 
942b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
9433972b0daSStefano Zampini   if (x) {
9443972b0daSStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
9453972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
94685c4d303SStefano Zampini     if (pcbddc->use_exact_dirichlet_trick) {
947b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
948b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
949b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
950b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
951b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
952b76ba322SStefano Zampini       if (ksp) {
953b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
954b76ba322SStefano Zampini       }
955b76ba322SStefano Zampini     }
9563972b0daSStefano Zampini   }
957b76ba322SStefano Zampini 
95892e3dcfbSStefano Zampini   /* prepare MatMult and rhs for solver */
959674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
960b76ba322SStefano Zampini     /* swap pointers for local matrices */
961b76ba322SStefano Zampini     temp_mat = matis->A;
962b76ba322SStefano Zampini     matis->A = pcbddc->local_mat;
963b76ba322SStefano Zampini     pcbddc->local_mat = temp_mat;
96492e3dcfbSStefano Zampini     if (rhs) {
965b76ba322SStefano Zampini       /* Get local rhs and apply transformation of basis */
966b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
967b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
968b76ba322SStefano Zampini       /* from original basis to modified basis */
969b76ba322SStefano Zampini       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
970b76ba322SStefano Zampini       /* put back modified values into the global vec using INSERT_VALUES copy mode */
971b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
972b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
973674ae819SStefano Zampini     }
97492e3dcfbSStefano Zampini   }
97592e3dcfbSStefano Zampini 
97692e3dcfbSStefano Zampini   /* remove nullspace if present */
9770bdf917eSStefano Zampini   if (ksp && pcbddc->NullSpace) {
978d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
979d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
980b76ba322SStefano Zampini   }
9810bdf917eSStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
982534831adSStefano Zampini   PetscFunctionReturn(0);
983534831adSStefano Zampini }
984534831adSStefano Zampini /* -------------------------------------------------------------------------- */
985534831adSStefano Zampini #undef __FUNCT__
986534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
987534831adSStefano Zampini /* -------------------------------------------------------------------------- */
988534831adSStefano Zampini /*
989534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
990534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
991534831adSStefano Zampini 
992534831adSStefano Zampini    Input Parameter:
993534831adSStefano Zampini +  pc - the preconditioner contex
994534831adSStefano Zampini 
995534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
996534831adSStefano Zampini 
997534831adSStefano Zampini    Notes:
998534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
999534831adSStefano Zampini    the user, but instead is called by KSPSolve().
1000534831adSStefano Zampini */
1001534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1002534831adSStefano Zampini {
1003534831adSStefano Zampini   PetscErrorCode ierr;
1004534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1005534831adSStefano Zampini   PC_IS          *pcis   = (PC_IS*)(pc->data);
1006534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
1007534831adSStefano Zampini   Mat            temp_mat;
1008534831adSStefano Zampini 
1009534831adSStefano Zampini   PetscFunctionBegin;
1010674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
1011534831adSStefano Zampini     /* swap pointers for local matrices */
1012534831adSStefano Zampini     temp_mat = matis->A;
1013534831adSStefano Zampini     matis->A = pcbddc->local_mat;
1014534831adSStefano Zampini     pcbddc->local_mat = temp_mat;
10153425bc38SStefano Zampini   }
10163308cffdSStefano Zampini   if (pcbddc->use_change_of_basis && x) {
1017534831adSStefano Zampini     /* Get Local boundary and apply transformation of basis to solution vector */
1018534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1019534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1020534831adSStefano Zampini     /* from modified basis to original basis */
1021534831adSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
1022534831adSStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
1023534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1024534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1025534831adSStefano Zampini   }
10263972b0daSStefano Zampini   /* add solution removed in presolve */
10273425bc38SStefano Zampini   if (x) {
10283425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
10293425bc38SStefano Zampini   }
1030fb223d50SStefano Zampini   /* restore rhs to its original state */
1031fb223d50SStefano Zampini   if (rhs) {
1032fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1033fb223d50SStefano Zampini   }
1034534831adSStefano Zampini   PetscFunctionReturn(0);
1035534831adSStefano Zampini }
1036534831adSStefano Zampini /* -------------------------------------------------------------------------- */
103753cdbc3dSStefano Zampini #undef __FUNCT__
103853cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
10390c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
10400c7d97c5SJed Brown /*
10410c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
10420c7d97c5SJed Brown                   by setting data structures and options.
10430c7d97c5SJed Brown 
10440c7d97c5SJed Brown    Input Parameter:
104553cdbc3dSStefano Zampini +  pc - the preconditioner context
10460c7d97c5SJed Brown 
10470c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
10480c7d97c5SJed Brown 
10490c7d97c5SJed Brown    Notes:
10500c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
10510c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
10520c7d97c5SJed Brown */
105353cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
10540c7d97c5SJed Brown {
10550c7d97c5SJed Brown   PetscErrorCode   ierr;
10560c7d97c5SJed Brown   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
1057f4ddd8eeSStefano Zampini   MatNullSpace     nearnullspace;
1058674ae819SStefano Zampini   MatStructure     flag;
1059674ae819SStefano Zampini   PetscBool        computeis,computetopography,computesolvers;
1060165b64e2SStefano Zampini   PetscBool        new_nearnullspace_provided;
10610c7d97c5SJed Brown 
10620c7d97c5SJed Brown   PetscFunctionBegin;
1063f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1064f4ddd8eeSStefano Zampini   /* PCIS does not support MatStructure flags different from SAME_PRECONDITIONER */
10653b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1066fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
1067f4ddd8eeSStefano Zampini 
1068f4ddd8eeSStefano Zampini   /* split work */
1069674ae819SStefano Zampini   if (pc->setupcalled) {
1070674ae819SStefano Zampini     computeis = PETSC_FALSE;
1071674ae819SStefano Zampini     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
1072674ae819SStefano Zampini     if (flag == SAME_PRECONDITIONER) {
1073f4ddd8eeSStefano Zampini       PetscFunctionReturn(0);
1074674ae819SStefano Zampini     } else if (flag == SAME_NONZERO_PATTERN) {
1075674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1076674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1077674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1078674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1079674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1080674ae819SStefano Zampini     }
1081674ae819SStefano Zampini   } else {
1082674ae819SStefano Zampini     computeis = PETSC_TRUE;
1083674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1084674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1085674ae819SStefano Zampini   }
1086fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
1087fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
1088fb180af4SStefano Zampini   }
1089f4ddd8eeSStefano Zampini 
1090f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
1091*70cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
1092*70cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
1093f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
1094f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1095*70cf5478SStefano Zampini     }
1096f4ddd8eeSStefano Zampini     if (pcbddc->current_level) {
1097f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
1098f4ddd8eeSStefano Zampini     }
1099f4ddd8eeSStefano Zampini   }
1100f4ddd8eeSStefano Zampini 
1101fcd409f5SStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
1102674ae819SStefano Zampini   if (computeis) {
1103fcd409f5SStefano Zampini     /* HACK INTO PCIS */
1104fcd409f5SStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
1105fcd409f5SStefano Zampini     pcis->computesolvers = PETSC_FALSE;
1106674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
110739e2fb2aSStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr);
1108674ae819SStefano Zampini   }
1109f4ddd8eeSStefano Zampini 
1110f4ddd8eeSStefano Zampini   /* Analyze interface */
1111674ae819SStefano Zampini   if (computetopography) {
1112674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1113fb8d54d4SStefano Zampini   }
1114fb8d54d4SStefano Zampini 
1115f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1116fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1117f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1118f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1119f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1120f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1121f4ddd8eeSStefano Zampini     } else {
1122f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1123f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1124f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1125165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1126f4ddd8eeSStefano Zampini         PetscInt         i;
1127165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1128165b64e2SStefano Zampini         PetscObjectState state;
1129165b64e2SStefano Zampini         PetscInt         nnsp_size;
1130165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1131f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1132f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1133165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1134f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1135f4ddd8eeSStefano Zampini             break;
1136f4ddd8eeSStefano Zampini           }
1137f4ddd8eeSStefano Zampini         }
1138f4ddd8eeSStefano Zampini       }
1139f4ddd8eeSStefano Zampini     }
1140f4ddd8eeSStefano Zampini   } else {
1141f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1142f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1143f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1144f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1145f4ddd8eeSStefano Zampini     }
1146f4ddd8eeSStefano Zampini   }
1147f4ddd8eeSStefano Zampini 
1148f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1149727cdba6SStefano Zampini   /* reset primal space flags */
1150f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1151727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
1152fb8d54d4SStefano Zampini   if (computetopography || new_nearnullspace_provided) {
1153727cdba6SStefano Zampini     /* It also sets the primal space flags */
1154674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1155e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1156f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
1157674ae819SStefano Zampini   }
1158fb8d54d4SStefano Zampini 
1159f4ddd8eeSStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
1160674ae819SStefano Zampini     /* reset data */
1161674ae819SStefano Zampini     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
116299cc7994SStefano Zampini     /* Create coarse and local stuffs */
116399cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1164674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
11650c7d97c5SJed Brown   }
1166*70cf5478SStefano Zampini 
11672b510759SStefano Zampini   if (pcbddc->dbg_flag && pcbddc->current_level) {
11682b510759SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
11692b510759SStefano Zampini   }
11700c7d97c5SJed Brown   PetscFunctionReturn(0);
11710c7d97c5SJed Brown }
11720c7d97c5SJed Brown 
11730c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
11740c7d97c5SJed Brown /*
11750c7d97c5SJed Brown    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
11760c7d97c5SJed Brown 
11770c7d97c5SJed Brown    Input Parameters:
11780c7d97c5SJed Brown .  pc - the preconditioner context
11790c7d97c5SJed Brown .  r - input vector (global)
11800c7d97c5SJed Brown 
11810c7d97c5SJed Brown    Output Parameter:
11820c7d97c5SJed Brown .  z - output vector (global)
11830c7d97c5SJed Brown 
11840c7d97c5SJed Brown    Application Interface Routine: PCApply()
11850c7d97c5SJed Brown  */
11860c7d97c5SJed Brown #undef __FUNCT__
11870c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
118853cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
11890c7d97c5SJed Brown {
11900c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
11910c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
11920c7d97c5SJed Brown   PetscErrorCode    ierr;
11933b03a366Sstefano_zampini   const PetscScalar one = 1.0;
11943b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
11952617d88aSStefano Zampini   const PetscScalar zero = 0.0;
11960c7d97c5SJed Brown 
11970c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
11980c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
11998eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
12000c7d97c5SJed Brown 
12010c7d97c5SJed Brown   PetscFunctionBegin;
120285c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
12030c7d97c5SJed Brown     /* First Dirichlet solve */
12040c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12050c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120653cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
12070c7d97c5SJed Brown     /*
12080c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1209674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1210674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
12110c7d97c5SJed Brown     */
12120c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
12130c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
12148eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
12150c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
12160c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
12170c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12180c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1219674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1220b76ba322SStefano Zampini   } else {
12210bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1222b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1223674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1224b76ba322SStefano Zampini   }
1225b76ba322SStefano Zampini 
12262617d88aSStefano Zampini   /* Apply interface preconditioner
12272617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
12282617d88aSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
12292617d88aSStefano Zampini 
1230674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1231674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
12320c7d97c5SJed Brown 
12333b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
12340c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12350c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12360c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
12378eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1238df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1239df187020SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1240df187020SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY (pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1241df187020SStefano Zampini   ierr = VecAXPY (pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
12420c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12430c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12440c7d97c5SJed Brown   PetscFunctionReturn(0);
12450c7d97c5SJed Brown }
1246da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1247674ae819SStefano Zampini 
1248da1bb401SStefano Zampini #undef __FUNCT__
1249da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1250da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1251da1bb401SStefano Zampini {
1252da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1253da1bb401SStefano Zampini   PetscErrorCode ierr;
1254da1bb401SStefano Zampini 
1255da1bb401SStefano Zampini   PetscFunctionBegin;
1256da1bb401SStefano Zampini   /* free data created by PCIS */
1257da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1258674ae819SStefano Zampini   /* free BDDC custom data  */
1259674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1260674ae819SStefano Zampini   /* destroy objects related to topography */
1261674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1262674ae819SStefano Zampini   /* free allocated graph structure */
1263da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1264674ae819SStefano Zampini   /* free data for scaling operator */
1265674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1266674ae819SStefano Zampini   /* free solvers stuff */
1267674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
126862a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
126962a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
127062a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
127139e2fb2aSStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr);
12723425bc38SStefano Zampini   /* remove functions */
1273674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1274bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
12752b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1276b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
12772b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1278bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1279bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
128082ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1281bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
128282ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1283bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
128482ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1285bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1286785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1287bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
128863602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1289bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1290bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1291bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1292bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1293674ae819SStefano Zampini   /* Free the private data structure */
1294674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1295da1bb401SStefano Zampini   PetscFunctionReturn(0);
1296da1bb401SStefano Zampini }
12973425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
12981e6b0712SBarry Smith 
12993425bc38SStefano Zampini #undef __FUNCT__
13003425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
13013425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
13023425bc38SStefano Zampini {
1303674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
13043425bc38SStefano Zampini   PC_IS*         pcis;
13053425bc38SStefano Zampini   PC_BDDC*       pcbddc;
13063425bc38SStefano Zampini   PetscErrorCode ierr;
13070c7d97c5SJed Brown 
13083425bc38SStefano Zampini   PetscFunctionBegin;
13093425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
13103425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
13113425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
13123425bc38SStefano Zampini 
13133425bc38SStefano Zampini   /* change of basis for physical rhs if needed
13143425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
13153308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
13163425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
13173425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13183425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1319fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1320fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
13213425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13223425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1323674ae819SStefano Zampini   /* Apply partition of unity */
13243425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1325674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
13268eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
13273425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
13283425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
13293425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
13303425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
13313425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
13323425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13333425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1334674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
13353425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13363425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13373425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
13383425bc38SStefano Zampini   }
13393425bc38SStefano Zampini   /* BDDC rhs */
13403425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
13418eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
13423425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
13433425bc38SStefano Zampini   }
13443425bc38SStefano Zampini   /* apply BDDC */
13453425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
13463425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
13473425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
13483425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
13493425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13503425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13513425bc38SStefano Zampini   /* restore original rhs */
13523425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
13533425bc38SStefano Zampini   PetscFunctionReturn(0);
13543425bc38SStefano Zampini }
13551e6b0712SBarry Smith 
13563425bc38SStefano Zampini #undef __FUNCT__
13573425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
13583425bc38SStefano Zampini /*@
135928509bceSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
13603425bc38SStefano Zampini 
13613425bc38SStefano Zampini    Collective
13623425bc38SStefano Zampini 
13633425bc38SStefano Zampini    Input Parameters:
136428509bceSStefano Zampini +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
136528509bceSStefano Zampini .  standard_rhs - the right-hand side for your linear system
13663425bc38SStefano Zampini 
13673425bc38SStefano Zampini    Output Parameters:
136828509bceSStefano Zampini -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
13693425bc38SStefano Zampini 
13703425bc38SStefano Zampini    Level: developer
13713425bc38SStefano Zampini 
13723425bc38SStefano Zampini    Notes:
13733425bc38SStefano Zampini 
137428509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
13753425bc38SStefano Zampini @*/
13763425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
13773425bc38SStefano Zampini {
1378674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
13793425bc38SStefano Zampini   PetscErrorCode ierr;
13803425bc38SStefano Zampini 
13813425bc38SStefano Zampini   PetscFunctionBegin;
13823425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
13833425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
13843425bc38SStefano Zampini   PetscFunctionReturn(0);
13853425bc38SStefano Zampini }
13863425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
13871e6b0712SBarry Smith 
13883425bc38SStefano Zampini #undef __FUNCT__
13893425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
13903425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
13913425bc38SStefano Zampini {
1392674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
13933425bc38SStefano Zampini   PC_IS*         pcis;
13943425bc38SStefano Zampini   PC_BDDC*       pcbddc;
13953425bc38SStefano Zampini   PetscErrorCode ierr;
13963425bc38SStefano Zampini 
13973425bc38SStefano Zampini   PetscFunctionBegin;
13983425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
13993425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
14003425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
14013425bc38SStefano Zampini 
14023425bc38SStefano Zampini   /* apply B_delta^T */
14033425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14043425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14053425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
14063425bc38SStefano Zampini   /* compute rhs for BDDC application */
14073425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
14088eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
14093425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
14103425bc38SStefano Zampini   }
14113425bc38SStefano Zampini   /* apply BDDC */
14123425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
14133425bc38SStefano Zampini   /* put values into standard global vector */
14143425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14153425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14168eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
14173425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
14183425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
14193425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
14203425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
14213425bc38SStefano Zampini   }
14223425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14233425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14243425bc38SStefano Zampini   /* final change of basis if needed
14253425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
14263308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
14273425bc38SStefano Zampini   PetscFunctionReturn(0);
14283425bc38SStefano Zampini }
14291e6b0712SBarry Smith 
14303425bc38SStefano Zampini #undef __FUNCT__
14313425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
14323425bc38SStefano Zampini /*@
143328509bceSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
14343425bc38SStefano Zampini 
14353425bc38SStefano Zampini    Collective
14363425bc38SStefano Zampini 
14373425bc38SStefano Zampini    Input Parameters:
143828509bceSStefano Zampini +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
143928509bceSStefano Zampini .  fetidp_flux_sol - the solution of the FETIDP linear system
14403425bc38SStefano Zampini 
14413425bc38SStefano Zampini    Output Parameters:
144228509bceSStefano Zampini -  standard_sol      - the solution defined on the physical domain
14433425bc38SStefano Zampini 
14443425bc38SStefano Zampini    Level: developer
14453425bc38SStefano Zampini 
14463425bc38SStefano Zampini    Notes:
14473425bc38SStefano Zampini 
144828509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
14493425bc38SStefano Zampini @*/
14503425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
14513425bc38SStefano Zampini {
1452674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
14533425bc38SStefano Zampini   PetscErrorCode ierr;
14543425bc38SStefano Zampini 
14553425bc38SStefano Zampini   PetscFunctionBegin;
14563425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
14573425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
14583425bc38SStefano Zampini   PetscFunctionReturn(0);
14593425bc38SStefano Zampini }
14603425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
14611e6b0712SBarry Smith 
1462f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1463f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1464f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1465f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1466674ae819SStefano Zampini 
14673425bc38SStefano Zampini #undef __FUNCT__
14683425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
14693425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
14703425bc38SStefano Zampini {
1471674ae819SStefano Zampini 
1472674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
14733425bc38SStefano Zampini   Mat            newmat;
1474674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
14753425bc38SStefano Zampini   PC             newpc;
1476ce94432eSBarry Smith   MPI_Comm       comm;
14773425bc38SStefano Zampini   PetscErrorCode ierr;
14783425bc38SStefano Zampini 
14793425bc38SStefano Zampini   PetscFunctionBegin;
1480ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
14813425bc38SStefano Zampini   /* FETIDP linear matrix */
14823425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
14833425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
14843425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
14853425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
14863425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
14873425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
14883425bc38SStefano Zampini   /* FETIDP preconditioner */
14893425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
14903425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
14913425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
14923425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
14933425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
14943425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
14953425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
14963425bc38SStefano Zampini   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
14973425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
14983425bc38SStefano Zampini   /* return pointers for objects created */
14993425bc38SStefano Zampini   *fetidp_mat=newmat;
15003425bc38SStefano Zampini   *fetidp_pc=newpc;
15013425bc38SStefano Zampini   PetscFunctionReturn(0);
15023425bc38SStefano Zampini }
15031e6b0712SBarry Smith 
15043425bc38SStefano Zampini #undef __FUNCT__
15053425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
15063425bc38SStefano Zampini /*@
150728509bceSStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
15083425bc38SStefano Zampini 
15093425bc38SStefano Zampini    Collective
15103425bc38SStefano Zampini 
15113425bc38SStefano Zampini    Input Parameters:
151228509bceSStefano Zampini +  pc - the BDDC preconditioning context already setup
151328509bceSStefano Zampini 
151428509bceSStefano Zampini    Output Parameters:
151528509bceSStefano Zampini .  fetidp_mat - shell FETIDP matrix object
151628509bceSStefano Zampini .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
151728509bceSStefano Zampini 
151828509bceSStefano Zampini    Options Database Keys:
151928509bceSStefano Zampini -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
15203425bc38SStefano Zampini 
15213425bc38SStefano Zampini    Level: developer
15223425bc38SStefano Zampini 
15233425bc38SStefano Zampini    Notes:
152428509bceSStefano Zampini      Currently the only operation provided for FETIDP matrix is MatMult
15253425bc38SStefano Zampini 
15263425bc38SStefano Zampini .seealso: PCBDDC
15273425bc38SStefano Zampini @*/
15283425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
15293425bc38SStefano Zampini {
15303425bc38SStefano Zampini   PetscErrorCode ierr;
15313425bc38SStefano Zampini 
15323425bc38SStefano Zampini   PetscFunctionBegin;
15333425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
15343425bc38SStefano Zampini   if (pc->setupcalled) {
1535516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1536f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
15373425bc38SStefano Zampini   PetscFunctionReturn(0);
15383425bc38SStefano Zampini }
15390c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1540da1bb401SStefano Zampini /*MC
1541da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
15420c7d97c5SJed Brown 
154328509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
154428509bceSStefano Zampini 
154528509bceSStefano Zampini .vb
154628509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
154728509bceSStefano 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
154828509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
154928509bceSStefano Zampini .ve
155028509bceSStefano Zampini 
155128509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
155228509bceSStefano Zampini 
1553b6fdb6dfSStefano Zampini    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
155428509bceSStefano Zampini 
155528509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
155628509bceSStefano Zampini 
1557b6fdb6dfSStefano 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.
1558b6fdb6dfSStefano Zampini 
155928509bceSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
156028509bceSStefano Zampini 
156128509bceSStefano 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
156228509bceSStefano Zampini 
156328509bceSStefano Zampini    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace.
156428509bceSStefano Zampini 
1565b6fdb6dfSStefano 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.
156628509bceSStefano Zampini 
156728509bceSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
156828509bceSStefano Zampini 
1569da1bb401SStefano Zampini    Options Database Keys:
157028509bceSStefano Zampini 
157128509bceSStefano Zampini .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
157228509bceSStefano Zampini .    -pc_bddc_use_edges <1> - use or not edges in primal space
1573b6fdb6dfSStefano Zampini .    -pc_bddc_use_faces <0> - use or not faces in primal space
157428509bceSStefano Zampini .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
157528509bceSStefano Zampini .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
157628509bceSStefano Zampini .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
157728509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
157828509bceSStefano Zampini .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
157928509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
158028509bceSStefano Zampini 
158128509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
158228509bceSStefano Zampini .vb
158328509bceSStefano Zampini       -pc_bddc_dirichlet_
158428509bceSStefano Zampini       -pc_bddc_neumann_
158528509bceSStefano Zampini       -pc_bddc_coarse_
158628509bceSStefano Zampini .ve
158728509bceSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
158828509bceSStefano Zampini 
158928509bceSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level can be specified as
159028509bceSStefano Zampini .vb
159128509bceSStefano Zampini       -pc_bddc_dirichlet_N_
159228509bceSStefano Zampini       -pc_bddc_neumann_N_
159328509bceSStefano Zampini       -pc_bddc_coarse_N_
159428509bceSStefano Zampini .ve
159528509bceSStefano Zampini    Note that level number ranges from the finest 0 to the coarsest N
1596da1bb401SStefano Zampini 
1597da1bb401SStefano Zampini    Level: intermediate
1598da1bb401SStefano Zampini 
1599b6fdb6dfSStefano Zampini    Developer notes:
160028509bceSStefano Zampini      Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose
1601da1bb401SStefano Zampini 
160228509bceSStefano Zampini      New deluxe scaling operator will be available soon.
1603da1bb401SStefano Zampini 
1604da1bb401SStefano Zampini    Contributed by Stefano Zampini
1605da1bb401SStefano Zampini 
1606da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1607da1bb401SStefano Zampini M*/
1608b2573a8aSBarry Smith 
1609da1bb401SStefano Zampini #undef __FUNCT__
1610da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
16118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1612da1bb401SStefano Zampini {
1613da1bb401SStefano Zampini   PetscErrorCode      ierr;
1614da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1615da1bb401SStefano Zampini 
1616da1bb401SStefano Zampini   PetscFunctionBegin;
1617da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1618da1bb401SStefano Zampini   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1619da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1620da1bb401SStefano Zampini 
1621da1bb401SStefano Zampini   /* create PCIS data structure */
1622da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1623da1bb401SStefano Zampini 
162447d04d0dSStefano Zampini   /* BDDC customization */
162547d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
162647d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
162747d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
162847d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
162947d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
163047d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
163147d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
163247d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
163347d04d0dSStefano Zampini 
163439e2fb2aSStefano Zampini   pcbddc->BtoNmap                    = 0;
1635727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
1636e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
1637e9189074SStefano Zampini   pcbddc->n_actual_vertices          = 0;
1638e9189074SStefano Zampini   pcbddc->n_constraints              = 0;
1639727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
1640fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
1641f4ddd8eeSStefano Zampini   pcbddc->coarse_size                = 0;
1642f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
1643727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
1644f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
1645f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
1646f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
1647674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
16480bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
16493972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1650534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1651534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1652534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1653da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1654da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1655da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1656da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1657da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
165815aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
165915aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1660da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1661da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1662da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1663da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1664da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1665da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1666da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1667da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1668da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1669da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
1670785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
1671785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
1672785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
167360d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
167460d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
167563602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal           = 0;
1676da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
167763602bcaSStefano Zampini   pcbddc->ISForDofsLocal             = 0;
1678da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
167985c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
168047d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
168147d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
16824fad6a16SStefano Zampini   pcbddc->current_level              = 0;
16832b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1684da1bb401SStefano Zampini 
1685674ae819SStefano Zampini   /* create local graph structure */
1686674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1687674ae819SStefano Zampini 
1688674ae819SStefano Zampini   /* scaling */
1689674ae819SStefano Zampini   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1690674ae819SStefano Zampini   pcbddc->work_scaling               = 0;
1691da1bb401SStefano Zampini 
1692da1bb401SStefano Zampini   /* function pointers */
1693da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
1694da1bb401SStefano Zampini   pc->ops->applytranspose      = 0;
1695da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1696da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1697da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1698da1bb401SStefano Zampini   pc->ops->view                = 0;
1699da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1700da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1701da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1702534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1703534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1704da1bb401SStefano Zampini 
1705da1bb401SStefano Zampini   /* composing function */
1706674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1707bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
17082b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1709b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
17102b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1711bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1712bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
171382ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1714bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
171582ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1716bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
171782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1718bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
171982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1720bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
172163602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
1722bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1723bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1724bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1725bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1726da1bb401SStefano Zampini   PetscFunctionReturn(0);
1727da1bb401SStefano Zampini }
17283425bc38SStefano Zampini 
1729