xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision d7de1dd987cf5fb7071be962e706ec2af58154b2)
153cdbc3dSStefano Zampini /* TODOLIST
2eb97c9d2SStefano Zampini 
3eb97c9d2SStefano Zampini    ConstraintsSetup
4eb97c9d2SStefano Zampini    - tolerances for constraints as an option (take care of single precision!)
5eb97c9d2SStefano Zampini 
6eb97c9d2SStefano Zampini    Solvers
7a0d3c3abSStefano Zampini    - Add support for cholesky for coarse solver (similar to local solvers)
8eb97c9d2SStefano Zampini    - Propagate ksp prefixes for solvers to mat objects?
9eb97c9d2SStefano Zampini    - Propagate nearnullspace info among levels
10eb97c9d2SStefano Zampini 
11eb97c9d2SStefano Zampini    User interface
12b9b85e73SStefano Zampini    - ** DofSplitting and DM attached to pc?
13eb97c9d2SStefano Zampini 
14eb97c9d2SStefano Zampini    Debugging output
15b9b85e73SStefano Zampini    - * Better management of verbosity levels of debugging output
16eb97c9d2SStefano Zampini 
17eb97c9d2SStefano Zampini    Build
18eb97c9d2SStefano Zampini    - make runexe59
19eb97c9d2SStefano Zampini 
20eb97c9d2SStefano Zampini    Extra
21b9b85e73SStefano Zampini    - ** GetRid of PCBDDCApplySchur, use MatSchur instead
22b9b85e73SStefano Zampini    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
23eb97c9d2SStefano Zampini    - add support for computing h,H and related using coordinates?
24c0b83709SStefano Zampini    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
25eb97c9d2SStefano Zampini    - Better management in PCIS code
26eb97c9d2SStefano Zampini    - BDDC with MG framework?
27eb97c9d2SStefano Zampini 
28eb97c9d2SStefano Zampini    FETIDP
29eb97c9d2SStefano Zampini    - Move FETIDP code to its own classes
30eb97c9d2SStefano Zampini 
31eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
32eb97c9d2SStefano Zampini    - Provide general case for subassembling
33b9b85e73SStefano Zampini    - *** Preallocation routines in MatISGetMPIAXAIJ
34eb97c9d2SStefano Zampini 
3553cdbc3dSStefano Zampini */
360c7d97c5SJed Brown 
3753cdbc3dSStefano Zampini /* ----------------------------------------------------------------------------------------------------------------------------------------------
380c7d97c5SJed Brown    Implementation of BDDC preconditioner based on:
390c7d97c5SJed Brown    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
4053cdbc3dSStefano Zampini    ---------------------------------------------------------------------------------------------------------------------------------------------- */
4153cdbc3dSStefano Zampini 
42ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
43ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
443b03a366Sstefano_zampini #include <petscblaslapack.h>
45674ae819SStefano Zampini 
460c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
470c7d97c5SJed Brown #undef __FUNCT__
480c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
490c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc)
500c7d97c5SJed Brown {
510c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
520c7d97c5SJed Brown   PetscErrorCode ierr;
530c7d97c5SJed Brown 
540c7d97c5SJed Brown   PetscFunctionBegin;
550c7d97c5SJed Brown   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
568eeda7d8SStefano Zampini   /* Verbose debugging */
578eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
588eeda7d8SStefano Zampini   /* Primal space cumstomization */
5908a5cf49SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);CHKERRQ(ierr);
608eeda7d8SStefano 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);
618eeda7d8SStefano 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);
628eeda7d8SStefano 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);
638eeda7d8SStefano Zampini   /* Change of basis */
64b9b85e73SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not internal change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr);
65b9b85e73SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not internal change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr);
66674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
67674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
68674ae819SStefano Zampini   }
698eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
708eeda7d8SStefano 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);
710298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
722b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
73323d395dSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);CHKERRQ(ierr);
74674ae819SStefano 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);
750c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
760c7d97c5SJed Brown   PetscFunctionReturn(0);
770c7d97c5SJed Brown }
780c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
79674ae819SStefano Zampini #undef __FUNCT__
80b9b85e73SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisLocalMat_BDDC"
81b9b85e73SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisLocalMat_BDDC(PC pc, Mat change)
82b9b85e73SStefano Zampini {
83b9b85e73SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
84b9b85e73SStefano Zampini   PetscErrorCode ierr;
85b9b85e73SStefano Zampini 
86b9b85e73SStefano Zampini   PetscFunctionBegin;
87b9b85e73SStefano Zampini   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
88b9b85e73SStefano Zampini   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
89b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix = change;
90b9b85e73SStefano Zampini   PetscFunctionReturn(0);
91b9b85e73SStefano Zampini }
92b9b85e73SStefano Zampini #undef __FUNCT__
93b9b85e73SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisLocalMat"
94b9b85e73SStefano Zampini /*@
95b9b85e73SStefano Zampini  PCBDDCSetChangeOfBasisLocalMat - Set user defined change of basis for local boundary dofs
96b9b85e73SStefano Zampini 
97b9b85e73SStefano Zampini    Collective on PC
98b9b85e73SStefano Zampini 
99b9b85e73SStefano Zampini    Input Parameters:
100b9b85e73SStefano Zampini +  pc - the preconditioning context
101b9b85e73SStefano Zampini -  change - the local change of basis matrix, either in local (internal + boundary) or in local boundary numbering
102b9b85e73SStefano Zampini 
103b9b85e73SStefano Zampini    Level: intermediate
104b9b85e73SStefano Zampini 
105b9b85e73SStefano Zampini    Notes:
106b9b85e73SStefano Zampini 
107b9b85e73SStefano Zampini .seealso: PCBDDC
108b9b85e73SStefano Zampini @*/
109b9b85e73SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisLocalMat(PC pc, Mat change)
110b9b85e73SStefano Zampini {
111b9b85e73SStefano Zampini   PetscErrorCode ierr;
112b9b85e73SStefano Zampini 
113b9b85e73SStefano Zampini   PetscFunctionBegin;
114b9b85e73SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
115b9b85e73SStefano Zampini   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
116b9b85e73SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisLocalMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
117b9b85e73SStefano Zampini   PetscFunctionReturn(0);
118b9b85e73SStefano Zampini }
119b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */
120b9b85e73SStefano Zampini #undef __FUNCT__
121674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
122674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
123674ae819SStefano Zampini {
124674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
125674ae819SStefano Zampini   PetscErrorCode ierr;
1261e6b0712SBarry Smith 
127674ae819SStefano Zampini   PetscFunctionBegin;
128674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
129674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
130674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
131674ae819SStefano Zampini   PetscFunctionReturn(0);
132674ae819SStefano Zampini }
133674ae819SStefano Zampini #undef __FUNCT__
134674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
135674ae819SStefano Zampini /*@
13628509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
137674ae819SStefano Zampini 
138674ae819SStefano Zampini    Not collective
139674ae819SStefano Zampini 
140674ae819SStefano Zampini    Input Parameters:
141674ae819SStefano Zampini +  pc - the preconditioning context
14228509bceSStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering
143674ae819SStefano Zampini 
144674ae819SStefano Zampini    Level: intermediate
145674ae819SStefano Zampini 
146674ae819SStefano Zampini    Notes:
147674ae819SStefano Zampini 
148674ae819SStefano Zampini .seealso: PCBDDC
149674ae819SStefano Zampini @*/
150674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
151674ae819SStefano Zampini {
152674ae819SStefano Zampini   PetscErrorCode ierr;
153674ae819SStefano Zampini 
154674ae819SStefano Zampini   PetscFunctionBegin;
155674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
156674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
157674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
158674ae819SStefano Zampini   PetscFunctionReturn(0);
159674ae819SStefano Zampini }
160674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1610c7d97c5SJed Brown #undef __FUNCT__
1624fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1634fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1644fad6a16SStefano Zampini {
1654fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1664fad6a16SStefano Zampini 
1674fad6a16SStefano Zampini   PetscFunctionBegin;
1684fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1694fad6a16SStefano Zampini   PetscFunctionReturn(0);
1704fad6a16SStefano Zampini }
1711e6b0712SBarry Smith 
1724fad6a16SStefano Zampini #undef __FUNCT__
1734fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1744fad6a16SStefano Zampini /*@
17528509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
1764fad6a16SStefano Zampini 
1774fad6a16SStefano Zampini    Logically collective on PC
1784fad6a16SStefano Zampini 
1794fad6a16SStefano Zampini    Input Parameters:
1804fad6a16SStefano Zampini +  pc - the preconditioning context
18128509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
1824fad6a16SStefano Zampini 
18328509bceSStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
1844fad6a16SStefano Zampini 
1854fad6a16SStefano Zampini    Level: intermediate
1864fad6a16SStefano Zampini 
1874fad6a16SStefano Zampini    Notes:
1884fad6a16SStefano Zampini 
1894fad6a16SStefano Zampini .seealso: PCBDDC
1904fad6a16SStefano Zampini @*/
1914fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
1924fad6a16SStefano Zampini {
1934fad6a16SStefano Zampini   PetscErrorCode ierr;
1944fad6a16SStefano Zampini 
1954fad6a16SStefano Zampini   PetscFunctionBegin;
1964fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1972b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
1984fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
1994fad6a16SStefano Zampini   PetscFunctionReturn(0);
2004fad6a16SStefano Zampini }
2012b510759SStefano Zampini 
202b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
2032b510759SStefano Zampini #undef __FUNCT__
204b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
205b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
206b8ffe317SStefano Zampini {
207b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
208b8ffe317SStefano Zampini 
209b8ffe317SStefano Zampini   PetscFunctionBegin;
21085c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
211b8ffe317SStefano Zampini   PetscFunctionReturn(0);
212b8ffe317SStefano Zampini }
213b8ffe317SStefano Zampini 
214b8ffe317SStefano Zampini #undef __FUNCT__
215b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
216b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
2172b510759SStefano Zampini {
2182b510759SStefano Zampini   PetscErrorCode ierr;
2192b510759SStefano Zampini 
2202b510759SStefano Zampini   PetscFunctionBegin;
2212b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
222b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
223b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
2242b510759SStefano Zampini   PetscFunctionReturn(0);
2252b510759SStefano Zampini }
2261e6b0712SBarry Smith 
2274fad6a16SStefano Zampini #undef __FUNCT__
2282b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
2292b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
2304fad6a16SStefano Zampini {
2314fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2324fad6a16SStefano Zampini 
2334fad6a16SStefano Zampini   PetscFunctionBegin;
2342b510759SStefano Zampini   pcbddc->current_level = level;
2354fad6a16SStefano Zampini   PetscFunctionReturn(0);
2364fad6a16SStefano Zampini }
2371e6b0712SBarry Smith 
2384fad6a16SStefano Zampini #undef __FUNCT__
239b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
240b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
241b8ffe317SStefano Zampini {
242b8ffe317SStefano Zampini   PetscErrorCode ierr;
243b8ffe317SStefano Zampini 
244b8ffe317SStefano Zampini   PetscFunctionBegin;
245b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
246b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
247b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
248b8ffe317SStefano Zampini   PetscFunctionReturn(0);
249b8ffe317SStefano Zampini }
250b8ffe317SStefano Zampini 
251b8ffe317SStefano Zampini #undef __FUNCT__
2522b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
2532b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
2542b510759SStefano Zampini {
2552b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2562b510759SStefano Zampini 
2572b510759SStefano Zampini   PetscFunctionBegin;
2582b510759SStefano Zampini   pcbddc->max_levels = levels;
2592b510759SStefano Zampini   PetscFunctionReturn(0);
2602b510759SStefano Zampini }
2612b510759SStefano Zampini 
2622b510759SStefano Zampini #undef __FUNCT__
2632b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2644fad6a16SStefano Zampini /*@
26528509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
2664fad6a16SStefano Zampini 
2674fad6a16SStefano Zampini    Logically collective on PC
2684fad6a16SStefano Zampini 
2694fad6a16SStefano Zampini    Input Parameters:
2704fad6a16SStefano Zampini +  pc - the preconditioning context
27128509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
2724fad6a16SStefano Zampini 
27328509bceSStefano Zampini    Default value is 0, i.e. traditional one-level BDDC
2744fad6a16SStefano Zampini 
2754fad6a16SStefano Zampini    Level: intermediate
2764fad6a16SStefano Zampini 
2774fad6a16SStefano Zampini    Notes:
2784fad6a16SStefano Zampini 
2794fad6a16SStefano Zampini .seealso: PCBDDC
2804fad6a16SStefano Zampini @*/
2812b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
2824fad6a16SStefano Zampini {
2834fad6a16SStefano Zampini   PetscErrorCode ierr;
2844fad6a16SStefano Zampini 
2854fad6a16SStefano Zampini   PetscFunctionBegin;
2864fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2872b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
288312be037SStefano Zampini   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
2892b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
2904fad6a16SStefano Zampini   PetscFunctionReturn(0);
2914fad6a16SStefano Zampini }
2924fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
2931e6b0712SBarry Smith 
2944fad6a16SStefano Zampini #undef __FUNCT__
2950bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
2960bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
2970bdf917eSStefano Zampini {
2980bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2990bdf917eSStefano Zampini   PetscErrorCode ierr;
3000bdf917eSStefano Zampini 
3010bdf917eSStefano Zampini   PetscFunctionBegin;
3020bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
3030bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
3040bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
3050bdf917eSStefano Zampini   PetscFunctionReturn(0);
3060bdf917eSStefano Zampini }
3071e6b0712SBarry Smith 
3080bdf917eSStefano Zampini #undef __FUNCT__
3090bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
3100bdf917eSStefano Zampini /*@
31128509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
3120bdf917eSStefano Zampini 
3130bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
3140bdf917eSStefano Zampini 
3150bdf917eSStefano Zampini    Input Parameters:
3160bdf917eSStefano Zampini +  pc - the preconditioning context
31728509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
3180bdf917eSStefano Zampini 
3190bdf917eSStefano Zampini    Level: intermediate
3200bdf917eSStefano Zampini 
3210bdf917eSStefano Zampini    Notes:
3220bdf917eSStefano Zampini 
3230bdf917eSStefano Zampini .seealso: PCBDDC
3240bdf917eSStefano Zampini @*/
3250bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
3260bdf917eSStefano Zampini {
3270bdf917eSStefano Zampini   PetscErrorCode ierr;
3280bdf917eSStefano Zampini 
3290bdf917eSStefano Zampini   PetscFunctionBegin;
3300bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
331674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
3322b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
3330bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
3340bdf917eSStefano Zampini   PetscFunctionReturn(0);
3350bdf917eSStefano Zampini }
3360bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
3371e6b0712SBarry Smith 
3380bdf917eSStefano Zampini #undef __FUNCT__
3393b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
3403b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
3413b03a366Sstefano_zampini {
3423b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3433b03a366Sstefano_zampini   PetscErrorCode ierr;
3443b03a366Sstefano_zampini 
3453b03a366Sstefano_zampini   PetscFunctionBegin;
346785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
347785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3483b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
34936e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
35036e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
351fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3523b03a366Sstefano_zampini   PetscFunctionReturn(0);
3533b03a366Sstefano_zampini }
3541e6b0712SBarry Smith 
3553b03a366Sstefano_zampini #undef __FUNCT__
3563b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
3573b03a366Sstefano_zampini /*@
35828509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
3593b03a366Sstefano_zampini 
360785d1243SStefano Zampini    Collective
3613b03a366Sstefano_zampini 
3623b03a366Sstefano_zampini    Input Parameters:
3633b03a366Sstefano_zampini +  pc - the preconditioning context
364785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
3653b03a366Sstefano_zampini 
3663b03a366Sstefano_zampini    Level: intermediate
3673b03a366Sstefano_zampini 
368785d1243SStefano Zampini    Notes: Any process can list any global node
3693b03a366Sstefano_zampini 
3703b03a366Sstefano_zampini .seealso: PCBDDC
3713b03a366Sstefano_zampini @*/
3723b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3733b03a366Sstefano_zampini {
3743b03a366Sstefano_zampini   PetscErrorCode ierr;
3753b03a366Sstefano_zampini 
3763b03a366Sstefano_zampini   PetscFunctionBegin;
3773b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
378674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
379785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
3803b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
3813b03a366Sstefano_zampini   PetscFunctionReturn(0);
3823b03a366Sstefano_zampini }
3833b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
3841e6b0712SBarry Smith 
3853b03a366Sstefano_zampini #undef __FUNCT__
38682ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
38782ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
3880c7d97c5SJed Brown {
3890c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3900c7d97c5SJed Brown   PetscErrorCode ierr;
3910c7d97c5SJed Brown 
3920c7d97c5SJed Brown   PetscFunctionBegin;
393785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
394785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3950c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
3960c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
397785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
3987d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3990c7d97c5SJed Brown   PetscFunctionReturn(0);
4000c7d97c5SJed Brown }
4010c7d97c5SJed Brown 
4020c7d97c5SJed Brown #undef __FUNCT__
40382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
4049c0446d6SStefano Zampini /*@
40582ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
4069c0446d6SStefano Zampini 
407785d1243SStefano Zampini    Collective
40853cdbc3dSStefano Zampini 
40953cdbc3dSStefano Zampini    Input Parameters:
41053cdbc3dSStefano Zampini +  pc - the preconditioning context
41182ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
41253cdbc3dSStefano Zampini 
41353cdbc3dSStefano Zampini    Level: intermediate
41453cdbc3dSStefano Zampini 
4159c0446d6SStefano Zampini    Notes:
41653cdbc3dSStefano Zampini 
41753cdbc3dSStefano Zampini .seealso: PCBDDC
41853cdbc3dSStefano Zampini @*/
41982ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
4200c7d97c5SJed Brown {
4210c7d97c5SJed Brown   PetscErrorCode ierr;
4220c7d97c5SJed Brown 
4230c7d97c5SJed Brown   PetscFunctionBegin;
4240c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4250c7d97c5SJed Brown   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
42682ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
42782ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4280c7d97c5SJed Brown   PetscFunctionReturn(0);
4290c7d97c5SJed Brown }
4300c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
4310c7d97c5SJed Brown 
4320c7d97c5SJed Brown #undef __FUNCT__
4330c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
43453cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
4350c7d97c5SJed Brown {
4360c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
43753cdbc3dSStefano Zampini   PetscErrorCode ierr;
4380c7d97c5SJed Brown 
4390c7d97c5SJed Brown   PetscFunctionBegin;
440785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
441785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
44253cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
44336e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
44436e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
445fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4460c7d97c5SJed Brown   PetscFunctionReturn(0);
4470c7d97c5SJed Brown }
4481e6b0712SBarry Smith 
4490c7d97c5SJed Brown #undef __FUNCT__
4500c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
45157527edcSJed Brown /*@
45228509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
45357527edcSJed Brown 
454785d1243SStefano Zampini    Collective
45557527edcSJed Brown 
45657527edcSJed Brown    Input Parameters:
45757527edcSJed Brown +  pc - the preconditioning context
458785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
45957527edcSJed Brown 
46057527edcSJed Brown    Level: intermediate
46157527edcSJed Brown 
462785d1243SStefano Zampini    Notes: Any process can list any global node
46357527edcSJed Brown 
46457527edcSJed Brown .seealso: PCBDDC
46557527edcSJed Brown @*/
46653cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
4670c7d97c5SJed Brown {
4680c7d97c5SJed Brown   PetscErrorCode ierr;
4690c7d97c5SJed Brown 
4700c7d97c5SJed Brown   PetscFunctionBegin;
4710c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
472674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
473785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
47453cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
47553cdbc3dSStefano Zampini   PetscFunctionReturn(0);
47653cdbc3dSStefano Zampini }
47753cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
4781e6b0712SBarry Smith 
47953cdbc3dSStefano Zampini #undef __FUNCT__
48082ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
48182ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
4820c7d97c5SJed Brown {
4830c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4840c7d97c5SJed Brown   PetscErrorCode ierr;
4850c7d97c5SJed Brown 
4860c7d97c5SJed Brown   PetscFunctionBegin;
487785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
488785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
4890c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
4900c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
491785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
4927d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4930c7d97c5SJed Brown   PetscFunctionReturn(0);
4940c7d97c5SJed Brown }
4950c7d97c5SJed Brown 
4960c7d97c5SJed Brown #undef __FUNCT__
49782ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
4980c7d97c5SJed Brown /*@
49982ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
5000c7d97c5SJed Brown 
501785d1243SStefano Zampini    Collective
5020c7d97c5SJed Brown 
5030c7d97c5SJed Brown    Input Parameters:
5040c7d97c5SJed Brown +  pc - the preconditioning context
50582ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
5060c7d97c5SJed Brown 
5070c7d97c5SJed Brown    Level: intermediate
5080c7d97c5SJed Brown 
5090c7d97c5SJed Brown    Notes:
5100c7d97c5SJed Brown 
5110c7d97c5SJed Brown .seealso: PCBDDC
5120c7d97c5SJed Brown @*/
51382ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
5140c7d97c5SJed Brown {
5150c7d97c5SJed Brown   PetscErrorCode ierr;
5160c7d97c5SJed Brown 
5170c7d97c5SJed Brown   PetscFunctionBegin;
5180c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5190c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
52082ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
52182ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
52253cdbc3dSStefano Zampini   PetscFunctionReturn(0);
52353cdbc3dSStefano Zampini }
52453cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
52553cdbc3dSStefano Zampini 
52653cdbc3dSStefano Zampini #undef __FUNCT__
527da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
528da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
529da1bb401SStefano Zampini {
530da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
531da1bb401SStefano Zampini 
532da1bb401SStefano Zampini   PetscFunctionBegin;
533da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
534da1bb401SStefano Zampini   PetscFunctionReturn(0);
535da1bb401SStefano Zampini }
5361e6b0712SBarry Smith 
537da1bb401SStefano Zampini #undef __FUNCT__
538da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
539da1bb401SStefano Zampini /*@
540785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
541da1bb401SStefano Zampini 
542785d1243SStefano Zampini    Collective
543785d1243SStefano Zampini 
544785d1243SStefano Zampini    Input Parameters:
545785d1243SStefano Zampini .  pc - the preconditioning context
546785d1243SStefano Zampini 
547785d1243SStefano Zampini    Output Parameters:
548785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
549785d1243SStefano Zampini 
550785d1243SStefano Zampini    Level: intermediate
551785d1243SStefano Zampini 
552785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
553785d1243SStefano Zampini 
554785d1243SStefano Zampini .seealso: PCBDDC
555785d1243SStefano Zampini @*/
556785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
557785d1243SStefano Zampini {
558785d1243SStefano Zampini   PetscErrorCode ierr;
559785d1243SStefano Zampini 
560785d1243SStefano Zampini   PetscFunctionBegin;
561785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
562785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
563785d1243SStefano Zampini   PetscFunctionReturn(0);
564785d1243SStefano Zampini }
565785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
566785d1243SStefano Zampini 
567785d1243SStefano Zampini #undef __FUNCT__
568785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
569785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
570785d1243SStefano Zampini {
571785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
572785d1243SStefano Zampini 
573785d1243SStefano Zampini   PetscFunctionBegin;
574785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
575785d1243SStefano Zampini   PetscFunctionReturn(0);
576785d1243SStefano Zampini }
577785d1243SStefano Zampini 
578785d1243SStefano Zampini #undef __FUNCT__
57982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
580da1bb401SStefano Zampini /*@
58182ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
582da1bb401SStefano Zampini 
583785d1243SStefano Zampini    Collective
584da1bb401SStefano Zampini 
585da1bb401SStefano Zampini    Input Parameters:
58628509bceSStefano Zampini .  pc - the preconditioning context
587da1bb401SStefano Zampini 
588da1bb401SStefano Zampini    Output Parameters:
58928509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
590da1bb401SStefano Zampini 
591da1bb401SStefano Zampini    Level: intermediate
592da1bb401SStefano Zampini 
593da1bb401SStefano Zampini    Notes:
594da1bb401SStefano Zampini 
595da1bb401SStefano Zampini .seealso: PCBDDC
596da1bb401SStefano Zampini @*/
59782ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
598da1bb401SStefano Zampini {
599da1bb401SStefano Zampini   PetscErrorCode ierr;
600da1bb401SStefano Zampini 
601da1bb401SStefano Zampini   PetscFunctionBegin;
602da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
60382ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
604da1bb401SStefano Zampini   PetscFunctionReturn(0);
605da1bb401SStefano Zampini }
606da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
6071e6b0712SBarry Smith 
608da1bb401SStefano Zampini #undef __FUNCT__
60953cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
61053cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
61153cdbc3dSStefano Zampini {
61253cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
61353cdbc3dSStefano Zampini 
61453cdbc3dSStefano Zampini   PetscFunctionBegin;
61553cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
61653cdbc3dSStefano Zampini   PetscFunctionReturn(0);
61753cdbc3dSStefano Zampini }
6181e6b0712SBarry Smith 
61953cdbc3dSStefano Zampini #undef __FUNCT__
62053cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
62153cdbc3dSStefano Zampini /*@
622785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
62353cdbc3dSStefano Zampini 
624785d1243SStefano Zampini    Collective
625785d1243SStefano Zampini 
626785d1243SStefano Zampini    Input Parameters:
627785d1243SStefano Zampini .  pc - the preconditioning context
628785d1243SStefano Zampini 
629785d1243SStefano Zampini    Output Parameters:
630785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
631785d1243SStefano Zampini 
632785d1243SStefano Zampini    Level: intermediate
633785d1243SStefano Zampini 
634785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
635785d1243SStefano Zampini 
636785d1243SStefano Zampini .seealso: PCBDDC
637785d1243SStefano Zampini @*/
638785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
639785d1243SStefano Zampini {
640785d1243SStefano Zampini   PetscErrorCode ierr;
641785d1243SStefano Zampini 
642785d1243SStefano Zampini   PetscFunctionBegin;
643785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
644785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
645785d1243SStefano Zampini   PetscFunctionReturn(0);
646785d1243SStefano Zampini }
647785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
648785d1243SStefano Zampini 
649785d1243SStefano Zampini #undef __FUNCT__
650785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
651785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
652785d1243SStefano Zampini {
653785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
654785d1243SStefano Zampini 
655785d1243SStefano Zampini   PetscFunctionBegin;
656785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
657785d1243SStefano Zampini   PetscFunctionReturn(0);
658785d1243SStefano Zampini }
659785d1243SStefano Zampini 
660785d1243SStefano Zampini #undef __FUNCT__
66182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
66253cdbc3dSStefano Zampini /*@
66382ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
66453cdbc3dSStefano Zampini 
665785d1243SStefano Zampini    Collective
66653cdbc3dSStefano Zampini 
66753cdbc3dSStefano Zampini    Input Parameters:
66828509bceSStefano Zampini .  pc - the preconditioning context
66953cdbc3dSStefano Zampini 
67053cdbc3dSStefano Zampini    Output Parameters:
67128509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
67253cdbc3dSStefano Zampini 
67353cdbc3dSStefano Zampini    Level: intermediate
67453cdbc3dSStefano Zampini 
67553cdbc3dSStefano Zampini    Notes:
67653cdbc3dSStefano Zampini 
67753cdbc3dSStefano Zampini .seealso: PCBDDC
67853cdbc3dSStefano Zampini @*/
67982ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
68053cdbc3dSStefano Zampini {
68153cdbc3dSStefano Zampini   PetscErrorCode ierr;
68253cdbc3dSStefano Zampini 
68353cdbc3dSStefano Zampini   PetscFunctionBegin;
68453cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
68582ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
6860c7d97c5SJed Brown   PetscFunctionReturn(0);
6870c7d97c5SJed Brown }
68836e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
6891e6b0712SBarry Smith 
69036e030ebSStefano Zampini #undef __FUNCT__
691da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
6921a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
69336e030ebSStefano Zampini {
69436e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
695da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
696da1bb401SStefano Zampini   PetscErrorCode ierr;
69736e030ebSStefano Zampini 
69836e030ebSStefano Zampini   PetscFunctionBegin;
699674ae819SStefano Zampini   /* free old CSR */
700674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
701fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
702674ae819SStefano Zampini   /* get CSR into graph structure */
703da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
704785e854fSJed Brown     ierr = PetscMalloc1((nvtxs+1),&mat_graph->xadj);CHKERRQ(ierr);
705785e854fSJed Brown     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
706674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
707674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
708da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
7091a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
7101a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
711674ae819SStefano Zampini   }
712575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
71336e030ebSStefano Zampini   PetscFunctionReturn(0);
71436e030ebSStefano Zampini }
7151e6b0712SBarry Smith 
71636e030ebSStefano Zampini #undef __FUNCT__
717da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
71836e030ebSStefano Zampini /*@
71928509bceSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
72036e030ebSStefano Zampini 
72136e030ebSStefano Zampini    Not collective
72236e030ebSStefano Zampini 
72336e030ebSStefano Zampini    Input Parameters:
72436e030ebSStefano Zampini +  pc - the preconditioning context
72528509bceSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
72628509bceSStefano Zampini .  xadj, adjncy - the CSR graph
72728509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
72836e030ebSStefano Zampini 
72936e030ebSStefano Zampini    Level: intermediate
73036e030ebSStefano Zampini 
73136e030ebSStefano Zampini    Notes:
73236e030ebSStefano Zampini 
73328509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
73436e030ebSStefano Zampini @*/
7351a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
73636e030ebSStefano Zampini {
737575ad6abSStefano Zampini   void (*f)(void) = 0;
73836e030ebSStefano Zampini   PetscErrorCode ierr;
73936e030ebSStefano Zampini 
74036e030ebSStefano Zampini   PetscFunctionBegin;
74136e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
742674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
743*d7de1dd9SStefano Zampini   PetscValidIntPointer(adjncy,4);
744674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
745674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
746da1bb401SStefano Zampini   }
74736e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
748575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
749575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
750575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
751575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
752575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
75336e030ebSStefano Zampini   }
75436e030ebSStefano Zampini   PetscFunctionReturn(0);
75536e030ebSStefano Zampini }
7569c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
7571e6b0712SBarry Smith 
7589c0446d6SStefano Zampini #undef __FUNCT__
75963602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
76063602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
76163602bcaSStefano Zampini {
76263602bcaSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
76363602bcaSStefano Zampini   PetscInt i;
76463602bcaSStefano Zampini   PetscErrorCode ierr;
76563602bcaSStefano Zampini 
76663602bcaSStefano Zampini   PetscFunctionBegin;
76763602bcaSStefano Zampini   /* Destroy ISes if they were already set */
76863602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
76963602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
77063602bcaSStefano Zampini   }
77163602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
77263602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
77363602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
77463602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
77563602bcaSStefano Zampini   }
77663602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
77763602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
77863602bcaSStefano Zampini   /* allocate space then set */
77963602bcaSStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
78063602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
78163602bcaSStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
78263602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
78363602bcaSStefano Zampini   }
78463602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal=n_is;
78563602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
7861a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
78763602bcaSStefano Zampini   PetscFunctionReturn(0);
78863602bcaSStefano Zampini }
78963602bcaSStefano Zampini 
79063602bcaSStefano Zampini #undef __FUNCT__
79163602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
79263602bcaSStefano Zampini /*@
79363602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
79463602bcaSStefano Zampini 
79563602bcaSStefano Zampini    Collective
79663602bcaSStefano Zampini 
79763602bcaSStefano Zampini    Input Parameters:
79863602bcaSStefano Zampini +  pc - the preconditioning context
79963602bcaSStefano Zampini -  n_is - number of index sets defining the fields
80063602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in local ordering
80163602bcaSStefano Zampini 
80263602bcaSStefano Zampini    Level: intermediate
80363602bcaSStefano Zampini 
80463602bcaSStefano 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.
80563602bcaSStefano Zampini 
80663602bcaSStefano Zampini .seealso: PCBDDC
80763602bcaSStefano Zampini @*/
80863602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
80963602bcaSStefano Zampini {
81063602bcaSStefano Zampini   PetscInt       i;
81163602bcaSStefano Zampini   PetscErrorCode ierr;
81263602bcaSStefano Zampini 
81363602bcaSStefano Zampini   PetscFunctionBegin;
81463602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
81563602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
81663602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
81763602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
81863602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
81963602bcaSStefano Zampini   }
82063602bcaSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
82163602bcaSStefano Zampini   PetscFunctionReturn(0);
82263602bcaSStefano Zampini }
82363602bcaSStefano Zampini /* -------------------------------------------------------------------------- */
82463602bcaSStefano Zampini 
82563602bcaSStefano Zampini #undef __FUNCT__
8269c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
8279c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
8289c0446d6SStefano Zampini {
8299c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
8309c0446d6SStefano Zampini   PetscInt i;
8319c0446d6SStefano Zampini   PetscErrorCode ierr;
8329c0446d6SStefano Zampini 
8339c0446d6SStefano Zampini   PetscFunctionBegin;
834da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
8359c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
8369c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
8379c0446d6SStefano Zampini   }
838d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
83963602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
84063602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
84163602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
84263602bcaSStefano Zampini   }
84363602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
84463602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
845da1bb401SStefano Zampini   /* allocate space then set */
846785e854fSJed Brown   ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
8479c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
848da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
849da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
8509c0446d6SStefano Zampini   }
8519c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
85263602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8531a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
8549c0446d6SStefano Zampini   PetscFunctionReturn(0);
8559c0446d6SStefano Zampini }
8561e6b0712SBarry Smith 
8579c0446d6SStefano Zampini #undef __FUNCT__
8589c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
8599c0446d6SStefano Zampini /*@
86063602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
8619c0446d6SStefano Zampini 
86263602bcaSStefano Zampini    Collective
8639c0446d6SStefano Zampini 
8649c0446d6SStefano Zampini    Input Parameters:
8659c0446d6SStefano Zampini +  pc - the preconditioning context
86628509bceSStefano Zampini -  n_is - number of index sets defining the fields
86763602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in global ordering
8689c0446d6SStefano Zampini 
8699c0446d6SStefano Zampini    Level: intermediate
8709c0446d6SStefano Zampini 
87163602bcaSStefano Zampini    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.
8729c0446d6SStefano Zampini 
8739c0446d6SStefano Zampini .seealso: PCBDDC
8749c0446d6SStefano Zampini @*/
8759c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
8769c0446d6SStefano Zampini {
8772b510759SStefano Zampini   PetscInt       i;
8789c0446d6SStefano Zampini   PetscErrorCode ierr;
8799c0446d6SStefano Zampini 
8809c0446d6SStefano Zampini   PetscFunctionBegin;
8819c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
88263602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
8832b510759SStefano Zampini   for (i=0;i<n_is;i++) {
88463602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
88563602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
8862b510759SStefano Zampini   }
8879c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
8889c0446d6SStefano Zampini   PetscFunctionReturn(0);
8899c0446d6SStefano Zampini }
890da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
891534831adSStefano Zampini #undef __FUNCT__
892534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
893534831adSStefano Zampini /* -------------------------------------------------------------------------- */
894534831adSStefano Zampini /*
895534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
896534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
8979c0446d6SStefano Zampini 
898534831adSStefano Zampini    Input Parameter:
899534831adSStefano Zampini +  pc - the preconditioner contex
900534831adSStefano Zampini 
901534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
902534831adSStefano Zampini 
903534831adSStefano Zampini    Notes:
904534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
905534831adSStefano Zampini    the user, but instead is called by KSPSolve().
906534831adSStefano Zampini */
907534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
908534831adSStefano Zampini {
909534831adSStefano Zampini   PetscErrorCode ierr;
910534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
911534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
912534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
913534831adSStefano Zampini   Mat            temp_mat;
9143972b0daSStefano Zampini   IS             dirIS;
9153972b0daSStefano Zampini   Vec            used_vec;
91682ba6b80SStefano Zampini   PetscBool      guess_nonzero;
917534831adSStefano Zampini 
918534831adSStefano Zampini   PetscFunctionBegin;
91985c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
92085c4d303SStefano Zampini   if (ksp) {
92185c4d303SStefano Zampini     PetscBool iscg;
92285c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
92385c4d303SStefano Zampini     if (!iscg) {
92485c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
92585c4d303SStefano Zampini     }
92685c4d303SStefano Zampini   }
92785c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
92862a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
92962a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
93062a6ff1dSStefano Zampini   }
93162a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
93262a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
93362a6ff1dSStefano Zampini   }
9343972b0daSStefano Zampini   if (x) {
9353972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
9363972b0daSStefano Zampini     used_vec = x;
9373972b0daSStefano Zampini   } else {
9383972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
9393972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
9403972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9413972b0daSStefano Zampini   }
9423972b0daSStefano Zampini   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
9433972b0daSStefano Zampini   if (ksp) {
9443972b0daSStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
9453972b0daSStefano Zampini     if (!guess_nonzero) {
9463972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9473972b0daSStefano Zampini     }
9483972b0daSStefano Zampini   }
9493308cffdSStefano Zampini 
9503972b0daSStefano Zampini   /* store the original rhs */
9513972b0daSStefano Zampini   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
9523972b0daSStefano Zampini 
9533972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
954785d1243SStefano Zampini   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
95582ba6b80SStefano Zampini   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
95682ba6b80SStefano Zampini   if (rhs && dirIS) {
957785d1243SStefano Zampini     PetscInt    dirsize,i,*is_indices;
958785d1243SStefano Zampini     PetscScalar *array_x,*array_diagonal;
959785d1243SStefano Zampini 
9603972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
9613972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
9623972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9633972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9643972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9653972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
96682ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
9673972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
9683972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
9693972b0daSStefano Zampini     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9702fa5cd67SKarl Rupp     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
9713972b0daSStefano Zampini     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9723972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
9733972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
9743972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9753972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
976b76ba322SStefano Zampini 
9773972b0daSStefano Zampini     /* remove the computed solution from the rhs */
9783972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
9793972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
9803972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
9813308cffdSStefano Zampini   }
982b76ba322SStefano Zampini 
983b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
9843972b0daSStefano Zampini   if (x) {
9853972b0daSStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
9863972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
98785c4d303SStefano Zampini     if (pcbddc->use_exact_dirichlet_trick) {
988b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
989b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
990b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
991b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
992b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
993b76ba322SStefano Zampini       if (ksp) {
994b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
995b76ba322SStefano Zampini       }
996b76ba322SStefano Zampini     }
9973972b0daSStefano Zampini   }
998b76ba322SStefano Zampini 
99992e3dcfbSStefano Zampini   /* prepare MatMult and rhs for solver */
1000b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1001b76ba322SStefano Zampini     /* swap pointers for local matrices */
1002b76ba322SStefano Zampini     temp_mat = matis->A;
1003b76ba322SStefano Zampini     matis->A = pcbddc->local_mat;
1004b76ba322SStefano Zampini     pcbddc->local_mat = temp_mat;
100592e3dcfbSStefano Zampini     if (rhs) {
1006b76ba322SStefano Zampini       /* Get local rhs and apply transformation of basis */
1007b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1008b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1009b76ba322SStefano Zampini       /* from original basis to modified basis */
1010b76ba322SStefano Zampini       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
1011b76ba322SStefano Zampini       /* put back modified values into the global vec using INSERT_VALUES copy mode */
1012b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1013b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1014674ae819SStefano Zampini     }
101592e3dcfbSStefano Zampini   }
101692e3dcfbSStefano Zampini 
101792e3dcfbSStefano Zampini   /* remove nullspace if present */
10180bdf917eSStefano Zampini   if (ksp && pcbddc->NullSpace) {
1019d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
1020d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1021b76ba322SStefano Zampini   }
10220bdf917eSStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1023534831adSStefano Zampini   PetscFunctionReturn(0);
1024534831adSStefano Zampini }
1025534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1026534831adSStefano Zampini #undef __FUNCT__
1027534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
1028534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1029534831adSStefano Zampini /*
1030534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1031534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
1032534831adSStefano Zampini 
1033534831adSStefano Zampini    Input Parameter:
1034534831adSStefano Zampini +  pc - the preconditioner contex
1035534831adSStefano Zampini 
1036534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
1037534831adSStefano Zampini 
1038534831adSStefano Zampini    Notes:
1039534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
1040534831adSStefano Zampini    the user, but instead is called by KSPSolve().
1041534831adSStefano Zampini */
1042534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1043534831adSStefano Zampini {
1044534831adSStefano Zampini   PetscErrorCode ierr;
1045534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1046534831adSStefano Zampini   PC_IS          *pcis   = (PC_IS*)(pc->data);
1047534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
1048534831adSStefano Zampini   Mat            temp_mat;
1049534831adSStefano Zampini 
1050534831adSStefano Zampini   PetscFunctionBegin;
1051b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1052534831adSStefano Zampini     /* swap pointers for local matrices */
1053534831adSStefano Zampini     temp_mat = matis->A;
1054534831adSStefano Zampini     matis->A = pcbddc->local_mat;
1055534831adSStefano Zampini     pcbddc->local_mat = temp_mat;
10563425bc38SStefano Zampini   }
1057b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix && x) {
1058534831adSStefano Zampini     /* Get Local boundary and apply transformation of basis to solution vector */
1059534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1060534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1061534831adSStefano Zampini     /* from modified basis to original basis */
1062534831adSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
1063534831adSStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
1064534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1065534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1066534831adSStefano Zampini   }
10673972b0daSStefano Zampini   /* add solution removed in presolve */
10683425bc38SStefano Zampini   if (x) {
10693425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
10703425bc38SStefano Zampini   }
1071fb223d50SStefano Zampini   /* restore rhs to its original state */
1072fb223d50SStefano Zampini   if (rhs) {
1073fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1074fb223d50SStefano Zampini   }
1075534831adSStefano Zampini   PetscFunctionReturn(0);
1076534831adSStefano Zampini }
1077534831adSStefano Zampini /* -------------------------------------------------------------------------- */
107853cdbc3dSStefano Zampini #undef __FUNCT__
107953cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
10800c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
10810c7d97c5SJed Brown /*
10820c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
10830c7d97c5SJed Brown                   by setting data structures and options.
10840c7d97c5SJed Brown 
10850c7d97c5SJed Brown    Input Parameter:
108653cdbc3dSStefano Zampini +  pc - the preconditioner context
10870c7d97c5SJed Brown 
10880c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
10890c7d97c5SJed Brown 
10900c7d97c5SJed Brown    Notes:
10910c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
10920c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
10930c7d97c5SJed Brown */
109453cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
10950c7d97c5SJed Brown {
10960c7d97c5SJed Brown   PetscErrorCode   ierr;
10970c7d97c5SJed Brown   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
1098f4ddd8eeSStefano Zampini   MatNullSpace     nearnullspace;
1099674ae819SStefano Zampini   PetscBool        computeis,computetopography,computesolvers;
1100165b64e2SStefano Zampini   PetscBool        new_nearnullspace_provided;
11010c7d97c5SJed Brown 
11020c7d97c5SJed Brown   PetscFunctionBegin;
1103f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
11043b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1105fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
1106f4ddd8eeSStefano Zampini   /* split work */
1107674ae819SStefano Zampini   if (pc->setupcalled) {
1108674ae819SStefano Zampini     computeis = PETSC_FALSE;
1109d1e9a80fSBarry Smith     if (pc->flag == SAME_NONZERO_PATTERN) {
1110674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1111674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1112674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1113674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1114674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1115674ae819SStefano Zampini     }
1116674ae819SStefano Zampini   } else {
1117674ae819SStefano Zampini     computeis = PETSC_TRUE;
1118674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1119674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1120674ae819SStefano Zampini   }
1121fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
1122fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
1123fb180af4SStefano Zampini   }
1124f4ddd8eeSStefano Zampini 
1125f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
112670cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
112770cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
112858a03d70SStefano Zampini       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1129f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1130f4ddd8eeSStefano Zampini     }
113158a03d70SStefano Zampini     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1132f4ddd8eeSStefano Zampini   }
1133f4ddd8eeSStefano Zampini 
1134fcd409f5SStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
1135674ae819SStefano Zampini   if (computeis) {
1136fcd409f5SStefano Zampini     /* HACK INTO PCIS */
1137fcd409f5SStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
1138fcd409f5SStefano Zampini     pcis->computesolvers = PETSC_FALSE;
1139674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
114039e2fb2aSStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr);
1141674ae819SStefano Zampini   }
1142f4ddd8eeSStefano Zampini 
1143b9b85e73SStefano Zampini   /* check user defined change of basis (if any) */
1144b9b85e73SStefano Zampini   if (pcbddc->user_ChangeOfBasisMatrix) {
1145b9b85e73SStefano Zampini     PC_IS* pcis= (PC_IS*)pc->data;
1146b9b85e73SStefano Zampini     PetscInt n,m;
1147b9b85e73SStefano Zampini     ierr = MatGetSize(pcbddc->user_ChangeOfBasisMatrix,&n,&m);CHKERRQ(ierr);
1148b9b85e73SStefano Zampini     if (n != m) {
1149b9b85e73SStefano Zampini       SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Change of basis matrix should be square %d != %d\n",n,m);
1150b9b85e73SStefano Zampini     } else if (n != pcis->n_B && n != pcis->n) {
1151b9b85e73SStefano Zampini       SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Size of change of basis matrix %d differs from boundary size %d and local size %d\n",n,pcis->n_B,pcis->n);
1152b9b85e73SStefano Zampini     }
1153b9b85e73SStefano Zampini     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1154b9b85e73SStefano Zampini     pcbddc->use_change_of_basis = PETSC_FALSE;
1155b9b85e73SStefano Zampini   }
1156f4ddd8eeSStefano Zampini   /* Analyze interface */
1157674ae819SStefano Zampini   if (computetopography) {
1158674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1159674ae819SStefano Zampini   }
1160fb8d54d4SStefano Zampini 
1161f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1162fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1163f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1164f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1165f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1166f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1167f4ddd8eeSStefano Zampini     } else {
1168f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1169f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1170f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1171165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1172f4ddd8eeSStefano Zampini         PetscInt         i;
1173165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1174165b64e2SStefano Zampini         PetscObjectState state;
1175165b64e2SStefano Zampini         PetscInt         nnsp_size;
1176165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1177f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1178f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1179165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1180f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1181f4ddd8eeSStefano Zampini             break;
1182f4ddd8eeSStefano Zampini           }
1183f4ddd8eeSStefano Zampini         }
1184f4ddd8eeSStefano Zampini       }
1185f4ddd8eeSStefano Zampini     }
1186f4ddd8eeSStefano Zampini   } else {
1187f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1188f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1189f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1190f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1191f4ddd8eeSStefano Zampini     }
1192f4ddd8eeSStefano Zampini   }
1193f4ddd8eeSStefano Zampini 
1194f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1195727cdba6SStefano Zampini   /* reset primal space flags */
1196f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1197727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
1198fb8d54d4SStefano Zampini   if (computetopography || new_nearnullspace_provided) {
1199727cdba6SStefano Zampini     /* It also sets the primal space flags */
1200674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1201e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1202f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
1203674ae819SStefano Zampini   }
1204fb8d54d4SStefano Zampini 
1205f4ddd8eeSStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
1206674ae819SStefano Zampini     /* reset data */
1207674ae819SStefano Zampini     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
120899cc7994SStefano Zampini     /* Create coarse and local stuffs */
120999cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1210674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
12110c7d97c5SJed Brown   }
121270cf5478SStefano Zampini 
121358a03d70SStefano Zampini   if (pcbddc->dbg_flag) {
121458a03d70SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
12152b510759SStefano Zampini   }
12160c7d97c5SJed Brown   PetscFunctionReturn(0);
12170c7d97c5SJed Brown }
12180c7d97c5SJed Brown 
12190c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
12200c7d97c5SJed Brown /*
122150efa1b5SStefano Zampini    PCApply_BDDC - Applies the BDDC operator to a vector.
12220c7d97c5SJed Brown 
12230c7d97c5SJed Brown    Input Parameters:
12240c7d97c5SJed Brown .  pc - the preconditioner context
12250c7d97c5SJed Brown .  r - input vector (global)
12260c7d97c5SJed Brown 
12270c7d97c5SJed Brown    Output Parameter:
12280c7d97c5SJed Brown .  z - output vector (global)
12290c7d97c5SJed Brown 
12300c7d97c5SJed Brown    Application Interface Routine: PCApply()
12310c7d97c5SJed Brown  */
12320c7d97c5SJed Brown #undef __FUNCT__
12330c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
123453cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
12350c7d97c5SJed Brown {
12360c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
12370c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
12380c7d97c5SJed Brown   PetscErrorCode    ierr;
12393b03a366Sstefano_zampini   const PetscScalar one = 1.0;
12403b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
12412617d88aSStefano Zampini   const PetscScalar zero = 0.0;
12420c7d97c5SJed Brown 
12430c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
12440c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
12458eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
12460c7d97c5SJed Brown 
12470c7d97c5SJed Brown   PetscFunctionBegin;
124885c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
12490c7d97c5SJed Brown     /* First Dirichlet solve */
12500c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12510c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
125253cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
12530c7d97c5SJed Brown     /*
12540c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1255674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1256674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
12570c7d97c5SJed Brown     */
12580c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
12590c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
12608eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
12610c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
12620c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
12630c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12640c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1265674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1266b76ba322SStefano Zampini   } else {
12670bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1268b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1269674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1270b76ba322SStefano Zampini   }
1271b76ba322SStefano Zampini 
12722617d88aSStefano Zampini   /* Apply interface preconditioner
12732617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1274dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
12752617d88aSStefano Zampini 
1276674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1277674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
12780c7d97c5SJed Brown 
12793b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
12800c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12810c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12820c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
12838eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1284df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1285df187020SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1286df187020SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1287df187020SStefano Zampini   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
12880c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12890c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12900c7d97c5SJed Brown   PetscFunctionReturn(0);
12910c7d97c5SJed Brown }
129250efa1b5SStefano Zampini 
129350efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */
129450efa1b5SStefano Zampini /*
129550efa1b5SStefano Zampini    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
129650efa1b5SStefano Zampini 
129750efa1b5SStefano Zampini    Input Parameters:
129850efa1b5SStefano Zampini .  pc - the preconditioner context
129950efa1b5SStefano Zampini .  r - input vector (global)
130050efa1b5SStefano Zampini 
130150efa1b5SStefano Zampini    Output Parameter:
130250efa1b5SStefano Zampini .  z - output vector (global)
130350efa1b5SStefano Zampini 
130450efa1b5SStefano Zampini    Application Interface Routine: PCApplyTranspose()
130550efa1b5SStefano Zampini  */
130650efa1b5SStefano Zampini #undef __FUNCT__
130750efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC"
130850efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
130950efa1b5SStefano Zampini {
131050efa1b5SStefano Zampini   PC_IS             *pcis = (PC_IS*)(pc->data);
131150efa1b5SStefano Zampini   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
131250efa1b5SStefano Zampini   PetscErrorCode    ierr;
131350efa1b5SStefano Zampini   const PetscScalar one = 1.0;
131450efa1b5SStefano Zampini   const PetscScalar m_one = -1.0;
131550efa1b5SStefano Zampini   const PetscScalar zero = 0.0;
131650efa1b5SStefano Zampini 
131750efa1b5SStefano Zampini   PetscFunctionBegin;
131850efa1b5SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
131950efa1b5SStefano Zampini     /* First Dirichlet solve */
132050efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
132150efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
132250efa1b5SStefano Zampini     ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
132350efa1b5SStefano Zampini     /*
132450efa1b5SStefano Zampini       Assembling right hand side for BDDC operator
132550efa1b5SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
132650efa1b5SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
132750efa1b5SStefano Zampini     */
132850efa1b5SStefano Zampini     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
132950efa1b5SStefano Zampini     ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
133050efa1b5SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
133150efa1b5SStefano Zampini     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
133250efa1b5SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
133350efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
133450efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
133550efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
133650efa1b5SStefano Zampini   } else {
133750efa1b5SStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
133850efa1b5SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
133950efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
134050efa1b5SStefano Zampini   }
134150efa1b5SStefano Zampini 
134250efa1b5SStefano Zampini   /* Apply interface preconditioner
134350efa1b5SStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1344dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
134550efa1b5SStefano Zampini 
134650efa1b5SStefano Zampini   /* Apply transpose of partition of unity operator */
134750efa1b5SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
134850efa1b5SStefano Zampini 
134950efa1b5SStefano Zampini   /* Second Dirichlet solve and assembling of output */
135050efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
135150efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
135250efa1b5SStefano Zampini   ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
135350efa1b5SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1354b0147a47SStefano Zampini   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1355b0147a47SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1356b0147a47SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1357b0147a47SStefano Zampini   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
135850efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
135950efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
136050efa1b5SStefano Zampini   PetscFunctionReturn(0);
136150efa1b5SStefano Zampini }
1362da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1363674ae819SStefano Zampini 
1364da1bb401SStefano Zampini #undef __FUNCT__
1365da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1366da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1367da1bb401SStefano Zampini {
1368da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1369da1bb401SStefano Zampini   PetscErrorCode ierr;
1370da1bb401SStefano Zampini 
1371da1bb401SStefano Zampini   PetscFunctionBegin;
1372da1bb401SStefano Zampini   /* free data created by PCIS */
1373da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1374674ae819SStefano Zampini   /* free BDDC custom data  */
1375674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1376674ae819SStefano Zampini   /* destroy objects related to topography */
1377674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1378674ae819SStefano Zampini   /* free allocated graph structure */
1379da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1380674ae819SStefano Zampini   /* free data for scaling operator */
1381674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1382674ae819SStefano Zampini   /* free solvers stuff */
1383674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
138462a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
138562a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
138662a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
138739e2fb2aSStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr);
13883425bc38SStefano Zampini   /* remove functions */
1389b9b85e73SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisLocalMat_C",NULL);CHKERRQ(ierr);
1390674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1391bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
13922b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1393b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
13942b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1395bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1396bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
139782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1398bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
139982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1400bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
140182ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1402bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1403785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1404bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
140563602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1406bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1407bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1408bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1409bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1410674ae819SStefano Zampini   /* Free the private data structure */
1411674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1412da1bb401SStefano Zampini   PetscFunctionReturn(0);
1413da1bb401SStefano Zampini }
14143425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
14151e6b0712SBarry Smith 
14163425bc38SStefano Zampini #undef __FUNCT__
14173425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
14183425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
14193425bc38SStefano Zampini {
1420674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
14213425bc38SStefano Zampini   PC_IS*         pcis;
14223425bc38SStefano Zampini   PC_BDDC*       pcbddc;
14233425bc38SStefano Zampini   PetscErrorCode ierr;
14240c7d97c5SJed Brown 
14253425bc38SStefano Zampini   PetscFunctionBegin;
14263425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
14273425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
14283425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
14293425bc38SStefano Zampini 
14303425bc38SStefano Zampini   /* change of basis for physical rhs if needed
14313425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
14323308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
14333425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
14343425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14353425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1436fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1437fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
14383425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14393425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1440674ae819SStefano Zampini   /* Apply partition of unity */
14413425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1442674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
14438eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
14443425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
14453425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
14463425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
14473425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
14483425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
14493425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14503425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1451674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
14523425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14533425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14543425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
14553425bc38SStefano Zampini   }
14563425bc38SStefano Zampini   /* BDDC rhs */
14573425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
14588eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
14593425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
14603425bc38SStefano Zampini   }
14613425bc38SStefano Zampini   /* apply BDDC */
1462dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
14633425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
14643425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
14653425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
14663425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14673425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14683425bc38SStefano Zampini   /* restore original rhs */
14693425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
14703425bc38SStefano Zampini   PetscFunctionReturn(0);
14713425bc38SStefano Zampini }
14721e6b0712SBarry Smith 
14733425bc38SStefano Zampini #undef __FUNCT__
14743425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
14753425bc38SStefano Zampini /*@
147628509bceSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
14773425bc38SStefano Zampini 
14783425bc38SStefano Zampini    Collective
14793425bc38SStefano Zampini 
14803425bc38SStefano Zampini    Input Parameters:
148128509bceSStefano Zampini +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
148228509bceSStefano Zampini .  standard_rhs - the right-hand side for your linear system
14833425bc38SStefano Zampini 
14843425bc38SStefano Zampini    Output Parameters:
148528509bceSStefano Zampini -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
14863425bc38SStefano Zampini 
14873425bc38SStefano Zampini    Level: developer
14883425bc38SStefano Zampini 
14893425bc38SStefano Zampini    Notes:
14903425bc38SStefano Zampini 
149128509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
14923425bc38SStefano Zampini @*/
14933425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
14943425bc38SStefano Zampini {
1495674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
14963425bc38SStefano Zampini   PetscErrorCode ierr;
14973425bc38SStefano Zampini 
14983425bc38SStefano Zampini   PetscFunctionBegin;
14993425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
15003425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
15013425bc38SStefano Zampini   PetscFunctionReturn(0);
15023425bc38SStefano Zampini }
15033425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
15041e6b0712SBarry Smith 
15053425bc38SStefano Zampini #undef __FUNCT__
15063425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
15073425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
15083425bc38SStefano Zampini {
1509674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
15103425bc38SStefano Zampini   PC_IS*         pcis;
15113425bc38SStefano Zampini   PC_BDDC*       pcbddc;
15123425bc38SStefano Zampini   PetscErrorCode ierr;
15133425bc38SStefano Zampini 
15143425bc38SStefano Zampini   PetscFunctionBegin;
15153425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
15163425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
15173425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
15183425bc38SStefano Zampini 
15193425bc38SStefano Zampini   /* apply B_delta^T */
15203425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15213425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15223425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
15233425bc38SStefano Zampini   /* compute rhs for BDDC application */
15243425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
15258eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
15263425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
15273425bc38SStefano Zampini   }
15283425bc38SStefano Zampini   /* apply BDDC */
1529dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
15303425bc38SStefano Zampini   /* put values into standard global vector */
15313425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15323425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15338eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
15343425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
15353425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
15363425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
15373425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
15383425bc38SStefano Zampini   }
15393425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15403425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15413425bc38SStefano Zampini   /* final change of basis if needed
15423425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
15433308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
15443425bc38SStefano Zampini   PetscFunctionReturn(0);
15453425bc38SStefano Zampini }
15461e6b0712SBarry Smith 
15473425bc38SStefano Zampini #undef __FUNCT__
15483425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
15493425bc38SStefano Zampini /*@
155028509bceSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
15513425bc38SStefano Zampini 
15523425bc38SStefano Zampini    Collective
15533425bc38SStefano Zampini 
15543425bc38SStefano Zampini    Input Parameters:
155528509bceSStefano Zampini +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
155628509bceSStefano Zampini .  fetidp_flux_sol - the solution of the FETIDP linear system
15573425bc38SStefano Zampini 
15583425bc38SStefano Zampini    Output Parameters:
155928509bceSStefano Zampini -  standard_sol      - the solution defined on the physical domain
15603425bc38SStefano Zampini 
15613425bc38SStefano Zampini    Level: developer
15623425bc38SStefano Zampini 
15633425bc38SStefano Zampini    Notes:
15643425bc38SStefano Zampini 
156528509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
15663425bc38SStefano Zampini @*/
15673425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
15683425bc38SStefano Zampini {
1569674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
15703425bc38SStefano Zampini   PetscErrorCode ierr;
15713425bc38SStefano Zampini 
15723425bc38SStefano Zampini   PetscFunctionBegin;
15733425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
15743425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
15753425bc38SStefano Zampini   PetscFunctionReturn(0);
15763425bc38SStefano Zampini }
15773425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
15781e6b0712SBarry Smith 
1579f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1580edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1581f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1582f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1583edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1584f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1585674ae819SStefano Zampini 
15863425bc38SStefano Zampini #undef __FUNCT__
15873425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
15883425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
15893425bc38SStefano Zampini {
1590674ae819SStefano Zampini 
1591674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
15923425bc38SStefano Zampini   Mat            newmat;
1593674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
15943425bc38SStefano Zampini   PC             newpc;
1595ce94432eSBarry Smith   MPI_Comm       comm;
15963425bc38SStefano Zampini   PetscErrorCode ierr;
15973425bc38SStefano Zampini 
15983425bc38SStefano Zampini   PetscFunctionBegin;
1599ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
16003425bc38SStefano Zampini   /* FETIDP linear matrix */
16013425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
16023425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
16033425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
16043425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1605edf7251bSStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
16063425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
16073425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
16083425bc38SStefano Zampini   /* FETIDP preconditioner */
16093425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
16103425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
16113425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
16123425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
16133425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
16143425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1615edf7251bSStefano Zampini   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
16163425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
161723ee1639SBarry Smith   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
16183425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
16193425bc38SStefano Zampini   /* return pointers for objects created */
16203425bc38SStefano Zampini   *fetidp_mat=newmat;
16213425bc38SStefano Zampini   *fetidp_pc=newpc;
16223425bc38SStefano Zampini   PetscFunctionReturn(0);
16233425bc38SStefano Zampini }
16241e6b0712SBarry Smith 
16253425bc38SStefano Zampini #undef __FUNCT__
16263425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
16273425bc38SStefano Zampini /*@
162828509bceSStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
16293425bc38SStefano Zampini 
16303425bc38SStefano Zampini    Collective
16313425bc38SStefano Zampini 
16323425bc38SStefano Zampini    Input Parameters:
163328509bceSStefano Zampini +  pc - the BDDC preconditioning context already setup
163428509bceSStefano Zampini 
163528509bceSStefano Zampini    Output Parameters:
163628509bceSStefano Zampini .  fetidp_mat - shell FETIDP matrix object
163728509bceSStefano Zampini .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
163828509bceSStefano Zampini 
163928509bceSStefano Zampini    Options Database Keys:
164028509bceSStefano Zampini -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
16413425bc38SStefano Zampini 
16423425bc38SStefano Zampini    Level: developer
16433425bc38SStefano Zampini 
16443425bc38SStefano Zampini    Notes:
164528509bceSStefano Zampini      Currently the only operation provided for FETIDP matrix is MatMult
16463425bc38SStefano Zampini 
16473425bc38SStefano Zampini .seealso: PCBDDC
16483425bc38SStefano Zampini @*/
16493425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
16503425bc38SStefano Zampini {
16513425bc38SStefano Zampini   PetscErrorCode ierr;
16523425bc38SStefano Zampini 
16533425bc38SStefano Zampini   PetscFunctionBegin;
16543425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
16553425bc38SStefano Zampini   if (pc->setupcalled) {
1656516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1657f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
16583425bc38SStefano Zampini   PetscFunctionReturn(0);
16593425bc38SStefano Zampini }
16600c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1661da1bb401SStefano Zampini /*MC
1662da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
16630c7d97c5SJed Brown 
166428509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
166528509bceSStefano Zampini 
166628509bceSStefano Zampini .vb
166728509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
166828509bceSStefano 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
166928509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
167028509bceSStefano Zampini .ve
167128509bceSStefano Zampini 
167228509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
167328509bceSStefano Zampini 
1674b6fdb6dfSStefano Zampini    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
167528509bceSStefano Zampini 
167628509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
167728509bceSStefano Zampini 
1678b6fdb6dfSStefano 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.
1679b6fdb6dfSStefano Zampini 
168028509bceSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
168128509bceSStefano Zampini 
16826a818285SBarry Smith    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()
168328509bceSStefano Zampini 
16846a818285SBarry Smith    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace().
168528509bceSStefano Zampini 
1686b6fdb6dfSStefano 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.
168728509bceSStefano Zampini 
168828509bceSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
168928509bceSStefano Zampini 
1690da1bb401SStefano Zampini    Options Database Keys:
169128509bceSStefano Zampini 
169228509bceSStefano Zampini .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
169328509bceSStefano Zampini .    -pc_bddc_use_edges <1> - use or not edges in primal space
1694b6fdb6dfSStefano Zampini .    -pc_bddc_use_faces <0> - use or not faces in primal space
169528509bceSStefano Zampini .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
169628509bceSStefano Zampini .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
169728509bceSStefano Zampini .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
169828509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
169928509bceSStefano Zampini .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
170028509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
170128509bceSStefano Zampini 
170228509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
170328509bceSStefano Zampini .vb
170428509bceSStefano Zampini       -pc_bddc_dirichlet_
170528509bceSStefano Zampini       -pc_bddc_neumann_
170628509bceSStefano Zampini       -pc_bddc_coarse_
170728509bceSStefano Zampini .ve
170828509bceSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
170928509bceSStefano Zampini 
171028509bceSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level can be specified as
171128509bceSStefano Zampini .vb
1712312be037SStefano Zampini       -pc_bddc_dirichlet_lN_
1713312be037SStefano Zampini       -pc_bddc_neumann_lN_
1714312be037SStefano Zampini       -pc_bddc_coarse_lN_
171528509bceSStefano Zampini .ve
1716312be037SStefano Zampini    Note that level number ranges from the finest 0 to the coarsest N.
1717da1bb401SStefano Zampini 
1718da1bb401SStefano Zampini    Level: intermediate
1719da1bb401SStefano Zampini 
1720b6fdb6dfSStefano Zampini    Developer notes:
1721da1bb401SStefano Zampini 
172228509bceSStefano Zampini    New deluxe scaling operator will be available soon.
1723da1bb401SStefano Zampini 
1724da1bb401SStefano Zampini    Contributed by Stefano Zampini
1725da1bb401SStefano Zampini 
1726da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1727da1bb401SStefano Zampini M*/
1728b2573a8aSBarry Smith 
1729da1bb401SStefano Zampini #undef __FUNCT__
1730da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
17318cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1732da1bb401SStefano Zampini {
1733da1bb401SStefano Zampini   PetscErrorCode      ierr;
1734da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1735da1bb401SStefano Zampini 
1736da1bb401SStefano Zampini   PetscFunctionBegin;
1737da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1738b00a9115SJed Brown   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
1739da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1740da1bb401SStefano Zampini 
1741da1bb401SStefano Zampini   /* create PCIS data structure */
1742da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1743da1bb401SStefano Zampini 
174447d04d0dSStefano Zampini   /* BDDC customization */
174508a5cf49SStefano Zampini   pcbddc->use_local_adj       = PETSC_TRUE;
174647d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
174747d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
174847d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
174947d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
175047d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
175147d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
175247d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
175347d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
1754b9d89cd5SStefano Zampini   /* private */
1755b9d89cd5SStefano Zampini   pcbddc->issym                      = PETSC_FALSE;
175639e2fb2aSStefano Zampini   pcbddc->BtoNmap                    = 0;
1757727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
1758e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
1759e9189074SStefano Zampini   pcbddc->n_actual_vertices          = 0;
1760e9189074SStefano Zampini   pcbddc->n_constraints              = 0;
1761727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
1762fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
176368457ee5SStefano Zampini   pcbddc->coarse_size                = -1;
1764f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
1765727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
1766f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
1767f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
1768f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
1769674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
17700bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
17713972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1772534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1773534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1774534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1775b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix   = 0;
1776da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1777da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1778da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1779da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1780da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
178115aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
178215aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1783da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1784da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1785da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1786da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1787da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1788da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1789da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1790da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1791da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1792da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
1793785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
1794785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
1795785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
179660d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
179760d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
179863602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal           = 0;
1799da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
180063602bcaSStefano Zampini   pcbddc->ISForDofsLocal             = 0;
1801da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
180285c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
180347d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
180447d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
18054fad6a16SStefano Zampini   pcbddc->current_level              = 0;
18062b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1807323d395dSStefano Zampini   pcbddc->use_coarse_estimates       = PETSC_FALSE;
1808f3bde8b3SStefano Zampini   pcbddc->coarse_subassembling       = 0;
1809323d395dSStefano Zampini   pcbddc->coarse_subassembling_init  = 0;
1810da1bb401SStefano Zampini 
1811674ae819SStefano Zampini   /* create local graph structure */
1812674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1813674ae819SStefano Zampini 
1814674ae819SStefano Zampini   /* scaling */
1815674ae819SStefano Zampini   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
181608a5cf49SStefano Zampini   pcbddc->deluxe_ctx                 = 0;
1817674ae819SStefano Zampini   pcbddc->work_scaling               = 0;
1818da1bb401SStefano Zampini 
1819da1bb401SStefano Zampini   /* function pointers */
1820da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
182193bd9ae7SStefano Zampini   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
1822da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1823da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1824da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1825da1bb401SStefano Zampini   pc->ops->view                = 0;
1826da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1827da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1828da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1829534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1830534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1831da1bb401SStefano Zampini 
1832da1bb401SStefano Zampini   /* composing function */
1833b9b85e73SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisLocalMat_C",PCBDDCSetChangeOfBasisLocalMat_BDDC);CHKERRQ(ierr);
1834674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1835bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
18362b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1837b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
18382b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1839bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1840bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
184182ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1842bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
184382ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1844bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
184582ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1846bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
184782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1848bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
184963602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
1850bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1851bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1852bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1853bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1854da1bb401SStefano Zampini   PetscFunctionReturn(0);
1855da1bb401SStefano Zampini }
18563425bc38SStefano Zampini 
1857