xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision be4bc59394babe1c31a7c6478ffb33b3e61f1848)
153cdbc3dSStefano Zampini /* TODOLIST
2*be4bc593SStefano Zampini    Don't use qr when number of primal dof per cc is 1
3*be4bc593SStefano Zampini    A_RR problem with bs>1!
485c4d303SStefano Zampini    Provide PCApplyTranpose
5*be4bc593SStefano Zampini    Man pages
6*be4bc593SStefano Zampini    Move FETIDP code
7*be4bc593SStefano Zampini    Provide general case for subassembling
8*be4bc593SStefano Zampini    Preallocation routines in MatConvert_IS_AIJ
9*be4bc593SStefano Zampini    Allow different customizations between solves
10*be4bc593SStefano Zampini    Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver?
11*be4bc593SStefano Zampini    Better management in PCIS code
1285c4d303SStefano Zampini    Is it possible working with PCBDDCGraph on boundary indices only?
13da1bb401SStefano Zampini    DofSplitting and DM attached to pc?
14da1bb401SStefano Zampini    Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet)
15b8ffe317SStefano Zampini    BDDC with MG framework?
1653cdbc3dSStefano Zampini */
170c7d97c5SJed Brown 
1853cdbc3dSStefano Zampini /* ----------------------------------------------------------------------------------------------------------------------------------------------
190c7d97c5SJed Brown    Implementation of BDDC preconditioner based on:
200c7d97c5SJed Brown    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
2153cdbc3dSStefano Zampini    ---------------------------------------------------------------------------------------------------------------------------------------------- */
2253cdbc3dSStefano Zampini 
23674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/  /* includes for fortran wrappers */
24674ae819SStefano Zampini #include "bddcprivate.h"
253b03a366Sstefano_zampini #include <petscblaslapack.h>
26674ae819SStefano Zampini 
270c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
280c7d97c5SJed Brown #undef __FUNCT__
290c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
300c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc)
310c7d97c5SJed Brown {
320c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
330c7d97c5SJed Brown   PetscErrorCode ierr;
340c7d97c5SJed Brown 
350c7d97c5SJed Brown   PetscFunctionBegin;
360c7d97c5SJed Brown   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
378eeda7d8SStefano Zampini   /* Verbose debugging */
388eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
398eeda7d8SStefano Zampini   /* Primal space cumstomization */
408eeda7d8SStefano 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);
418eeda7d8SStefano 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);
428eeda7d8SStefano 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);
438eeda7d8SStefano Zampini   /* Change of basis */
448eeda7d8SStefano 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);
458eeda7d8SStefano 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);
46674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
47674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
48674ae819SStefano Zampini   }
498eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
508eeda7d8SStefano 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);
510298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
522b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
53674ae819SStefano 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);
540c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
550c7d97c5SJed Brown   PetscFunctionReturn(0);
560c7d97c5SJed Brown }
570c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
58674ae819SStefano Zampini #undef __FUNCT__
59674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
60674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
61674ae819SStefano Zampini {
62674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
63674ae819SStefano Zampini   PetscErrorCode ierr;
641e6b0712SBarry Smith 
65674ae819SStefano Zampini   PetscFunctionBegin;
66674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
67674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
68674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
69674ae819SStefano Zampini   PetscFunctionReturn(0);
70674ae819SStefano Zampini }
71674ae819SStefano Zampini #undef __FUNCT__
72674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
73674ae819SStefano Zampini /*@
74674ae819SStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set user defined primal vertices in PCBDDC.
75674ae819SStefano Zampini 
76674ae819SStefano Zampini    Not collective
77674ae819SStefano Zampini 
78674ae819SStefano Zampini    Input Parameters:
79674ae819SStefano Zampini +  pc - the preconditioning context
80674ae819SStefano Zampini -  PrimalVertices - index sets of primal vertices in local numbering
81674ae819SStefano Zampini 
82674ae819SStefano Zampini    Level: intermediate
83674ae819SStefano Zampini 
84674ae819SStefano Zampini    Notes:
85674ae819SStefano Zampini 
86674ae819SStefano Zampini .seealso: PCBDDC
87674ae819SStefano Zampini @*/
88674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
89674ae819SStefano Zampini {
90674ae819SStefano Zampini   PetscErrorCode ierr;
91674ae819SStefano Zampini 
92674ae819SStefano Zampini   PetscFunctionBegin;
93674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
94674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
95674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
96674ae819SStefano Zampini   PetscFunctionReturn(0);
97674ae819SStefano Zampini }
98674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
990c7d97c5SJed Brown #undef __FUNCT__
1004fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1014fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1024fad6a16SStefano Zampini {
1034fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1044fad6a16SStefano Zampini 
1054fad6a16SStefano Zampini   PetscFunctionBegin;
1064fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1074fad6a16SStefano Zampini   PetscFunctionReturn(0);
1084fad6a16SStefano Zampini }
1091e6b0712SBarry Smith 
1104fad6a16SStefano Zampini #undef __FUNCT__
1114fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1124fad6a16SStefano Zampini /*@
1134fad6a16SStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel coarsening
1144fad6a16SStefano Zampini 
1154fad6a16SStefano Zampini    Logically collective on PC
1164fad6a16SStefano Zampini 
1174fad6a16SStefano Zampini    Input Parameters:
1184fad6a16SStefano Zampini +  pc - the preconditioning context
1194fad6a16SStefano Zampini -  k - coarsening ratio
1204fad6a16SStefano Zampini 
1214fad6a16SStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level.
1224fad6a16SStefano Zampini 
1234fad6a16SStefano Zampini    Level: intermediate
1244fad6a16SStefano Zampini 
1254fad6a16SStefano Zampini    Notes:
1264fad6a16SStefano Zampini 
1274fad6a16SStefano Zampini .seealso: PCBDDC
1284fad6a16SStefano Zampini @*/
1294fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
1304fad6a16SStefano Zampini {
1314fad6a16SStefano Zampini   PetscErrorCode ierr;
1324fad6a16SStefano Zampini 
1334fad6a16SStefano Zampini   PetscFunctionBegin;
1344fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1352b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
1364fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
1374fad6a16SStefano Zampini   PetscFunctionReturn(0);
1384fad6a16SStefano Zampini }
1392b510759SStefano Zampini 
140b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
1412b510759SStefano Zampini #undef __FUNCT__
142b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
143b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
144b8ffe317SStefano Zampini {
145b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
146b8ffe317SStefano Zampini 
147b8ffe317SStefano Zampini   PetscFunctionBegin;
14885c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
149b8ffe317SStefano Zampini   PetscFunctionReturn(0);
150b8ffe317SStefano Zampini }
151b8ffe317SStefano Zampini 
152b8ffe317SStefano Zampini #undef __FUNCT__
153b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
154b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
1552b510759SStefano Zampini {
1562b510759SStefano Zampini   PetscErrorCode ierr;
1572b510759SStefano Zampini 
1582b510759SStefano Zampini   PetscFunctionBegin;
1592b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
160b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
161b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1622b510759SStefano Zampini   PetscFunctionReturn(0);
1632b510759SStefano Zampini }
1641e6b0712SBarry Smith 
1654fad6a16SStefano Zampini #undef __FUNCT__
1662b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
1672b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
1684fad6a16SStefano Zampini {
1694fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1704fad6a16SStefano Zampini 
1714fad6a16SStefano Zampini   PetscFunctionBegin;
1722b510759SStefano Zampini   pcbddc->current_level = level;
1734fad6a16SStefano Zampini   PetscFunctionReturn(0);
1744fad6a16SStefano Zampini }
1751e6b0712SBarry Smith 
1764fad6a16SStefano Zampini #undef __FUNCT__
177b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
178b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
179b8ffe317SStefano Zampini {
180b8ffe317SStefano Zampini   PetscErrorCode ierr;
181b8ffe317SStefano Zampini 
182b8ffe317SStefano Zampini   PetscFunctionBegin;
183b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
184b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
185b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
186b8ffe317SStefano Zampini   PetscFunctionReturn(0);
187b8ffe317SStefano Zampini }
188b8ffe317SStefano Zampini 
189b8ffe317SStefano Zampini #undef __FUNCT__
1902b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
1912b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
1922b510759SStefano Zampini {
1932b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1942b510759SStefano Zampini 
1952b510759SStefano Zampini   PetscFunctionBegin;
1962b510759SStefano Zampini   pcbddc->max_levels = levels;
1972b510759SStefano Zampini   PetscFunctionReturn(0);
1982b510759SStefano Zampini }
1992b510759SStefano Zampini 
2002b510759SStefano Zampini #undef __FUNCT__
2012b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2024fad6a16SStefano Zampini /*@
2032b510759SStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels within the multilevel approach.
2044fad6a16SStefano Zampini 
2054fad6a16SStefano Zampini    Logically collective on PC
2064fad6a16SStefano Zampini 
2074fad6a16SStefano Zampini    Input Parameters:
2084fad6a16SStefano Zampini +  pc - the preconditioning context
2092b510759SStefano Zampini -  levels - the maximum number of levels
2104fad6a16SStefano Zampini 
2112b510759SStefano Zampini    Default value is 0, i.e. coarse problem will be solved exactly
2124fad6a16SStefano Zampini 
2134fad6a16SStefano Zampini    Level: intermediate
2144fad6a16SStefano Zampini 
2154fad6a16SStefano Zampini    Notes:
2164fad6a16SStefano Zampini 
2174fad6a16SStefano Zampini .seealso: PCBDDC
2184fad6a16SStefano Zampini @*/
2192b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
2204fad6a16SStefano Zampini {
2214fad6a16SStefano Zampini   PetscErrorCode ierr;
2224fad6a16SStefano Zampini 
2234fad6a16SStefano Zampini   PetscFunctionBegin;
2244fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2252b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
2262b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
2274fad6a16SStefano Zampini   PetscFunctionReturn(0);
2284fad6a16SStefano Zampini }
2294fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
2301e6b0712SBarry Smith 
2314fad6a16SStefano Zampini #undef __FUNCT__
2320bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
2330bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
2340bdf917eSStefano Zampini {
2350bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2360bdf917eSStefano Zampini   PetscErrorCode ierr;
2370bdf917eSStefano Zampini 
2380bdf917eSStefano Zampini   PetscFunctionBegin;
2390bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
2400bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
2410bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
2420bdf917eSStefano Zampini   PetscFunctionReturn(0);
2430bdf917eSStefano Zampini }
2441e6b0712SBarry Smith 
2450bdf917eSStefano Zampini #undef __FUNCT__
2460bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
2470bdf917eSStefano Zampini /*@
2480bdf917eSStefano Zampini  PCBDDCSetNullSpace - Set NullSpace of global operator of BDDC preconditioned mat.
2490bdf917eSStefano Zampini 
2500bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
2510bdf917eSStefano Zampini 
2520bdf917eSStefano Zampini    Input Parameters:
2530bdf917eSStefano Zampini +  pc - the preconditioning context
2540bdf917eSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned.
2550bdf917eSStefano Zampini 
2560bdf917eSStefano Zampini    Level: intermediate
2570bdf917eSStefano Zampini 
2580bdf917eSStefano Zampini    Notes:
2590bdf917eSStefano Zampini 
2600bdf917eSStefano Zampini .seealso: PCBDDC
2610bdf917eSStefano Zampini @*/
2620bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
2630bdf917eSStefano Zampini {
2640bdf917eSStefano Zampini   PetscErrorCode ierr;
2650bdf917eSStefano Zampini 
2660bdf917eSStefano Zampini   PetscFunctionBegin;
2670bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
268674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
2692b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
2700bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
2710bdf917eSStefano Zampini   PetscFunctionReturn(0);
2720bdf917eSStefano Zampini }
2730bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
2741e6b0712SBarry Smith 
2750bdf917eSStefano Zampini #undef __FUNCT__
2763b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
2773b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
2783b03a366Sstefano_zampini {
2793b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2803b03a366Sstefano_zampini   PetscErrorCode ierr;
2813b03a366Sstefano_zampini 
2823b03a366Sstefano_zampini   PetscFunctionBegin;
2833b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
28436e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
28536e030ebSStefano Zampini   pcbddc->DirichletBoundaries=DirichletBoundaries;
2863b03a366Sstefano_zampini   PetscFunctionReturn(0);
2873b03a366Sstefano_zampini }
2881e6b0712SBarry Smith 
2893b03a366Sstefano_zampini #undef __FUNCT__
2903b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
2913b03a366Sstefano_zampini /*@
292da1bb401SStefano Zampini  PCBDDCSetDirichletBoundaries - Set index set defining subdomain part (in local ordering)
293da1bb401SStefano Zampini                               of Dirichlet boundaries for the global problem.
2943b03a366Sstefano_zampini 
2953b03a366Sstefano_zampini    Not collective
2963b03a366Sstefano_zampini 
2973b03a366Sstefano_zampini    Input Parameters:
2983b03a366Sstefano_zampini +  pc - the preconditioning context
2990298fd71SBarry Smith -  DirichletBoundaries - sequential index set defining the subdomain part of Dirichlet boundaries (can be NULL)
3003b03a366Sstefano_zampini 
3013b03a366Sstefano_zampini    Level: intermediate
3023b03a366Sstefano_zampini 
3033b03a366Sstefano_zampini    Notes:
3043b03a366Sstefano_zampini 
3053b03a366Sstefano_zampini .seealso: PCBDDC
3063b03a366Sstefano_zampini @*/
3073b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3083b03a366Sstefano_zampini {
3093b03a366Sstefano_zampini   PetscErrorCode ierr;
3103b03a366Sstefano_zampini 
3113b03a366Sstefano_zampini   PetscFunctionBegin;
3123b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
313674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
3143b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
3153b03a366Sstefano_zampini   PetscFunctionReturn(0);
3163b03a366Sstefano_zampini }
3173b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
3181e6b0712SBarry Smith 
3193b03a366Sstefano_zampini #undef __FUNCT__
3200c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
32153cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
3220c7d97c5SJed Brown {
3230c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
32453cdbc3dSStefano Zampini   PetscErrorCode ierr;
3250c7d97c5SJed Brown 
3260c7d97c5SJed Brown   PetscFunctionBegin;
32753cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
32836e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
32936e030ebSStefano Zampini   pcbddc->NeumannBoundaries=NeumannBoundaries;
3300c7d97c5SJed Brown   PetscFunctionReturn(0);
3310c7d97c5SJed Brown }
3321e6b0712SBarry Smith 
3330c7d97c5SJed Brown #undef __FUNCT__
3340c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
33557527edcSJed Brown /*@
336da1bb401SStefano Zampini  PCBDDCSetNeumannBoundaries - Set index set defining subdomain part (in local ordering)
337da1bb401SStefano Zampini                               of Neumann boundaries for the global problem.
33857527edcSJed Brown 
3399c0446d6SStefano Zampini    Not collective
34057527edcSJed Brown 
34157527edcSJed Brown    Input Parameters:
34257527edcSJed Brown +  pc - the preconditioning context
3430298fd71SBarry Smith -  NeumannBoundaries - sequential index set defining the subdomain part of Neumann boundaries (can be NULL)
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 /*@
377da1bb401SStefano Zampini  PCBDDCGetDirichletBoundaries - Get index set defining subdomain part (in local ordering)
378da1bb401SStefano Zampini                                 of Dirichlet boundaries for the global problem.
379da1bb401SStefano Zampini 
380da1bb401SStefano Zampini    Not collective
381da1bb401SStefano Zampini 
382da1bb401SStefano Zampini    Input Parameters:
383da1bb401SStefano Zampini +  pc - the preconditioning context
384da1bb401SStefano Zampini 
385da1bb401SStefano Zampini    Output Parameters:
386da1bb401SStefano Zampini +  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
387da1bb401SStefano Zampini 
388da1bb401SStefano Zampini    Level: intermediate
389da1bb401SStefano Zampini 
390da1bb401SStefano Zampini    Notes:
391da1bb401SStefano Zampini 
392da1bb401SStefano Zampini .seealso: PCBDDC
393da1bb401SStefano Zampini @*/
394da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
395da1bb401SStefano Zampini {
396da1bb401SStefano Zampini   PetscErrorCode ierr;
397da1bb401SStefano Zampini 
398da1bb401SStefano Zampini   PetscFunctionBegin;
399da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
400da1bb401SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
401da1bb401SStefano Zampini   PetscFunctionReturn(0);
402da1bb401SStefano Zampini }
403da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
4041e6b0712SBarry Smith 
405da1bb401SStefano Zampini #undef __FUNCT__
40653cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
40753cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
40853cdbc3dSStefano Zampini {
40953cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
41053cdbc3dSStefano Zampini 
41153cdbc3dSStefano Zampini   PetscFunctionBegin;
41253cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
41353cdbc3dSStefano Zampini   PetscFunctionReturn(0);
41453cdbc3dSStefano Zampini }
4151e6b0712SBarry Smith 
41653cdbc3dSStefano Zampini #undef __FUNCT__
41753cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
41853cdbc3dSStefano Zampini /*@
419da1bb401SStefano Zampini  PCBDDCGetNeumannBoundaries - Get index set defining subdomain part (in local ordering)
420da1bb401SStefano Zampini                               of Neumann boundaries for the global problem.
42153cdbc3dSStefano Zampini 
4229c0446d6SStefano Zampini    Not collective
42353cdbc3dSStefano Zampini 
42453cdbc3dSStefano Zampini    Input Parameters:
42553cdbc3dSStefano Zampini +  pc - the preconditioning context
42653cdbc3dSStefano Zampini 
42753cdbc3dSStefano Zampini    Output Parameters:
42853cdbc3dSStefano Zampini +  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
42953cdbc3dSStefano Zampini 
43053cdbc3dSStefano Zampini    Level: intermediate
43153cdbc3dSStefano Zampini 
43253cdbc3dSStefano Zampini    Notes:
43353cdbc3dSStefano Zampini 
43453cdbc3dSStefano Zampini .seealso: PCBDDC
43553cdbc3dSStefano Zampini @*/
43653cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
43753cdbc3dSStefano Zampini {
43853cdbc3dSStefano Zampini   PetscErrorCode ierr;
43953cdbc3dSStefano Zampini 
44053cdbc3dSStefano Zampini   PetscFunctionBegin;
44153cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
44253cdbc3dSStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
4430c7d97c5SJed Brown   PetscFunctionReturn(0);
4440c7d97c5SJed Brown }
44536e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
4461e6b0712SBarry Smith 
44736e030ebSStefano Zampini #undef __FUNCT__
448da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
4491a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
45036e030ebSStefano Zampini {
45136e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
452da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
453da1bb401SStefano Zampini   PetscErrorCode ierr;
45436e030ebSStefano Zampini 
45536e030ebSStefano Zampini   PetscFunctionBegin;
456674ae819SStefano Zampini   /* free old CSR */
457674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
458fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
459674ae819SStefano Zampini   /* get CSR into graph structure */
460da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
461674ae819SStefano Zampini     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
462674ae819SStefano Zampini     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
463674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
464674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
465da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
4661a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
4671a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
468674ae819SStefano Zampini   }
469575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
47036e030ebSStefano Zampini   PetscFunctionReturn(0);
47136e030ebSStefano Zampini }
4721e6b0712SBarry Smith 
47336e030ebSStefano Zampini #undef __FUNCT__
474da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
47536e030ebSStefano Zampini /*@
476da1bb401SStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set CSR graph of local matrix for use of PCBDDC.
47736e030ebSStefano Zampini 
47836e030ebSStefano Zampini    Not collective
47936e030ebSStefano Zampini 
48036e030ebSStefano Zampini    Input Parameters:
48136e030ebSStefano Zampini +  pc - the preconditioning context
482da1bb401SStefano Zampini -  nvtxs - number of local vertices of the graph
483da1bb401SStefano Zampini -  xadj, adjncy - the CSR graph
484da1bb401SStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. In the former case the user must free the array passed in;
485da1bb401SStefano Zampini                                                              in the latter case, memory must be obtained with PetscMalloc.
48636e030ebSStefano Zampini 
48736e030ebSStefano Zampini    Level: intermediate
48836e030ebSStefano Zampini 
48936e030ebSStefano Zampini    Notes:
49036e030ebSStefano Zampini 
49136e030ebSStefano Zampini .seealso: PCBDDC
49236e030ebSStefano Zampini @*/
4931a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
49436e030ebSStefano Zampini {
495575ad6abSStefano Zampini   void (*f)(void) = 0;
49636e030ebSStefano Zampini   PetscErrorCode ierr;
49736e030ebSStefano Zampini 
49836e030ebSStefano Zampini   PetscFunctionBegin;
49936e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
500674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
501674ae819SStefano Zampini   PetscValidIntPointer(xadj,4);
502674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
503674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
504da1bb401SStefano Zampini   }
50536e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
506575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
507575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
508575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
509575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
510575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
51136e030ebSStefano Zampini   }
51236e030ebSStefano Zampini   PetscFunctionReturn(0);
51336e030ebSStefano Zampini }
5149c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
5151e6b0712SBarry Smith 
5169c0446d6SStefano Zampini #undef __FUNCT__
5179c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
5189c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
5199c0446d6SStefano Zampini {
5209c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5219c0446d6SStefano Zampini   PetscInt i;
5229c0446d6SStefano Zampini   PetscErrorCode ierr;
5239c0446d6SStefano Zampini 
5249c0446d6SStefano Zampini   PetscFunctionBegin;
525da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
5269c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
5279c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
5289c0446d6SStefano Zampini   }
529d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
530da1bb401SStefano Zampini   /* allocate space then set */
5319c0446d6SStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
5329c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
533da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
534da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
5359c0446d6SStefano Zampini   }
5369c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
5379c0446d6SStefano Zampini   PetscFunctionReturn(0);
5389c0446d6SStefano Zampini }
5391e6b0712SBarry Smith 
5409c0446d6SStefano Zampini #undef __FUNCT__
5419c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
5429c0446d6SStefano Zampini /*@
543da1bb401SStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of local mat.
5449c0446d6SStefano Zampini 
5459c0446d6SStefano Zampini    Not collective
5469c0446d6SStefano Zampini 
5479c0446d6SStefano Zampini    Input Parameters:
5489c0446d6SStefano Zampini +  pc - the preconditioning context
549da1bb401SStefano Zampini -  n - number of index sets defining the fields
550da1bb401SStefano Zampini -  IS[] - array of IS describing the fields
5519c0446d6SStefano Zampini 
5529c0446d6SStefano Zampini    Level: intermediate
5539c0446d6SStefano Zampini 
5549c0446d6SStefano Zampini    Notes:
5559c0446d6SStefano Zampini 
5569c0446d6SStefano Zampini .seealso: PCBDDC
5579c0446d6SStefano Zampini @*/
5589c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
5599c0446d6SStefano Zampini {
5602b510759SStefano Zampini   PetscInt       i;
5619c0446d6SStefano Zampini   PetscErrorCode ierr;
5629c0446d6SStefano Zampini 
5639c0446d6SStefano Zampini   PetscFunctionBegin;
5649c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5652b510759SStefano Zampini   for (i=0;i<n_is;i++) {
5662b510759SStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2);
5672b510759SStefano Zampini   }
5689c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
5699c0446d6SStefano Zampini   PetscFunctionReturn(0);
5709c0446d6SStefano Zampini }
571da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
572534831adSStefano Zampini #undef __FUNCT__
573534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
574534831adSStefano Zampini /* -------------------------------------------------------------------------- */
575534831adSStefano Zampini /*
576534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
577534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
5789c0446d6SStefano Zampini 
579534831adSStefano Zampini    Input Parameter:
580534831adSStefano Zampini +  pc - the preconditioner contex
581534831adSStefano Zampini 
582534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
583534831adSStefano Zampini 
584534831adSStefano Zampini    Notes:
585534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
586534831adSStefano Zampini    the user, but instead is called by KSPSolve().
587534831adSStefano Zampini */
588534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
589534831adSStefano Zampini {
590534831adSStefano Zampini   PetscErrorCode ierr;
591534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
592534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
593534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
594534831adSStefano Zampini   Mat            temp_mat;
5953972b0daSStefano Zampini   IS             dirIS;
5963972b0daSStefano Zampini   PetscInt       dirsize,i,*is_indices;
5973972b0daSStefano Zampini   PetscScalar    *array_x,*array_diagonal;
5983972b0daSStefano Zampini   Vec            used_vec;
59992e3dcfbSStefano Zampini   PetscBool      guess_nonzero,flg,bddc_has_dirichlet_boundaries;
600534831adSStefano Zampini 
601534831adSStefano Zampini   PetscFunctionBegin;
60285c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
60385c4d303SStefano Zampini   if (ksp) {
60485c4d303SStefano Zampini     PetscBool iscg;
60585c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
60685c4d303SStefano Zampini     if (!iscg) {
60785c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
60885c4d303SStefano Zampini     }
60985c4d303SStefano Zampini   }
61085c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
61162a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
61262a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
61362a6ff1dSStefano Zampini   }
61462a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
61562a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
61662a6ff1dSStefano Zampini   }
6173972b0daSStefano Zampini   if (x) {
6183972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
6193972b0daSStefano Zampini     used_vec = x;
6203972b0daSStefano Zampini   } else {
6213972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
6223972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
6233972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6243972b0daSStefano Zampini   }
6253972b0daSStefano Zampini   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
6263972b0daSStefano Zampini   if (ksp) {
6273972b0daSStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
6283972b0daSStefano Zampini     if (!guess_nonzero) {
6293972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6303972b0daSStefano Zampini     }
6313972b0daSStefano Zampini   }
6323308cffdSStefano Zampini 
63392e3dcfbSStefano Zampini   /* TODO: remove when Dirichlet boundaries will be shared */
63492e3dcfbSStefano Zampini   ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr);
63592e3dcfbSStefano Zampini   flg = PETSC_FALSE;
63692e3dcfbSStefano Zampini   if (dirIS) flg = PETSC_TRUE;
63792e3dcfbSStefano Zampini   ierr = MPI_Allreduce(&flg,&bddc_has_dirichlet_boundaries,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
63892e3dcfbSStefano Zampini 
6393972b0daSStefano Zampini   /* store the original rhs */
6403972b0daSStefano Zampini   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
6413972b0daSStefano Zampini 
6423972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
64392e3dcfbSStefano Zampini   if (rhs && bddc_has_dirichlet_boundaries) {
6443972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
6453972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
6463972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6473972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6483972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6493972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6503972b0daSStefano Zampini     if (dirIS) {
6513972b0daSStefano Zampini       ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr);
6523972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6533972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6543972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6552fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
6563972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6573972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6583972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6593972b0daSStefano Zampini     }
6603972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6613972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
662b76ba322SStefano Zampini 
6633972b0daSStefano Zampini     /* remove the computed solution from the rhs */
6643972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6653972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
6663972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6673308cffdSStefano Zampini   }
668b76ba322SStefano Zampini 
669b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
6703972b0daSStefano Zampini   if (x) {
6713972b0daSStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
6723972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
67385c4d303SStefano Zampini     if (pcbddc->use_exact_dirichlet_trick) {
674b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
675b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
676b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
677b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
678b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
679b76ba322SStefano Zampini       if (ksp) {
680b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
681b76ba322SStefano Zampini       }
682b76ba322SStefano Zampini     }
6833972b0daSStefano Zampini   }
684b76ba322SStefano Zampini 
68592e3dcfbSStefano Zampini   /* prepare MatMult and rhs for solver */
686674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
687b76ba322SStefano Zampini     /* swap pointers for local matrices */
688b76ba322SStefano Zampini     temp_mat = matis->A;
689b76ba322SStefano Zampini     matis->A = pcbddc->local_mat;
690b76ba322SStefano Zampini     pcbddc->local_mat = temp_mat;
69192e3dcfbSStefano Zampini     if (rhs) {
692b76ba322SStefano Zampini       /* Get local rhs and apply transformation of basis */
693b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
694b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
695b76ba322SStefano Zampini       /* from original basis to modified basis */
696b76ba322SStefano Zampini       ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
697b76ba322SStefano Zampini       /* put back modified values into the global vec using INSERT_VALUES copy mode */
698b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
699b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
700674ae819SStefano Zampini     }
70192e3dcfbSStefano Zampini   }
70292e3dcfbSStefano Zampini 
70392e3dcfbSStefano Zampini   /* remove nullspace if present */
7040bdf917eSStefano Zampini   if (ksp && pcbddc->NullSpace) {
705d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
706d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
707b76ba322SStefano Zampini   }
7080bdf917eSStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
709534831adSStefano Zampini   PetscFunctionReturn(0);
710534831adSStefano Zampini }
711534831adSStefano Zampini /* -------------------------------------------------------------------------- */
712534831adSStefano Zampini #undef __FUNCT__
713534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
714534831adSStefano Zampini /* -------------------------------------------------------------------------- */
715534831adSStefano Zampini /*
716534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
717534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
718534831adSStefano Zampini 
719534831adSStefano Zampini    Input Parameter:
720534831adSStefano Zampini +  pc - the preconditioner contex
721534831adSStefano Zampini 
722534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
723534831adSStefano Zampini 
724534831adSStefano Zampini    Notes:
725534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
726534831adSStefano Zampini    the user, but instead is called by KSPSolve().
727534831adSStefano Zampini */
728534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
729534831adSStefano Zampini {
730534831adSStefano Zampini   PetscErrorCode ierr;
731534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
732534831adSStefano Zampini   PC_IS          *pcis   = (PC_IS*)(pc->data);
733534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
734534831adSStefano Zampini   Mat            temp_mat;
735534831adSStefano Zampini 
736534831adSStefano Zampini   PetscFunctionBegin;
737674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
738534831adSStefano Zampini     /* swap pointers for local matrices */
739534831adSStefano Zampini     temp_mat = matis->A;
740534831adSStefano Zampini     matis->A = pcbddc->local_mat;
741534831adSStefano Zampini     pcbddc->local_mat = temp_mat;
7423425bc38SStefano Zampini   }
7433308cffdSStefano Zampini   if (pcbddc->use_change_of_basis && x) {
744534831adSStefano Zampini     /* Get Local boundary and apply transformation of basis to solution vector */
745534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
746534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
747534831adSStefano Zampini     /* from modified basis to original basis */
748534831adSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
749534831adSStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
750534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
751534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
752534831adSStefano Zampini   }
7533972b0daSStefano Zampini   /* add solution removed in presolve */
7543425bc38SStefano Zampini   if (x) {
7553425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
7563425bc38SStefano Zampini   }
757fb223d50SStefano Zampini   /* restore rhs to its original state */
758fb223d50SStefano Zampini   if (rhs) {
759fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
760fb223d50SStefano Zampini   }
761534831adSStefano Zampini   PetscFunctionReturn(0);
762534831adSStefano Zampini }
763534831adSStefano Zampini /* -------------------------------------------------------------------------- */
76453cdbc3dSStefano Zampini #undef __FUNCT__
76553cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
7660c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
7670c7d97c5SJed Brown /*
7680c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
7690c7d97c5SJed Brown                   by setting data structures and options.
7700c7d97c5SJed Brown 
7710c7d97c5SJed Brown    Input Parameter:
77253cdbc3dSStefano Zampini +  pc - the preconditioner context
7730c7d97c5SJed Brown 
7740c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
7750c7d97c5SJed Brown 
7760c7d97c5SJed Brown    Notes:
7770c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
7780c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
7790c7d97c5SJed Brown */
78053cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
7810c7d97c5SJed Brown {
7820c7d97c5SJed Brown   PetscErrorCode ierr;
7830c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
784674ae819SStefano Zampini   MatStructure   flag;
785674ae819SStefano Zampini   PetscBool      computeis,computetopography,computesolvers;
7860c7d97c5SJed Brown 
7870c7d97c5SJed Brown   PetscFunctionBegin;
788674ae819SStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */
789fcd409f5SStefano Zampini   /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */
7903b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
791fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
7923b03a366Sstefano_zampini   /* Get stdout for dbg */
793674ae819SStefano Zampini   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
794ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
795e269702eSStefano Zampini     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
7962b510759SStefano Zampini     if (pcbddc->current_level) {
7972b510759SStefano Zampini       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
7982b510759SStefano Zampini     }
799e269702eSStefano Zampini   }
800674ae819SStefano Zampini   /* first attempt to split work */
801674ae819SStefano Zampini   if (pc->setupcalled) {
802674ae819SStefano Zampini     computeis = PETSC_FALSE;
803674ae819SStefano Zampini     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
804674ae819SStefano Zampini     if (flag == SAME_PRECONDITIONER) {
805674ae819SStefano Zampini       computetopography = PETSC_FALSE;
806674ae819SStefano Zampini       computesolvers = PETSC_FALSE;
807674ae819SStefano Zampini     } else if (flag == SAME_NONZERO_PATTERN) {
808674ae819SStefano Zampini       computetopography = PETSC_FALSE;
809674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
810674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
811674ae819SStefano Zampini       computetopography = PETSC_TRUE;
812674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
813674ae819SStefano Zampini     }
814674ae819SStefano Zampini   } else {
815674ae819SStefano Zampini     computeis = PETSC_TRUE;
816674ae819SStefano Zampini     computetopography = PETSC_TRUE;
817674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
818674ae819SStefano Zampini   }
819fcd409f5SStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
820674ae819SStefano Zampini   if (computeis) {
821fcd409f5SStefano Zampini     /* HACK INTO PCIS */
822fcd409f5SStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
823fcd409f5SStefano Zampini     pcis->computesolvers = PETSC_FALSE;
824674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
825674ae819SStefano Zampini   }
826674ae819SStefano Zampini   /* Analyze interface and set up local constraint and change of basis matrices */
827674ae819SStefano Zampini   if (computetopography) {
828674ae819SStefano Zampini     /* reset data */
829674ae819SStefano Zampini     ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
830674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
831674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
832674ae819SStefano Zampini   }
833674ae819SStefano Zampini   if (computesolvers) {
834674ae819SStefano Zampini     /* reset data */
835674ae819SStefano Zampini     ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
836674ae819SStefano Zampini     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
83799cc7994SStefano Zampini     /* Create coarse and local stuffs */
83899cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
839674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
8400c7d97c5SJed Brown   }
8412b510759SStefano Zampini   if (pcbddc->dbg_flag && pcbddc->current_level) {
8422b510759SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
8432b510759SStefano Zampini   }
8440c7d97c5SJed Brown   PetscFunctionReturn(0);
8450c7d97c5SJed Brown }
8460c7d97c5SJed Brown 
8470c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
8480c7d97c5SJed Brown /*
8490c7d97c5SJed Brown    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
8500c7d97c5SJed Brown 
8510c7d97c5SJed Brown    Input Parameters:
8520c7d97c5SJed Brown .  pc - the preconditioner context
8530c7d97c5SJed Brown .  r - input vector (global)
8540c7d97c5SJed Brown 
8550c7d97c5SJed Brown    Output Parameter:
8560c7d97c5SJed Brown .  z - output vector (global)
8570c7d97c5SJed Brown 
8580c7d97c5SJed Brown    Application Interface Routine: PCApply()
8590c7d97c5SJed Brown  */
8600c7d97c5SJed Brown #undef __FUNCT__
8610c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
86253cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
8630c7d97c5SJed Brown {
8640c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
8650c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
8660c7d97c5SJed Brown   PetscErrorCode    ierr;
8673b03a366Sstefano_zampini   const PetscScalar one = 1.0;
8683b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
8692617d88aSStefano Zampini   const PetscScalar zero = 0.0;
8700c7d97c5SJed Brown 
8710c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
8720c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
8738eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
8740c7d97c5SJed Brown 
8750c7d97c5SJed Brown   PetscFunctionBegin;
87685c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
8770c7d97c5SJed Brown     /* First Dirichlet solve */
8780c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8790c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
88053cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
8810c7d97c5SJed Brown     /*
8820c7d97c5SJed Brown       Assembling right hand side for BDDC operator
883674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
884674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
8850c7d97c5SJed Brown     */
8860c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
8870c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
8888eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
8890c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
8900c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
8910c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8920c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
893674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
894b76ba322SStefano Zampini   } else {
8950bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
896b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
897674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
898b76ba322SStefano Zampini   }
899b76ba322SStefano Zampini 
9002617d88aSStefano Zampini   /* Apply interface preconditioner
9012617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
9022617d88aSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
9032617d88aSStefano Zampini 
904674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
905674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
9060c7d97c5SJed Brown 
9073b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
9080c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9090c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9100c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
9118eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
91253cdbc3dSStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr);
9130c7d97c5SJed Brown   ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr);
9148eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
9150c7d97c5SJed Brown   ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr);
9160c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9170c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9180c7d97c5SJed Brown   PetscFunctionReturn(0);
9190c7d97c5SJed Brown }
920da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
921674ae819SStefano Zampini 
922da1bb401SStefano Zampini #undef __FUNCT__
923da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
924da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
925da1bb401SStefano Zampini {
926da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
927da1bb401SStefano Zampini   PetscErrorCode ierr;
928da1bb401SStefano Zampini 
929da1bb401SStefano Zampini   PetscFunctionBegin;
930da1bb401SStefano Zampini   /* free data created by PCIS */
931da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
932674ae819SStefano Zampini   /* free BDDC custom data  */
933674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
934674ae819SStefano Zampini   /* destroy objects related to topography */
935674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
936674ae819SStefano Zampini   /* free allocated graph structure */
937da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
938674ae819SStefano Zampini   /* free data for scaling operator */
939674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
940674ae819SStefano Zampini   /* free solvers stuff */
941674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
94233bc96a4SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
94333bc96a4SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
9449881197aSStefano Zampini   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
945ac78edfcSStefano Zampini   ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
94662a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
94762a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
94862a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
9493425bc38SStefano Zampini   /* remove functions */
950674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
951bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
9522b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
953b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
9542b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
955bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
956bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
957bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
958bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
959bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
960bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
961bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
962bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
963bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
964bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
965674ae819SStefano Zampini   /* Free the private data structure */
966674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
967da1bb401SStefano Zampini   PetscFunctionReturn(0);
968da1bb401SStefano Zampini }
9693425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
9701e6b0712SBarry Smith 
9713425bc38SStefano Zampini #undef __FUNCT__
9723425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
9733425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
9743425bc38SStefano Zampini {
975674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
9763425bc38SStefano Zampini   PC_IS*         pcis;
9773425bc38SStefano Zampini   PC_BDDC*       pcbddc;
9783425bc38SStefano Zampini   PetscErrorCode ierr;
9790c7d97c5SJed Brown 
9803425bc38SStefano Zampini   PetscFunctionBegin;
9813425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
9823425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
9833425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
9843425bc38SStefano Zampini 
9853425bc38SStefano Zampini   /* change of basis for physical rhs if needed
9863425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
9873308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
9883425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
9893425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9903425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
991fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
992fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
9933425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9943425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
995674ae819SStefano Zampini   /* Apply partition of unity */
9963425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
997674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
9988eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
9993425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
10003425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10013425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
10023425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
10033425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
10043425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10053425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1006674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
10073425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10083425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10093425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
10103425bc38SStefano Zampini   }
10113425bc38SStefano Zampini   /* BDDC rhs */
10123425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
10138eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
10143425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10153425bc38SStefano Zampini   }
10163425bc38SStefano Zampini   /* apply BDDC */
10173425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
10183425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
10193425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
10203425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
10213425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10223425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10233425bc38SStefano Zampini   /* restore original rhs */
10243425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
10253425bc38SStefano Zampini   PetscFunctionReturn(0);
10263425bc38SStefano Zampini }
10271e6b0712SBarry Smith 
10283425bc38SStefano Zampini #undef __FUNCT__
10293425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
10303425bc38SStefano Zampini /*@
10313425bc38SStefano Zampini  PCBDDCMatFETIDPGetRHS - Get rhs for FETIDP linear system.
10323425bc38SStefano Zampini 
10333425bc38SStefano Zampini    Collective
10343425bc38SStefano Zampini 
10353425bc38SStefano Zampini    Input Parameters:
10363425bc38SStefano Zampini +  fetidp_mat   - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators
10373425bc38SStefano Zampini +  standard_rhs - the rhs of your linear system
10383425bc38SStefano Zampini 
10393425bc38SStefano Zampini    Output Parameters:
10403425bc38SStefano Zampini +  fetidp_flux_rhs   - the rhs of the FETIDP linear system
10413425bc38SStefano Zampini 
10423425bc38SStefano Zampini    Level: developer
10433425bc38SStefano Zampini 
10443425bc38SStefano Zampini    Notes:
10453425bc38SStefano Zampini 
10463425bc38SStefano Zampini .seealso: PCBDDC
10473425bc38SStefano Zampini @*/
10483425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
10493425bc38SStefano Zampini {
1050674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10513425bc38SStefano Zampini   PetscErrorCode ierr;
10523425bc38SStefano Zampini 
10533425bc38SStefano Zampini   PetscFunctionBegin;
10543425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10553425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
10563425bc38SStefano Zampini   PetscFunctionReturn(0);
10573425bc38SStefano Zampini }
10583425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
10591e6b0712SBarry Smith 
10603425bc38SStefano Zampini #undef __FUNCT__
10613425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
10623425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
10633425bc38SStefano Zampini {
1064674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10653425bc38SStefano Zampini   PC_IS*         pcis;
10663425bc38SStefano Zampini   PC_BDDC*       pcbddc;
10673425bc38SStefano Zampini   PetscErrorCode ierr;
10683425bc38SStefano Zampini 
10693425bc38SStefano Zampini   PetscFunctionBegin;
10703425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10713425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
10723425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
10733425bc38SStefano Zampini 
10743425bc38SStefano Zampini   /* apply B_delta^T */
10753425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10763425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10773425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
10783425bc38SStefano Zampini   /* compute rhs for BDDC application */
10793425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
10808eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
10813425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10823425bc38SStefano Zampini   }
10833425bc38SStefano Zampini   /* apply BDDC */
10843425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
10853425bc38SStefano Zampini   /* put values into standard global vector */
10863425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10873425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10888eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
10893425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
10903425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
10913425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
10923425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10933425bc38SStefano Zampini   }
10943425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10953425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10963425bc38SStefano Zampini   /* final change of basis if needed
10973425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
10983308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
10993425bc38SStefano Zampini   PetscFunctionReturn(0);
11003425bc38SStefano Zampini }
11011e6b0712SBarry Smith 
11023425bc38SStefano Zampini #undef __FUNCT__
11033425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
11043425bc38SStefano Zampini /*@
11053425bc38SStefano Zampini  PCBDDCMatFETIDPGetSolution - Get Solution for FETIDP linear system.
11063425bc38SStefano Zampini 
11073425bc38SStefano Zampini    Collective
11083425bc38SStefano Zampini 
11093425bc38SStefano Zampini    Input Parameters:
11103425bc38SStefano Zampini +  fetidp_mat        - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators
11113425bc38SStefano Zampini +  fetidp_flux_sol - the solution of the FETIDP linear system
11123425bc38SStefano Zampini 
11133425bc38SStefano Zampini    Output Parameters:
11143425bc38SStefano Zampini +  standard_sol      - the solution on the global domain
11153425bc38SStefano Zampini 
11163425bc38SStefano Zampini    Level: developer
11173425bc38SStefano Zampini 
11183425bc38SStefano Zampini    Notes:
11193425bc38SStefano Zampini 
11203425bc38SStefano Zampini .seealso: PCBDDC
11213425bc38SStefano Zampini @*/
11223425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
11233425bc38SStefano Zampini {
1124674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
11253425bc38SStefano Zampini   PetscErrorCode ierr;
11263425bc38SStefano Zampini 
11273425bc38SStefano Zampini   PetscFunctionBegin;
11283425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
11293425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
11303425bc38SStefano Zampini   PetscFunctionReturn(0);
11313425bc38SStefano Zampini }
11323425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
11331e6b0712SBarry Smith 
1134f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1135f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1136f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1137f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1138674ae819SStefano Zampini 
11393425bc38SStefano Zampini #undef __FUNCT__
11403425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
11413425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
11423425bc38SStefano Zampini {
1143674ae819SStefano Zampini 
1144674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
11453425bc38SStefano Zampini   Mat            newmat;
1146674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
11473425bc38SStefano Zampini   PC             newpc;
1148ce94432eSBarry Smith   MPI_Comm       comm;
11493425bc38SStefano Zampini   PetscErrorCode ierr;
11503425bc38SStefano Zampini 
11513425bc38SStefano Zampini   PetscFunctionBegin;
1152ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
11533425bc38SStefano Zampini   /* FETIDP linear matrix */
11543425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
11553425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
11563425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
11573425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
11583425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
11593425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
11603425bc38SStefano Zampini   /* FETIDP preconditioner */
11613425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
11623425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
11633425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
11643425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
11653425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
11663425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
11673425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
11683425bc38SStefano Zampini   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
11693425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
11703425bc38SStefano Zampini   /* return pointers for objects created */
11713425bc38SStefano Zampini   *fetidp_mat=newmat;
11723425bc38SStefano Zampini   *fetidp_pc=newpc;
11733425bc38SStefano Zampini   PetscFunctionReturn(0);
11743425bc38SStefano Zampini }
11751e6b0712SBarry Smith 
11763425bc38SStefano Zampini #undef __FUNCT__
11773425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
11783425bc38SStefano Zampini /*@
11793425bc38SStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP.
11803425bc38SStefano Zampini 
11813425bc38SStefano Zampini    Collective
11823425bc38SStefano Zampini 
11833425bc38SStefano Zampini    Input Parameters:
11843425bc38SStefano Zampini +  pc - the BDDC preconditioning context (setup must be already called)
11853425bc38SStefano Zampini 
11863425bc38SStefano Zampini    Level: developer
11873425bc38SStefano Zampini 
11883425bc38SStefano Zampini    Notes:
11893425bc38SStefano Zampini 
11903425bc38SStefano Zampini .seealso: PCBDDC
11913425bc38SStefano Zampini @*/
11923425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
11933425bc38SStefano Zampini {
11943425bc38SStefano Zampini   PetscErrorCode ierr;
11953425bc38SStefano Zampini 
11963425bc38SStefano Zampini   PetscFunctionBegin;
11973425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11983425bc38SStefano Zampini   if (pc->setupcalled) {
1199516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1200f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
12013425bc38SStefano Zampini   PetscFunctionReturn(0);
12023425bc38SStefano Zampini }
12030c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1204da1bb401SStefano Zampini /*MC
1205da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
12060c7d97c5SJed Brown 
1207da1bb401SStefano Zampini    Options Database Keys:
1208da1bb401SStefano Zampini .    -pcbddc ??? -
1209da1bb401SStefano Zampini 
1210da1bb401SStefano Zampini    Level: intermediate
1211da1bb401SStefano Zampini 
1212da1bb401SStefano Zampini    Notes: The matrix used with this preconditioner must be of type MATIS
1213da1bb401SStefano Zampini 
1214da1bb401SStefano Zampini           Unlike more 'conventional' interface preconditioners, this iterates over ALL the
1215da1bb401SStefano Zampini           degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers
1216da1bb401SStefano Zampini           on the subdomains).
1217da1bb401SStefano Zampini 
1218da1bb401SStefano Zampini           Options for the coarse grid preconditioner can be set with -
1219da1bb401SStefano Zampini           Options for the Dirichlet subproblem can be set with -
1220da1bb401SStefano Zampini           Options for the Neumann subproblem can be set with -
1221da1bb401SStefano Zampini 
1222da1bb401SStefano Zampini    Contributed by Stefano Zampini
1223da1bb401SStefano Zampini 
1224da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1225da1bb401SStefano Zampini M*/
1226b2573a8aSBarry Smith 
1227da1bb401SStefano Zampini #undef __FUNCT__
1228da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
12298cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1230da1bb401SStefano Zampini {
1231da1bb401SStefano Zampini   PetscErrorCode      ierr;
1232da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1233da1bb401SStefano Zampini 
1234da1bb401SStefano Zampini   PetscFunctionBegin;
1235da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1236da1bb401SStefano Zampini   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1237da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1238da1bb401SStefano Zampini 
1239da1bb401SStefano Zampini   /* create PCIS data structure */
1240da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1241da1bb401SStefano Zampini 
124247d04d0dSStefano Zampini   /* BDDC customization */
124347d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
124447d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
124547d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
124647d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
124747d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
124847d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
124947d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
125047d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
125147d04d0dSStefano Zampini 
1252674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
12530bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
12543972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1255534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1256534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1257534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1258da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1259da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1260da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1261da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1262da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
126315aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
126415aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1265da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1266da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1267da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1268da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1269da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1270da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1271da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1272da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1273da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1274da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
1275da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
1276da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
127785c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
127847d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
127947d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
12804fad6a16SStefano Zampini   pcbddc->current_level              = 0;
12812b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1282da1bb401SStefano Zampini 
1283674ae819SStefano Zampini   /* create local graph structure */
1284674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1285674ae819SStefano Zampini 
1286674ae819SStefano Zampini   /* scaling */
1287674ae819SStefano Zampini   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1288674ae819SStefano Zampini   pcbddc->work_scaling               = 0;
1289da1bb401SStefano Zampini 
1290da1bb401SStefano Zampini   /* function pointers */
1291da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
1292da1bb401SStefano Zampini   pc->ops->applytranspose      = 0;
1293da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1294da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1295da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1296da1bb401SStefano Zampini   pc->ops->view                = 0;
1297da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1298da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1299da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1300534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1301534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1302da1bb401SStefano Zampini 
1303da1bb401SStefano Zampini   /* composing function */
1304674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1305bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
13062b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1307b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
13082b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1309bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1310bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1311bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1312bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1313bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1314bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1315bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1316bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1317bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1318bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1319da1bb401SStefano Zampini   PetscFunctionReturn(0);
1320da1bb401SStefano Zampini }
13213425bc38SStefano Zampini 
1322