xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision b6fdb6df294eae8bbbbd30da17bc256f98e2e680)
153cdbc3dSStefano Zampini /* TODOLIST
260d44989SStefano Zampini    Better management for block size > 1 for idx_R_local and others
385c4d303SStefano Zampini    Provide PCApplyTranpose
4af732b37SStefano Zampini    make runexe59
5be4bc593SStefano Zampini    Man pages
628509bceSStefano Zampini    Propagate nearnullspace info among levels
760d44989SStefano Zampini    Change of basis approach does not work with my nonlinear mechanics example. why? maybe an issue with l2gmap?
8be4bc593SStefano Zampini    Move FETIDP code
9be4bc593SStefano Zampini    Provide general case for subassembling
10be4bc593SStefano Zampini    Preallocation routines in MatConvert_IS_AIJ
1128509bceSStefano Zampini    Allow different customizations among different linear solves (requires also reset/destroy of ksp_R and coarse_ksp)
12be4bc593SStefano Zampini    Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver?
13be4bc593SStefano Zampini    Better management in PCIS code
1485c4d303SStefano Zampini    Is it possible working with PCBDDCGraph on boundary indices only?
15da1bb401SStefano Zampini    DofSplitting and DM attached to pc?
16da1bb401SStefano Zampini    Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet)
17b8ffe317SStefano Zampini    BDDC with MG framework?
1853cdbc3dSStefano Zampini */
190c7d97c5SJed Brown 
2053cdbc3dSStefano Zampini /* ----------------------------------------------------------------------------------------------------------------------------------------------
210c7d97c5SJed Brown    Implementation of BDDC preconditioner based on:
220c7d97c5SJed Brown    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
2353cdbc3dSStefano Zampini    ---------------------------------------------------------------------------------------------------------------------------------------------- */
2453cdbc3dSStefano Zampini 
25674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/  /* includes for fortran wrappers */
26674ae819SStefano Zampini #include "bddcprivate.h"
273b03a366Sstefano_zampini #include <petscblaslapack.h>
28674ae819SStefano Zampini 
290c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
300c7d97c5SJed Brown #undef __FUNCT__
310c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
320c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc)
330c7d97c5SJed Brown {
340c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
350c7d97c5SJed Brown   PetscErrorCode ierr;
360c7d97c5SJed Brown 
370c7d97c5SJed Brown   PetscFunctionBegin;
380c7d97c5SJed Brown   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
398eeda7d8SStefano Zampini   /* Verbose debugging */
408eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
418eeda7d8SStefano Zampini   /* Primal space cumstomization */
428eeda7d8SStefano 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);
438eeda7d8SStefano 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);
448eeda7d8SStefano 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);
458eeda7d8SStefano Zampini   /* Change of basis */
468eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr);
478eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr);
48674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
49674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
50674ae819SStefano Zampini   }
518eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
528eeda7d8SStefano 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);
530298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
542b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
55674ae819SStefano 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);
560c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
570c7d97c5SJed Brown   PetscFunctionReturn(0);
580c7d97c5SJed Brown }
590c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
60674ae819SStefano Zampini #undef __FUNCT__
61674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
62674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
63674ae819SStefano Zampini {
64674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
65674ae819SStefano Zampini   PetscErrorCode ierr;
661e6b0712SBarry Smith 
67674ae819SStefano Zampini   PetscFunctionBegin;
68674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
69674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
70674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
71674ae819SStefano Zampini   PetscFunctionReturn(0);
72674ae819SStefano Zampini }
73674ae819SStefano Zampini #undef __FUNCT__
74674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
75674ae819SStefano Zampini /*@
7628509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
77674ae819SStefano Zampini 
78674ae819SStefano Zampini    Not collective
79674ae819SStefano Zampini 
80674ae819SStefano Zampini    Input Parameters:
81674ae819SStefano Zampini +  pc - the preconditioning context
8228509bceSStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering
83674ae819SStefano Zampini 
84674ae819SStefano Zampini    Level: intermediate
85674ae819SStefano Zampini 
86674ae819SStefano Zampini    Notes:
87674ae819SStefano Zampini 
88674ae819SStefano Zampini .seealso: PCBDDC
89674ae819SStefano Zampini @*/
90674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
91674ae819SStefano Zampini {
92674ae819SStefano Zampini   PetscErrorCode ierr;
93674ae819SStefano Zampini 
94674ae819SStefano Zampini   PetscFunctionBegin;
95674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
96674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
97674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
98674ae819SStefano Zampini   PetscFunctionReturn(0);
99674ae819SStefano Zampini }
100674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1010c7d97c5SJed Brown #undef __FUNCT__
1024fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1034fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1044fad6a16SStefano Zampini {
1054fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1064fad6a16SStefano Zampini 
1074fad6a16SStefano Zampini   PetscFunctionBegin;
1084fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1094fad6a16SStefano Zampini   PetscFunctionReturn(0);
1104fad6a16SStefano Zampini }
1111e6b0712SBarry Smith 
1124fad6a16SStefano Zampini #undef __FUNCT__
1134fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1144fad6a16SStefano Zampini /*@
11528509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
1164fad6a16SStefano Zampini 
1174fad6a16SStefano Zampini    Logically collective on PC
1184fad6a16SStefano Zampini 
1194fad6a16SStefano Zampini    Input Parameters:
1204fad6a16SStefano Zampini +  pc - the preconditioning context
12128509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
1224fad6a16SStefano Zampini 
12328509bceSStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
1244fad6a16SStefano Zampini 
1254fad6a16SStefano Zampini    Level: intermediate
1264fad6a16SStefano Zampini 
1274fad6a16SStefano Zampini    Notes:
1284fad6a16SStefano Zampini 
1294fad6a16SStefano Zampini .seealso: PCBDDC
1304fad6a16SStefano Zampini @*/
1314fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
1324fad6a16SStefano Zampini {
1334fad6a16SStefano Zampini   PetscErrorCode ierr;
1344fad6a16SStefano Zampini 
1354fad6a16SStefano Zampini   PetscFunctionBegin;
1364fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1372b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
1384fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
1394fad6a16SStefano Zampini   PetscFunctionReturn(0);
1404fad6a16SStefano Zampini }
1412b510759SStefano Zampini 
142b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
1432b510759SStefano Zampini #undef __FUNCT__
144b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
145b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
146b8ffe317SStefano Zampini {
147b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
148b8ffe317SStefano Zampini 
149b8ffe317SStefano Zampini   PetscFunctionBegin;
15085c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
151b8ffe317SStefano Zampini   PetscFunctionReturn(0);
152b8ffe317SStefano Zampini }
153b8ffe317SStefano Zampini 
154b8ffe317SStefano Zampini #undef __FUNCT__
155b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
156b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
1572b510759SStefano Zampini {
1582b510759SStefano Zampini   PetscErrorCode ierr;
1592b510759SStefano Zampini 
1602b510759SStefano Zampini   PetscFunctionBegin;
1612b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
162b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
163b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1642b510759SStefano Zampini   PetscFunctionReturn(0);
1652b510759SStefano Zampini }
1661e6b0712SBarry Smith 
1674fad6a16SStefano Zampini #undef __FUNCT__
1682b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
1692b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
1704fad6a16SStefano Zampini {
1714fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1724fad6a16SStefano Zampini 
1734fad6a16SStefano Zampini   PetscFunctionBegin;
1742b510759SStefano Zampini   pcbddc->current_level = level;
1754fad6a16SStefano Zampini   PetscFunctionReturn(0);
1764fad6a16SStefano Zampini }
1771e6b0712SBarry Smith 
1784fad6a16SStefano Zampini #undef __FUNCT__
179b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
180b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
181b8ffe317SStefano Zampini {
182b8ffe317SStefano Zampini   PetscErrorCode ierr;
183b8ffe317SStefano Zampini 
184b8ffe317SStefano Zampini   PetscFunctionBegin;
185b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
186b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
187b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
188b8ffe317SStefano Zampini   PetscFunctionReturn(0);
189b8ffe317SStefano Zampini }
190b8ffe317SStefano Zampini 
191b8ffe317SStefano Zampini #undef __FUNCT__
1922b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
1932b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
1942b510759SStefano Zampini {
1952b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1962b510759SStefano Zampini 
1972b510759SStefano Zampini   PetscFunctionBegin;
1982b510759SStefano Zampini   pcbddc->max_levels = levels;
1992b510759SStefano Zampini   PetscFunctionReturn(0);
2002b510759SStefano Zampini }
2012b510759SStefano Zampini 
2022b510759SStefano Zampini #undef __FUNCT__
2032b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2044fad6a16SStefano Zampini /*@
20528509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
2064fad6a16SStefano Zampini 
2074fad6a16SStefano Zampini    Logically collective on PC
2084fad6a16SStefano Zampini 
2094fad6a16SStefano Zampini    Input Parameters:
2104fad6a16SStefano Zampini +  pc - the preconditioning context
21128509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
2124fad6a16SStefano Zampini 
21328509bceSStefano Zampini    Default value is 0, i.e. traditional one-level BDDC
2144fad6a16SStefano Zampini 
2154fad6a16SStefano Zampini    Level: intermediate
2164fad6a16SStefano Zampini 
2174fad6a16SStefano Zampini    Notes:
2184fad6a16SStefano Zampini 
2194fad6a16SStefano Zampini .seealso: PCBDDC
2204fad6a16SStefano Zampini @*/
2212b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
2224fad6a16SStefano Zampini {
2234fad6a16SStefano Zampini   PetscErrorCode ierr;
2244fad6a16SStefano Zampini 
2254fad6a16SStefano Zampini   PetscFunctionBegin;
2264fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2272b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
2282b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
2294fad6a16SStefano Zampini   PetscFunctionReturn(0);
2304fad6a16SStefano Zampini }
2314fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
2321e6b0712SBarry Smith 
2334fad6a16SStefano Zampini #undef __FUNCT__
2340bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
2350bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
2360bdf917eSStefano Zampini {
2370bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2380bdf917eSStefano Zampini   PetscErrorCode ierr;
2390bdf917eSStefano Zampini 
2400bdf917eSStefano Zampini   PetscFunctionBegin;
2410bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
2420bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
2430bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
2440bdf917eSStefano Zampini   PetscFunctionReturn(0);
2450bdf917eSStefano Zampini }
2461e6b0712SBarry Smith 
2470bdf917eSStefano Zampini #undef __FUNCT__
2480bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
2490bdf917eSStefano Zampini /*@
25028509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
2510bdf917eSStefano Zampini 
2520bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
2530bdf917eSStefano Zampini 
2540bdf917eSStefano Zampini    Input Parameters:
2550bdf917eSStefano Zampini +  pc - the preconditioning context
25628509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
2570bdf917eSStefano Zampini 
2580bdf917eSStefano Zampini    Level: intermediate
2590bdf917eSStefano Zampini 
2600bdf917eSStefano Zampini    Notes:
2610bdf917eSStefano Zampini 
2620bdf917eSStefano Zampini .seealso: PCBDDC
2630bdf917eSStefano Zampini @*/
2640bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
2650bdf917eSStefano Zampini {
2660bdf917eSStefano Zampini   PetscErrorCode ierr;
2670bdf917eSStefano Zampini 
2680bdf917eSStefano Zampini   PetscFunctionBegin;
2690bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
270674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
2712b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
2720bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
2730bdf917eSStefano Zampini   PetscFunctionReturn(0);
2740bdf917eSStefano Zampini }
2750bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
2761e6b0712SBarry Smith 
2770bdf917eSStefano Zampini #undef __FUNCT__
2783b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
2793b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
2803b03a366Sstefano_zampini {
2813b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2823b03a366Sstefano_zampini   PetscErrorCode ierr;
2833b03a366Sstefano_zampini 
2843b03a366Sstefano_zampini   PetscFunctionBegin;
2853b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
28636e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
28736e030ebSStefano Zampini   pcbddc->DirichletBoundaries=DirichletBoundaries;
2883b03a366Sstefano_zampini   PetscFunctionReturn(0);
2893b03a366Sstefano_zampini }
2901e6b0712SBarry Smith 
2913b03a366Sstefano_zampini #undef __FUNCT__
2923b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
2933b03a366Sstefano_zampini /*@
29428509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
2953b03a366Sstefano_zampini 
2963b03a366Sstefano_zampini    Not collective
2973b03a366Sstefano_zampini 
2983b03a366Sstefano_zampini    Input Parameters:
2993b03a366Sstefano_zampini +  pc - the preconditioning context
30028509bceSStefano Zampini -  DirichletBoundaries - sequential IS defining the subdomain part of Dirichlet boundaries (in local ordering)
3013b03a366Sstefano_zampini 
3023b03a366Sstefano_zampini    Level: intermediate
3033b03a366Sstefano_zampini 
3043b03a366Sstefano_zampini    Notes:
3053b03a366Sstefano_zampini 
3063b03a366Sstefano_zampini .seealso: PCBDDC
3073b03a366Sstefano_zampini @*/
3083b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3093b03a366Sstefano_zampini {
3103b03a366Sstefano_zampini   PetscErrorCode ierr;
3113b03a366Sstefano_zampini 
3123b03a366Sstefano_zampini   PetscFunctionBegin;
3133b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
314674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
3153b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
3163b03a366Sstefano_zampini   PetscFunctionReturn(0);
3173b03a366Sstefano_zampini }
3183b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
3191e6b0712SBarry Smith 
3203b03a366Sstefano_zampini #undef __FUNCT__
3210c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
32253cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
3230c7d97c5SJed Brown {
3240c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
32553cdbc3dSStefano Zampini   PetscErrorCode ierr;
3260c7d97c5SJed Brown 
3270c7d97c5SJed Brown   PetscFunctionBegin;
32853cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
32936e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
33036e030ebSStefano Zampini   pcbddc->NeumannBoundaries=NeumannBoundaries;
3310c7d97c5SJed Brown   PetscFunctionReturn(0);
3320c7d97c5SJed Brown }
3331e6b0712SBarry Smith 
3340c7d97c5SJed Brown #undef __FUNCT__
3350c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
33657527edcSJed Brown /*@
33728509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
33857527edcSJed Brown 
3399c0446d6SStefano Zampini    Not collective
34057527edcSJed Brown 
34157527edcSJed Brown    Input Parameters:
34257527edcSJed Brown +  pc - the preconditioning context
34328509bceSStefano Zampini -  NeumannBoundaries - sequential IS defining the subdomain part of Neumann boundaries (in local ordering)
34457527edcSJed Brown 
34557527edcSJed Brown    Level: intermediate
34657527edcSJed Brown 
34757527edcSJed Brown    Notes:
34857527edcSJed Brown 
34957527edcSJed Brown .seealso: PCBDDC
35057527edcSJed Brown @*/
35153cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
3520c7d97c5SJed Brown {
3530c7d97c5SJed Brown   PetscErrorCode ierr;
3540c7d97c5SJed Brown 
3550c7d97c5SJed Brown   PetscFunctionBegin;
3560c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
357674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
35853cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
35953cdbc3dSStefano Zampini   PetscFunctionReturn(0);
36053cdbc3dSStefano Zampini }
36153cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
3621e6b0712SBarry Smith 
36353cdbc3dSStefano Zampini #undef __FUNCT__
364da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
365da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
366da1bb401SStefano Zampini {
367da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
368da1bb401SStefano Zampini 
369da1bb401SStefano Zampini   PetscFunctionBegin;
370da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
371da1bb401SStefano Zampini   PetscFunctionReturn(0);
372da1bb401SStefano Zampini }
3731e6b0712SBarry Smith 
374da1bb401SStefano Zampini #undef __FUNCT__
375da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
376da1bb401SStefano Zampini /*@
37728509bceSStefano Zampini  PCBDDCGetDirichletBoundaries - Get IS for local Dirichlet boundaries
378da1bb401SStefano Zampini 
379da1bb401SStefano Zampini    Not collective
380da1bb401SStefano Zampini 
381da1bb401SStefano Zampini    Input Parameters:
38228509bceSStefano Zampini .  pc - the preconditioning context
383da1bb401SStefano Zampini 
384da1bb401SStefano Zampini    Output Parameters:
38528509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
386da1bb401SStefano Zampini 
387da1bb401SStefano Zampini    Level: intermediate
388da1bb401SStefano Zampini 
389da1bb401SStefano Zampini    Notes:
390da1bb401SStefano Zampini 
391da1bb401SStefano Zampini .seealso: PCBDDC
392da1bb401SStefano Zampini @*/
393da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
394da1bb401SStefano Zampini {
395da1bb401SStefano Zampini   PetscErrorCode ierr;
396da1bb401SStefano Zampini 
397da1bb401SStefano Zampini   PetscFunctionBegin;
398da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
399da1bb401SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
400da1bb401SStefano Zampini   PetscFunctionReturn(0);
401da1bb401SStefano Zampini }
402da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
4031e6b0712SBarry Smith 
404da1bb401SStefano Zampini #undef __FUNCT__
40553cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
40653cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
40753cdbc3dSStefano Zampini {
40853cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
40953cdbc3dSStefano Zampini 
41053cdbc3dSStefano Zampini   PetscFunctionBegin;
41153cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
41253cdbc3dSStefano Zampini   PetscFunctionReturn(0);
41353cdbc3dSStefano Zampini }
4141e6b0712SBarry Smith 
41553cdbc3dSStefano Zampini #undef __FUNCT__
41653cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
41753cdbc3dSStefano Zampini /*@
41828509bceSStefano Zampini  PCBDDCGetNeumannBoundaries - Get IS for local Neumann boundaries
41953cdbc3dSStefano Zampini 
4209c0446d6SStefano Zampini    Not collective
42153cdbc3dSStefano Zampini 
42253cdbc3dSStefano Zampini    Input Parameters:
42328509bceSStefano Zampini .  pc - the preconditioning context
42453cdbc3dSStefano Zampini 
42553cdbc3dSStefano Zampini    Output Parameters:
42628509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
42753cdbc3dSStefano Zampini 
42853cdbc3dSStefano Zampini    Level: intermediate
42953cdbc3dSStefano Zampini 
43053cdbc3dSStefano Zampini    Notes:
43153cdbc3dSStefano Zampini 
43253cdbc3dSStefano Zampini .seealso: PCBDDC
43353cdbc3dSStefano Zampini @*/
43453cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
43553cdbc3dSStefano Zampini {
43653cdbc3dSStefano Zampini   PetscErrorCode ierr;
43753cdbc3dSStefano Zampini 
43853cdbc3dSStefano Zampini   PetscFunctionBegin;
43953cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
44053cdbc3dSStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
4410c7d97c5SJed Brown   PetscFunctionReturn(0);
4420c7d97c5SJed Brown }
44336e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
4441e6b0712SBarry Smith 
44536e030ebSStefano Zampini #undef __FUNCT__
446da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
4471a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
44836e030ebSStefano Zampini {
44936e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
450da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
451da1bb401SStefano Zampini   PetscErrorCode ierr;
45236e030ebSStefano Zampini 
45336e030ebSStefano Zampini   PetscFunctionBegin;
454674ae819SStefano Zampini   /* free old CSR */
455674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
456fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
457674ae819SStefano Zampini   /* get CSR into graph structure */
458da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
459674ae819SStefano Zampini     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
460674ae819SStefano Zampini     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
461674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
462674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
463da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
4641a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
4651a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
466674ae819SStefano Zampini   }
467575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
46836e030ebSStefano Zampini   PetscFunctionReturn(0);
46936e030ebSStefano Zampini }
4701e6b0712SBarry Smith 
47136e030ebSStefano Zampini #undef __FUNCT__
472da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
47336e030ebSStefano Zampini /*@
47428509bceSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
47536e030ebSStefano Zampini 
47636e030ebSStefano Zampini    Not collective
47736e030ebSStefano Zampini 
47836e030ebSStefano Zampini    Input Parameters:
47936e030ebSStefano Zampini +  pc - the preconditioning context
48028509bceSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
48128509bceSStefano Zampini .  xadj, adjncy - the CSR graph
48228509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
48336e030ebSStefano Zampini 
48436e030ebSStefano Zampini    Level: intermediate
48536e030ebSStefano Zampini 
48636e030ebSStefano Zampini    Notes:
48736e030ebSStefano Zampini 
48828509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
48936e030ebSStefano Zampini @*/
4901a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
49136e030ebSStefano Zampini {
492575ad6abSStefano Zampini   void (*f)(void) = 0;
49336e030ebSStefano Zampini   PetscErrorCode ierr;
49436e030ebSStefano Zampini 
49536e030ebSStefano Zampini   PetscFunctionBegin;
49636e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
497674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
498674ae819SStefano Zampini   PetscValidIntPointer(xadj,4);
499674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
500674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
501da1bb401SStefano Zampini   }
50236e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
503575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
504575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
505575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
506575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
507575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
50836e030ebSStefano Zampini   }
50936e030ebSStefano Zampini   PetscFunctionReturn(0);
51036e030ebSStefano Zampini }
5119c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
5121e6b0712SBarry Smith 
5139c0446d6SStefano Zampini #undef __FUNCT__
5149c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
5159c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
5169c0446d6SStefano Zampini {
5179c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5189c0446d6SStefano Zampini   PetscInt i;
5199c0446d6SStefano Zampini   PetscErrorCode ierr;
5209c0446d6SStefano Zampini 
5219c0446d6SStefano Zampini   PetscFunctionBegin;
522da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
5239c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
5249c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
5259c0446d6SStefano Zampini   }
526d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
527da1bb401SStefano Zampini   /* allocate space then set */
5289c0446d6SStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
5299c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
530da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
531da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
5329c0446d6SStefano Zampini   }
5339c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
53460d44989SStefano Zampini   pcbddc->user_provided_isfordofs = PETSC_TRUE;
5359c0446d6SStefano Zampini   PetscFunctionReturn(0);
5369c0446d6SStefano Zampini }
5371e6b0712SBarry Smith 
5389c0446d6SStefano Zampini #undef __FUNCT__
5399c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
5409c0446d6SStefano Zampini /*@
54128509bceSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the local Neumann matrix
5429c0446d6SStefano Zampini 
5439c0446d6SStefano Zampini    Not collective
5449c0446d6SStefano Zampini 
5459c0446d6SStefano Zampini    Input Parameters:
5469c0446d6SStefano Zampini +  pc - the preconditioning context
54728509bceSStefano Zampini -  n_is - number of index sets defining the fields
54828509bceSStefano Zampini .  ISForDofs - array of IS describing the fields
5499c0446d6SStefano Zampini 
5509c0446d6SStefano Zampini    Level: intermediate
5519c0446d6SStefano Zampini 
5529c0446d6SStefano Zampini    Notes:
5539c0446d6SStefano Zampini 
5549c0446d6SStefano Zampini .seealso: PCBDDC
5559c0446d6SStefano Zampini @*/
5569c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
5579c0446d6SStefano Zampini {
5582b510759SStefano Zampini   PetscInt       i;
5599c0446d6SStefano Zampini   PetscErrorCode ierr;
5609c0446d6SStefano Zampini 
5619c0446d6SStefano Zampini   PetscFunctionBegin;
5629c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5632b510759SStefano Zampini   for (i=0;i<n_is;i++) {
5642b510759SStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2);
5652b510759SStefano Zampini   }
5669c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
5679c0446d6SStefano Zampini   PetscFunctionReturn(0);
5689c0446d6SStefano Zampini }
569da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
570534831adSStefano Zampini #undef __FUNCT__
571534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
572534831adSStefano Zampini /* -------------------------------------------------------------------------- */
573534831adSStefano Zampini /*
574534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
575534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
5769c0446d6SStefano Zampini 
577534831adSStefano Zampini    Input Parameter:
578534831adSStefano Zampini +  pc - the preconditioner contex
579534831adSStefano Zampini 
580534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
581534831adSStefano Zampini 
582534831adSStefano Zampini    Notes:
583534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
584534831adSStefano Zampini    the user, but instead is called by KSPSolve().
585534831adSStefano Zampini */
586534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
587534831adSStefano Zampini {
588534831adSStefano Zampini   PetscErrorCode ierr;
589534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
590534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
591534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
592534831adSStefano Zampini   Mat            temp_mat;
5933972b0daSStefano Zampini   IS             dirIS;
5943972b0daSStefano Zampini   PetscInt       dirsize,i,*is_indices;
5953972b0daSStefano Zampini   PetscScalar    *array_x,*array_diagonal;
5963972b0daSStefano Zampini   Vec            used_vec;
59792e3dcfbSStefano Zampini   PetscBool      guess_nonzero,flg,bddc_has_dirichlet_boundaries;
598534831adSStefano Zampini 
599534831adSStefano Zampini   PetscFunctionBegin;
60085c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
60185c4d303SStefano Zampini   if (ksp) {
60285c4d303SStefano Zampini     PetscBool iscg;
60385c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
60485c4d303SStefano Zampini     if (!iscg) {
60585c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
60685c4d303SStefano Zampini     }
60785c4d303SStefano Zampini   }
60885c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
60962a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
61062a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
61162a6ff1dSStefano Zampini   }
61262a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
61362a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
61462a6ff1dSStefano Zampini   }
6153972b0daSStefano Zampini   if (x) {
6163972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
6173972b0daSStefano Zampini     used_vec = x;
6183972b0daSStefano Zampini   } else {
6193972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
6203972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
6213972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6223972b0daSStefano Zampini   }
6233972b0daSStefano Zampini   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
6243972b0daSStefano Zampini   if (ksp) {
6253972b0daSStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
6263972b0daSStefano Zampini     if (!guess_nonzero) {
6273972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6283972b0daSStefano Zampini     }
6293972b0daSStefano Zampini   }
6303308cffdSStefano Zampini 
63192e3dcfbSStefano Zampini   /* TODO: remove when Dirichlet boundaries will be shared */
63292e3dcfbSStefano Zampini   ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr);
63392e3dcfbSStefano Zampini   flg = PETSC_FALSE;
63492e3dcfbSStefano Zampini   if (dirIS) flg = PETSC_TRUE;
63592e3dcfbSStefano Zampini   ierr = MPI_Allreduce(&flg,&bddc_has_dirichlet_boundaries,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
63692e3dcfbSStefano Zampini 
6373972b0daSStefano Zampini   /* store the original rhs */
6383972b0daSStefano Zampini   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
6393972b0daSStefano Zampini 
6403972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
64192e3dcfbSStefano Zampini   if (rhs && bddc_has_dirichlet_boundaries) {
6423972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
6433972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
6443972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6453972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6463972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6473972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6483972b0daSStefano Zampini     if (dirIS) {
6493972b0daSStefano Zampini       ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr);
6503972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6513972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6523972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6532fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
6543972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6553972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6563972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6573972b0daSStefano Zampini     }
6583972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6593972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
660b76ba322SStefano Zampini 
6613972b0daSStefano Zampini     /* remove the computed solution from the rhs */
6623972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6633972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
6643972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6653308cffdSStefano Zampini   }
666b76ba322SStefano Zampini 
667b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
6683972b0daSStefano Zampini   if (x) {
6693972b0daSStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
6703972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
67185c4d303SStefano Zampini     if (pcbddc->use_exact_dirichlet_trick) {
672b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
673b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
674b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
675b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
676b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
677b76ba322SStefano Zampini       if (ksp) {
678b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
679b76ba322SStefano Zampini       }
680b76ba322SStefano Zampini     }
6813972b0daSStefano Zampini   }
682b76ba322SStefano Zampini 
68392e3dcfbSStefano Zampini   /* prepare MatMult and rhs for solver */
684674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
685b76ba322SStefano Zampini     /* swap pointers for local matrices */
686b76ba322SStefano Zampini     temp_mat = matis->A;
687b76ba322SStefano Zampini     matis->A = pcbddc->local_mat;
688b76ba322SStefano Zampini     pcbddc->local_mat = temp_mat;
68992e3dcfbSStefano Zampini     if (rhs) {
690b76ba322SStefano Zampini       /* Get local rhs and apply transformation of basis */
691b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
692b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
693b76ba322SStefano Zampini       /* from original basis to modified basis */
694b76ba322SStefano Zampini       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
695b76ba322SStefano Zampini       /* put back modified values into the global vec using INSERT_VALUES copy mode */
696b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
697b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
698674ae819SStefano Zampini     }
69992e3dcfbSStefano Zampini   }
70092e3dcfbSStefano Zampini 
70192e3dcfbSStefano Zampini   /* remove nullspace if present */
7020bdf917eSStefano Zampini   if (ksp && pcbddc->NullSpace) {
703d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
704d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
705b76ba322SStefano Zampini   }
7060bdf917eSStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
707534831adSStefano Zampini   PetscFunctionReturn(0);
708534831adSStefano Zampini }
709534831adSStefano Zampini /* -------------------------------------------------------------------------- */
710534831adSStefano Zampini #undef __FUNCT__
711534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
712534831adSStefano Zampini /* -------------------------------------------------------------------------- */
713534831adSStefano Zampini /*
714534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
715534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
716534831adSStefano Zampini 
717534831adSStefano Zampini    Input Parameter:
718534831adSStefano Zampini +  pc - the preconditioner contex
719534831adSStefano Zampini 
720534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
721534831adSStefano Zampini 
722534831adSStefano Zampini    Notes:
723534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
724534831adSStefano Zampini    the user, but instead is called by KSPSolve().
725534831adSStefano Zampini */
726534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
727534831adSStefano Zampini {
728534831adSStefano Zampini   PetscErrorCode ierr;
729534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
730534831adSStefano Zampini   PC_IS          *pcis   = (PC_IS*)(pc->data);
731534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
732534831adSStefano Zampini   Mat            temp_mat;
733534831adSStefano Zampini 
734534831adSStefano Zampini   PetscFunctionBegin;
735674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
736534831adSStefano Zampini     /* swap pointers for local matrices */
737534831adSStefano Zampini     temp_mat = matis->A;
738534831adSStefano Zampini     matis->A = pcbddc->local_mat;
739534831adSStefano Zampini     pcbddc->local_mat = temp_mat;
7403425bc38SStefano Zampini   }
7413308cffdSStefano Zampini   if (pcbddc->use_change_of_basis && x) {
742534831adSStefano Zampini     /* Get Local boundary and apply transformation of basis to solution vector */
743534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
744534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
745534831adSStefano Zampini     /* from modified basis to original basis */
746534831adSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
747534831adSStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
748534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
749534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
750534831adSStefano Zampini   }
7513972b0daSStefano Zampini   /* add solution removed in presolve */
7523425bc38SStefano Zampini   if (x) {
7533425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
7543425bc38SStefano Zampini   }
755fb223d50SStefano Zampini   /* restore rhs to its original state */
756fb223d50SStefano Zampini   if (rhs) {
757fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
758fb223d50SStefano Zampini   }
759534831adSStefano Zampini   PetscFunctionReturn(0);
760534831adSStefano Zampini }
761534831adSStefano Zampini /* -------------------------------------------------------------------------- */
76253cdbc3dSStefano Zampini #undef __FUNCT__
76353cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
7640c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
7650c7d97c5SJed Brown /*
7660c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
7670c7d97c5SJed Brown                   by setting data structures and options.
7680c7d97c5SJed Brown 
7690c7d97c5SJed Brown    Input Parameter:
77053cdbc3dSStefano Zampini +  pc - the preconditioner context
7710c7d97c5SJed Brown 
7720c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
7730c7d97c5SJed Brown 
7740c7d97c5SJed Brown    Notes:
7750c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
7760c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
7770c7d97c5SJed Brown */
77853cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
7790c7d97c5SJed Brown {
7800c7d97c5SJed Brown   PetscErrorCode ierr;
7810c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
782674ae819SStefano Zampini   MatStructure   flag;
783674ae819SStefano Zampini   PetscBool      computeis,computetopography,computesolvers;
7840c7d97c5SJed Brown 
7850c7d97c5SJed Brown   PetscFunctionBegin;
786674ae819SStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */
787fcd409f5SStefano Zampini   /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */
7883b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
789fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
7903b03a366Sstefano_zampini   /* Get stdout for dbg */
791674ae819SStefano Zampini   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
792ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
793e269702eSStefano Zampini     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
7942b510759SStefano Zampini     if (pcbddc->current_level) {
7952b510759SStefano Zampini       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
7962b510759SStefano Zampini     }
797e269702eSStefano Zampini   }
798674ae819SStefano Zampini   /* first attempt to split work */
799674ae819SStefano Zampini   if (pc->setupcalled) {
800674ae819SStefano Zampini     computeis = PETSC_FALSE;
801674ae819SStefano Zampini     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
802674ae819SStefano Zampini     if (flag == SAME_PRECONDITIONER) {
803674ae819SStefano Zampini       computetopography = PETSC_FALSE;
804674ae819SStefano Zampini       computesolvers = PETSC_FALSE;
805674ae819SStefano Zampini     } else if (flag == SAME_NONZERO_PATTERN) {
806674ae819SStefano Zampini       computetopography = PETSC_FALSE;
807674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
808674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
809674ae819SStefano Zampini       computetopography = PETSC_TRUE;
810674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
811674ae819SStefano Zampini     }
812674ae819SStefano Zampini   } else {
813674ae819SStefano Zampini     computeis = PETSC_TRUE;
814674ae819SStefano Zampini     computetopography = PETSC_TRUE;
815674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
816674ae819SStefano Zampini   }
817fcd409f5SStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
818674ae819SStefano Zampini   if (computeis) {
819fcd409f5SStefano Zampini     /* HACK INTO PCIS */
820fcd409f5SStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
821fcd409f5SStefano Zampini     pcis->computesolvers = PETSC_FALSE;
822674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
823674ae819SStefano Zampini   }
824674ae819SStefano Zampini   /* Analyze interface and set up local constraint and change of basis matrices */
825674ae819SStefano Zampini   if (computetopography) {
826674ae819SStefano Zampini     /* reset data */
827674ae819SStefano Zampini     ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
828674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
829674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
830674ae819SStefano Zampini   }
831674ae819SStefano Zampini   if (computesolvers) {
832674ae819SStefano Zampini     /* reset data */
833674ae819SStefano Zampini     ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
834674ae819SStefano Zampini     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
83599cc7994SStefano Zampini     /* Create coarse and local stuffs */
83699cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
837674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
8380c7d97c5SJed Brown   }
8392b510759SStefano Zampini   if (pcbddc->dbg_flag && pcbddc->current_level) {
8402b510759SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
8412b510759SStefano Zampini   }
8420c7d97c5SJed Brown   PetscFunctionReturn(0);
8430c7d97c5SJed Brown }
8440c7d97c5SJed Brown 
8450c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
8460c7d97c5SJed Brown /*
8470c7d97c5SJed Brown    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
8480c7d97c5SJed Brown 
8490c7d97c5SJed Brown    Input Parameters:
8500c7d97c5SJed Brown .  pc - the preconditioner context
8510c7d97c5SJed Brown .  r - input vector (global)
8520c7d97c5SJed Brown 
8530c7d97c5SJed Brown    Output Parameter:
8540c7d97c5SJed Brown .  z - output vector (global)
8550c7d97c5SJed Brown 
8560c7d97c5SJed Brown    Application Interface Routine: PCApply()
8570c7d97c5SJed Brown  */
8580c7d97c5SJed Brown #undef __FUNCT__
8590c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
86053cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
8610c7d97c5SJed Brown {
8620c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
8630c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
8640c7d97c5SJed Brown   PetscErrorCode    ierr;
8653b03a366Sstefano_zampini   const PetscScalar one = 1.0;
8663b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
8672617d88aSStefano Zampini   const PetscScalar zero = 0.0;
8680c7d97c5SJed Brown 
8690c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
8700c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
8718eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
8720c7d97c5SJed Brown 
8730c7d97c5SJed Brown   PetscFunctionBegin;
87485c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
8750c7d97c5SJed Brown     /* First Dirichlet solve */
8760c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8770c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
87853cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
8790c7d97c5SJed Brown     /*
8800c7d97c5SJed Brown       Assembling right hand side for BDDC operator
881674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
882674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
8830c7d97c5SJed Brown     */
8840c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
8850c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
8868eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
8870c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
8880c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
8890c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8900c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
891674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
892b76ba322SStefano Zampini   } else {
8930bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
894b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
895674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
896b76ba322SStefano Zampini   }
897b76ba322SStefano Zampini 
8982617d88aSStefano Zampini   /* Apply interface preconditioner
8992617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
9002617d88aSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
9012617d88aSStefano Zampini 
902674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
903674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
9040c7d97c5SJed Brown 
9053b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
9060c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9070c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9080c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
9098eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
91053cdbc3dSStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr);
9110c7d97c5SJed Brown   ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr);
9128eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
9130c7d97c5SJed Brown   ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr);
9140c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9150c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9160c7d97c5SJed Brown   PetscFunctionReturn(0);
9170c7d97c5SJed Brown }
918da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
919674ae819SStefano Zampini 
920da1bb401SStefano Zampini #undef __FUNCT__
921da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
922da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
923da1bb401SStefano Zampini {
924da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
925da1bb401SStefano Zampini   PetscErrorCode ierr;
926da1bb401SStefano Zampini 
927da1bb401SStefano Zampini   PetscFunctionBegin;
928da1bb401SStefano Zampini   /* free data created by PCIS */
929da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
930674ae819SStefano Zampini   /* free BDDC custom data  */
931674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
932674ae819SStefano Zampini   /* destroy objects related to topography */
933674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
934674ae819SStefano Zampini   /* free allocated graph structure */
935da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
936674ae819SStefano Zampini   /* free data for scaling operator */
937674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
938674ae819SStefano Zampini   /* free solvers stuff */
939674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
94033bc96a4SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
94133bc96a4SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
9429881197aSStefano Zampini   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
943ac78edfcSStefano Zampini   ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
94462a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
94562a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
94662a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
9473425bc38SStefano Zampini   /* remove functions */
948674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
949bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
9502b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
951b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
9522b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
953bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
954bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
955bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
956bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
957bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
958bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
959bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
960bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
961bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
962bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
963674ae819SStefano Zampini   /* Free the private data structure */
964674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
965da1bb401SStefano Zampini   PetscFunctionReturn(0);
966da1bb401SStefano Zampini }
9673425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
9681e6b0712SBarry Smith 
9693425bc38SStefano Zampini #undef __FUNCT__
9703425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
9713425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
9723425bc38SStefano Zampini {
973674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
9743425bc38SStefano Zampini   PC_IS*         pcis;
9753425bc38SStefano Zampini   PC_BDDC*       pcbddc;
9763425bc38SStefano Zampini   PetscErrorCode ierr;
9770c7d97c5SJed Brown 
9783425bc38SStefano Zampini   PetscFunctionBegin;
9793425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
9803425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
9813425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
9823425bc38SStefano Zampini 
9833425bc38SStefano Zampini   /* change of basis for physical rhs if needed
9843425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
9853308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
9863425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
9873425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9883425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
989fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
990fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
9913425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9923425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
993674ae819SStefano Zampini   /* Apply partition of unity */
9943425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
995674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
9968eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
9973425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
9983425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
9993425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
10003425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
10013425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
10023425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10033425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1004674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
10053425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10063425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10073425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
10083425bc38SStefano Zampini   }
10093425bc38SStefano Zampini   /* BDDC rhs */
10103425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
10118eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
10123425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10133425bc38SStefano Zampini   }
10143425bc38SStefano Zampini   /* apply BDDC */
10153425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
10163425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
10173425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
10183425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
10193425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10203425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10213425bc38SStefano Zampini   /* restore original rhs */
10223425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
10233425bc38SStefano Zampini   PetscFunctionReturn(0);
10243425bc38SStefano Zampini }
10251e6b0712SBarry Smith 
10263425bc38SStefano Zampini #undef __FUNCT__
10273425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
10283425bc38SStefano Zampini /*@
102928509bceSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
10303425bc38SStefano Zampini 
10313425bc38SStefano Zampini    Collective
10323425bc38SStefano Zampini 
10333425bc38SStefano Zampini    Input Parameters:
103428509bceSStefano Zampini +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
103528509bceSStefano Zampini .  standard_rhs - the right-hand side for your linear system
10363425bc38SStefano Zampini 
10373425bc38SStefano Zampini    Output Parameters:
103828509bceSStefano Zampini -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
10393425bc38SStefano Zampini 
10403425bc38SStefano Zampini    Level: developer
10413425bc38SStefano Zampini 
10423425bc38SStefano Zampini    Notes:
10433425bc38SStefano Zampini 
104428509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
10453425bc38SStefano Zampini @*/
10463425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
10473425bc38SStefano Zampini {
1048674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10493425bc38SStefano Zampini   PetscErrorCode ierr;
10503425bc38SStefano Zampini 
10513425bc38SStefano Zampini   PetscFunctionBegin;
10523425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10533425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
10543425bc38SStefano Zampini   PetscFunctionReturn(0);
10553425bc38SStefano Zampini }
10563425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
10571e6b0712SBarry Smith 
10583425bc38SStefano Zampini #undef __FUNCT__
10593425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
10603425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
10613425bc38SStefano Zampini {
1062674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10633425bc38SStefano Zampini   PC_IS*         pcis;
10643425bc38SStefano Zampini   PC_BDDC*       pcbddc;
10653425bc38SStefano Zampini   PetscErrorCode ierr;
10663425bc38SStefano Zampini 
10673425bc38SStefano Zampini   PetscFunctionBegin;
10683425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10693425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
10703425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
10713425bc38SStefano Zampini 
10723425bc38SStefano Zampini   /* apply B_delta^T */
10733425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10743425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10753425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
10763425bc38SStefano Zampini   /* compute rhs for BDDC application */
10773425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
10788eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
10793425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10803425bc38SStefano Zampini   }
10813425bc38SStefano Zampini   /* apply BDDC */
10823425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
10833425bc38SStefano Zampini   /* put values into standard global vector */
10843425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10853425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10868eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
10873425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
10883425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
10893425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
10903425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10913425bc38SStefano Zampini   }
10923425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10933425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10943425bc38SStefano Zampini   /* final change of basis if needed
10953425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
10963308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
10973425bc38SStefano Zampini   PetscFunctionReturn(0);
10983425bc38SStefano Zampini }
10991e6b0712SBarry Smith 
11003425bc38SStefano Zampini #undef __FUNCT__
11013425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
11023425bc38SStefano Zampini /*@
110328509bceSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
11043425bc38SStefano Zampini 
11053425bc38SStefano Zampini    Collective
11063425bc38SStefano Zampini 
11073425bc38SStefano Zampini    Input Parameters:
110828509bceSStefano Zampini +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
110928509bceSStefano Zampini .  fetidp_flux_sol - the solution of the FETIDP linear system
11103425bc38SStefano Zampini 
11113425bc38SStefano Zampini    Output Parameters:
111228509bceSStefano Zampini -  standard_sol      - the solution defined on the physical domain
11133425bc38SStefano Zampini 
11143425bc38SStefano Zampini    Level: developer
11153425bc38SStefano Zampini 
11163425bc38SStefano Zampini    Notes:
11173425bc38SStefano Zampini 
111828509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
11193425bc38SStefano Zampini @*/
11203425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
11213425bc38SStefano Zampini {
1122674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
11233425bc38SStefano Zampini   PetscErrorCode ierr;
11243425bc38SStefano Zampini 
11253425bc38SStefano Zampini   PetscFunctionBegin;
11263425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
11273425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
11283425bc38SStefano Zampini   PetscFunctionReturn(0);
11293425bc38SStefano Zampini }
11303425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
11311e6b0712SBarry Smith 
1132f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1133f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1134f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1135f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1136674ae819SStefano Zampini 
11373425bc38SStefano Zampini #undef __FUNCT__
11383425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
11393425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
11403425bc38SStefano Zampini {
1141674ae819SStefano Zampini 
1142674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
11433425bc38SStefano Zampini   Mat            newmat;
1144674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
11453425bc38SStefano Zampini   PC             newpc;
1146ce94432eSBarry Smith   MPI_Comm       comm;
11473425bc38SStefano Zampini   PetscErrorCode ierr;
11483425bc38SStefano Zampini 
11493425bc38SStefano Zampini   PetscFunctionBegin;
1150ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
11513425bc38SStefano Zampini   /* FETIDP linear matrix */
11523425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
11533425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
11543425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
11553425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
11563425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
11573425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
11583425bc38SStefano Zampini   /* FETIDP preconditioner */
11593425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
11603425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
11613425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
11623425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
11633425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
11643425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
11653425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
11663425bc38SStefano Zampini   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
11673425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
11683425bc38SStefano Zampini   /* return pointers for objects created */
11693425bc38SStefano Zampini   *fetidp_mat=newmat;
11703425bc38SStefano Zampini   *fetidp_pc=newpc;
11713425bc38SStefano Zampini   PetscFunctionReturn(0);
11723425bc38SStefano Zampini }
11731e6b0712SBarry Smith 
11743425bc38SStefano Zampini #undef __FUNCT__
11753425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
11763425bc38SStefano Zampini /*@
117728509bceSStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
11783425bc38SStefano Zampini 
11793425bc38SStefano Zampini    Collective
11803425bc38SStefano Zampini 
11813425bc38SStefano Zampini    Input Parameters:
118228509bceSStefano Zampini +  pc - the BDDC preconditioning context already setup
118328509bceSStefano Zampini 
118428509bceSStefano Zampini    Output Parameters:
118528509bceSStefano Zampini .  fetidp_mat - shell FETIDP matrix object
118628509bceSStefano Zampini .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
118728509bceSStefano Zampini 
118828509bceSStefano Zampini    Options Database Keys:
118928509bceSStefano Zampini -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
11903425bc38SStefano Zampini 
11913425bc38SStefano Zampini    Level: developer
11923425bc38SStefano Zampini 
11933425bc38SStefano Zampini    Notes:
119428509bceSStefano Zampini      Currently the only operation provided for FETIDP matrix is MatMult
11953425bc38SStefano Zampini 
11963425bc38SStefano Zampini .seealso: PCBDDC
11973425bc38SStefano Zampini @*/
11983425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
11993425bc38SStefano Zampini {
12003425bc38SStefano Zampini   PetscErrorCode ierr;
12013425bc38SStefano Zampini 
12023425bc38SStefano Zampini   PetscFunctionBegin;
12033425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12043425bc38SStefano Zampini   if (pc->setupcalled) {
1205516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1206f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
12073425bc38SStefano Zampini   PetscFunctionReturn(0);
12083425bc38SStefano Zampini }
12090c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1210da1bb401SStefano Zampini /*MC
1211da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
12120c7d97c5SJed Brown 
121328509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
121428509bceSStefano Zampini 
121528509bceSStefano Zampini .vb
121628509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
121728509bceSStefano 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
121828509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
121928509bceSStefano Zampini .ve
122028509bceSStefano Zampini 
122128509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
122228509bceSStefano Zampini 
1223*b6fdb6dfSStefano Zampini    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
122428509bceSStefano Zampini 
122528509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
122628509bceSStefano Zampini 
1227*b6fdb6dfSStefano 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.
1228*b6fdb6dfSStefano Zampini 
122928509bceSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
123028509bceSStefano Zampini 
123128509bceSStefano Zampini    Boundary nodes are split in vertices, edges and faces using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph
123228509bceSStefano Zampini 
123328509bceSStefano Zampini    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace.
123428509bceSStefano Zampini 
1235*b6fdb6dfSStefano 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.
123628509bceSStefano Zampini 
123728509bceSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
123828509bceSStefano Zampini 
1239da1bb401SStefano Zampini    Options Database Keys:
124028509bceSStefano Zampini 
124128509bceSStefano Zampini .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
124228509bceSStefano Zampini .    -pc_bddc_use_edges <1> - use or not edges in primal space
1243*b6fdb6dfSStefano Zampini .    -pc_bddc_use_faces <0> - use or not faces in primal space
124428509bceSStefano Zampini .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
124528509bceSStefano Zampini .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
124628509bceSStefano Zampini .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
124728509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
124828509bceSStefano Zampini .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
124928509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
125028509bceSStefano Zampini 
125128509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
125228509bceSStefano Zampini .vb
125328509bceSStefano Zampini       -pc_bddc_dirichlet_
125428509bceSStefano Zampini       -pc_bddc_neumann_
125528509bceSStefano Zampini       -pc_bddc_coarse_
125628509bceSStefano Zampini .ve
125728509bceSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
125828509bceSStefano Zampini 
125928509bceSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level can be specified as
126028509bceSStefano Zampini .vb
126128509bceSStefano Zampini       -pc_bddc_dirichlet_N_
126228509bceSStefano Zampini       -pc_bddc_neumann_N_
126328509bceSStefano Zampini       -pc_bddc_coarse_N_
126428509bceSStefano Zampini .ve
126528509bceSStefano Zampini    Note that level number ranges from the finest 0 to the coarsest N
1266da1bb401SStefano Zampini 
1267da1bb401SStefano Zampini    Level: intermediate
1268da1bb401SStefano Zampini 
1269*b6fdb6dfSStefano Zampini    Developer notes:
127028509bceSStefano Zampini      Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose
1271da1bb401SStefano Zampini 
127228509bceSStefano Zampini      New deluxe scaling operator will be available soon.
1273da1bb401SStefano Zampini 
1274da1bb401SStefano Zampini    Contributed by Stefano Zampini
1275da1bb401SStefano Zampini 
1276da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1277da1bb401SStefano Zampini M*/
1278b2573a8aSBarry Smith 
1279da1bb401SStefano Zampini #undef __FUNCT__
1280da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
12818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1282da1bb401SStefano Zampini {
1283da1bb401SStefano Zampini   PetscErrorCode      ierr;
1284da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1285da1bb401SStefano Zampini 
1286da1bb401SStefano Zampini   PetscFunctionBegin;
1287da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1288da1bb401SStefano Zampini   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1289da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1290da1bb401SStefano Zampini 
1291da1bb401SStefano Zampini   /* create PCIS data structure */
1292da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1293da1bb401SStefano Zampini 
129447d04d0dSStefano Zampini   /* BDDC customization */
129547d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
129647d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
129747d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
129847d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
129947d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
130047d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
130147d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
130247d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
130347d04d0dSStefano Zampini 
1304674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
13050bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
13063972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1307534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1308534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1309534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1310da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1311da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1312da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1313da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1314da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
131515aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
131615aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1317da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1318da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1319da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1320da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1321da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1322da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1323da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1324da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1325da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1326da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
132760d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
132860d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
1329da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
1330da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
133185c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
133247d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
133347d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
13344fad6a16SStefano Zampini   pcbddc->current_level              = 0;
13352b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1336da1bb401SStefano Zampini 
1337674ae819SStefano Zampini   /* create local graph structure */
1338674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1339674ae819SStefano Zampini 
1340674ae819SStefano Zampini   /* scaling */
1341674ae819SStefano Zampini   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1342674ae819SStefano Zampini   pcbddc->work_scaling               = 0;
1343da1bb401SStefano Zampini 
1344da1bb401SStefano Zampini   /* function pointers */
1345da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
1346da1bb401SStefano Zampini   pc->ops->applytranspose      = 0;
1347da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1348da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1349da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1350da1bb401SStefano Zampini   pc->ops->view                = 0;
1351da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1352da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1353da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1354534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1355534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1356da1bb401SStefano Zampini 
1357da1bb401SStefano Zampini   /* composing function */
1358674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1359bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
13602b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1361b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
13622b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1363bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1364bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1365bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1366bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1367bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1368bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1369bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1370bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1371bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1372bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1373da1bb401SStefano Zampini   PetscFunctionReturn(0);
1374da1bb401SStefano Zampini }
13753425bc38SStefano Zampini 
1376