xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision fcd409f55db887dc8ef901cc2123e1395f10ef8d)
153cdbc3dSStefano Zampini /* TODOLIST
2b8ffe317SStefano Zampini    Remove PetscOptionsSetValue for pcis!
3da1bb401SStefano Zampini    DofSplitting and DM attached to pc?
4da1bb401SStefano Zampini    Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet)
5b8ffe317SStefano Zampini    BDDC with MG framework?
6a0ba757dSStefano Zampini    provide other ops? Ask to developers
7a0ba757dSStefano Zampini    man pages
853cdbc3dSStefano Zampini */
90c7d97c5SJed Brown 
1053cdbc3dSStefano Zampini /* ----------------------------------------------------------------------------------------------------------------------------------------------
110c7d97c5SJed Brown    Implementation of BDDC preconditioner based on:
120c7d97c5SJed Brown    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
1353cdbc3dSStefano Zampini    ---------------------------------------------------------------------------------------------------------------------------------------------- */
1453cdbc3dSStefano Zampini 
15674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/  /* includes for fortran wrappers */
16674ae819SStefano Zampini #include "bddcprivate.h"
173b03a366Sstefano_zampini #include <petscblaslapack.h>
18674ae819SStefano Zampini 
190c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
200c7d97c5SJed Brown #undef __FUNCT__
210c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
220c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc)
230c7d97c5SJed Brown {
240c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
250c7d97c5SJed Brown   PetscErrorCode ierr;
260c7d97c5SJed Brown 
270c7d97c5SJed Brown   PetscFunctionBegin;
280c7d97c5SJed Brown   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
298eeda7d8SStefano Zampini   /* Verbose debugging */
308eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
318eeda7d8SStefano Zampini   /* Primal space cumstomization */
328eeda7d8SStefano 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);
338eeda7d8SStefano 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);
348eeda7d8SStefano 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);
358eeda7d8SStefano Zampini   /* Change of basis */
368eeda7d8SStefano 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);
378eeda7d8SStefano 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);
38674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
39674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
40674ae819SStefano Zampini   }
418eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
428eeda7d8SStefano 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);
430298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
442b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
45674ae819SStefano 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);
460c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
470c7d97c5SJed Brown   PetscFunctionReturn(0);
480c7d97c5SJed Brown }
490c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
50674ae819SStefano Zampini #undef __FUNCT__
51674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
52674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
53674ae819SStefano Zampini {
54674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
55674ae819SStefano Zampini   PetscErrorCode ierr;
561e6b0712SBarry Smith 
57674ae819SStefano Zampini   PetscFunctionBegin;
58674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
59674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
60674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
61674ae819SStefano Zampini   PetscFunctionReturn(0);
62674ae819SStefano Zampini }
63674ae819SStefano Zampini #undef __FUNCT__
64674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
65674ae819SStefano Zampini /*@
66674ae819SStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set user defined primal vertices in PCBDDC.
67674ae819SStefano Zampini 
68674ae819SStefano Zampini    Not collective
69674ae819SStefano Zampini 
70674ae819SStefano Zampini    Input Parameters:
71674ae819SStefano Zampini +  pc - the preconditioning context
72674ae819SStefano Zampini -  PrimalVertices - index sets of primal vertices in local numbering
73674ae819SStefano Zampini 
74674ae819SStefano Zampini    Level: intermediate
75674ae819SStefano Zampini 
76674ae819SStefano Zampini    Notes:
77674ae819SStefano Zampini 
78674ae819SStefano Zampini .seealso: PCBDDC
79674ae819SStefano Zampini @*/
80674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
81674ae819SStefano Zampini {
82674ae819SStefano Zampini   PetscErrorCode ierr;
83674ae819SStefano Zampini 
84674ae819SStefano Zampini   PetscFunctionBegin;
85674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
86674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
87674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
88674ae819SStefano Zampini   PetscFunctionReturn(0);
89674ae819SStefano Zampini }
90674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
910c7d97c5SJed Brown #undef __FUNCT__
924fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
934fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
944fad6a16SStefano Zampini {
954fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
964fad6a16SStefano Zampini 
974fad6a16SStefano Zampini   PetscFunctionBegin;
984fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
994fad6a16SStefano Zampini   PetscFunctionReturn(0);
1004fad6a16SStefano Zampini }
1011e6b0712SBarry Smith 
1024fad6a16SStefano Zampini #undef __FUNCT__
1034fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1044fad6a16SStefano Zampini /*@
1054fad6a16SStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel coarsening
1064fad6a16SStefano Zampini 
1074fad6a16SStefano Zampini    Logically collective on PC
1084fad6a16SStefano Zampini 
1094fad6a16SStefano Zampini    Input Parameters:
1104fad6a16SStefano Zampini +  pc - the preconditioning context
1114fad6a16SStefano Zampini -  k - coarsening ratio
1124fad6a16SStefano Zampini 
1134fad6a16SStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level.
1144fad6a16SStefano Zampini 
1154fad6a16SStefano Zampini    Level: intermediate
1164fad6a16SStefano Zampini 
1174fad6a16SStefano Zampini    Notes:
1184fad6a16SStefano Zampini 
1194fad6a16SStefano Zampini .seealso: PCBDDC
1204fad6a16SStefano Zampini @*/
1214fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
1224fad6a16SStefano Zampini {
1234fad6a16SStefano Zampini   PetscErrorCode ierr;
1244fad6a16SStefano Zampini 
1254fad6a16SStefano Zampini   PetscFunctionBegin;
1264fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1272b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
1284fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
1294fad6a16SStefano Zampini   PetscFunctionReturn(0);
1304fad6a16SStefano Zampini }
1312b510759SStefano Zampini 
132b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
1332b510759SStefano Zampini #undef __FUNCT__
134b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
135b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
136b8ffe317SStefano Zampini {
137b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
138b8ffe317SStefano Zampini 
139b8ffe317SStefano Zampini   PetscFunctionBegin;
140b8ffe317SStefano Zampini   pcbddc->use_exact_dirichlet = flg;
141b8ffe317SStefano Zampini   PetscFunctionReturn(0);
142b8ffe317SStefano Zampini }
143b8ffe317SStefano Zampini 
144b8ffe317SStefano Zampini #undef __FUNCT__
145b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
146b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
1472b510759SStefano Zampini {
1482b510759SStefano Zampini   PetscErrorCode ierr;
1492b510759SStefano Zampini 
1502b510759SStefano Zampini   PetscFunctionBegin;
1512b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
152b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
153b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
1542b510759SStefano Zampini   PetscFunctionReturn(0);
1552b510759SStefano Zampini }
1561e6b0712SBarry Smith 
1574fad6a16SStefano Zampini #undef __FUNCT__
1582b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
1592b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
1604fad6a16SStefano Zampini {
1614fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1624fad6a16SStefano Zampini 
1634fad6a16SStefano Zampini   PetscFunctionBegin;
1642b510759SStefano Zampini   pcbddc->current_level = level;
1654fad6a16SStefano Zampini   PetscFunctionReturn(0);
1664fad6a16SStefano Zampini }
1671e6b0712SBarry Smith 
1684fad6a16SStefano Zampini #undef __FUNCT__
169b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
170b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
171b8ffe317SStefano Zampini {
172b8ffe317SStefano Zampini   PetscErrorCode ierr;
173b8ffe317SStefano Zampini 
174b8ffe317SStefano Zampini   PetscFunctionBegin;
175b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
176b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
177b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
178b8ffe317SStefano Zampini   PetscFunctionReturn(0);
179b8ffe317SStefano Zampini }
180b8ffe317SStefano Zampini 
181b8ffe317SStefano Zampini #undef __FUNCT__
1822b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
1832b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
1842b510759SStefano Zampini {
1852b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1862b510759SStefano Zampini 
1872b510759SStefano Zampini   PetscFunctionBegin;
1882b510759SStefano Zampini   pcbddc->max_levels = levels;
1892b510759SStefano Zampini   PetscFunctionReturn(0);
1902b510759SStefano Zampini }
1912b510759SStefano Zampini 
1922b510759SStefano Zampini #undef __FUNCT__
1932b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
1944fad6a16SStefano Zampini /*@
1952b510759SStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels within the multilevel approach.
1964fad6a16SStefano Zampini 
1974fad6a16SStefano Zampini    Logically collective on PC
1984fad6a16SStefano Zampini 
1994fad6a16SStefano Zampini    Input Parameters:
2004fad6a16SStefano Zampini +  pc - the preconditioning context
2012b510759SStefano Zampini -  levels - the maximum number of levels
2024fad6a16SStefano Zampini 
2032b510759SStefano Zampini    Default value is 0, i.e. coarse problem will be solved exactly
2044fad6a16SStefano Zampini 
2054fad6a16SStefano Zampini    Level: intermediate
2064fad6a16SStefano Zampini 
2074fad6a16SStefano Zampini    Notes:
2084fad6a16SStefano Zampini 
2094fad6a16SStefano Zampini .seealso: PCBDDC
2104fad6a16SStefano Zampini @*/
2112b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
2124fad6a16SStefano Zampini {
2134fad6a16SStefano Zampini   PetscErrorCode ierr;
2144fad6a16SStefano Zampini 
2154fad6a16SStefano Zampini   PetscFunctionBegin;
2164fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2172b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
2182b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
2194fad6a16SStefano Zampini   PetscFunctionReturn(0);
2204fad6a16SStefano Zampini }
2214fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
2221e6b0712SBarry Smith 
2234fad6a16SStefano Zampini #undef __FUNCT__
2240bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
2250bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
2260bdf917eSStefano Zampini {
2270bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2280bdf917eSStefano Zampini   PetscErrorCode ierr;
2290bdf917eSStefano Zampini 
2300bdf917eSStefano Zampini   PetscFunctionBegin;
2310bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
2320bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
2330bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
2340bdf917eSStefano Zampini   PetscFunctionReturn(0);
2350bdf917eSStefano Zampini }
2361e6b0712SBarry Smith 
2370bdf917eSStefano Zampini #undef __FUNCT__
2380bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
2390bdf917eSStefano Zampini /*@
2400bdf917eSStefano Zampini  PCBDDCSetNullSpace - Set NullSpace of global operator of BDDC preconditioned mat.
2410bdf917eSStefano Zampini 
2420bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
2430bdf917eSStefano Zampini 
2440bdf917eSStefano Zampini    Input Parameters:
2450bdf917eSStefano Zampini +  pc - the preconditioning context
2460bdf917eSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned.
2470bdf917eSStefano Zampini 
2480bdf917eSStefano Zampini    Level: intermediate
2490bdf917eSStefano Zampini 
2500bdf917eSStefano Zampini    Notes:
2510bdf917eSStefano Zampini 
2520bdf917eSStefano Zampini .seealso: PCBDDC
2530bdf917eSStefano Zampini @*/
2540bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
2550bdf917eSStefano Zampini {
2560bdf917eSStefano Zampini   PetscErrorCode ierr;
2570bdf917eSStefano Zampini 
2580bdf917eSStefano Zampini   PetscFunctionBegin;
2590bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
260674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
2612b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
2620bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
2630bdf917eSStefano Zampini   PetscFunctionReturn(0);
2640bdf917eSStefano Zampini }
2650bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
2661e6b0712SBarry Smith 
2670bdf917eSStefano Zampini #undef __FUNCT__
2683b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
2693b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
2703b03a366Sstefano_zampini {
2713b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2723b03a366Sstefano_zampini   PetscErrorCode ierr;
2733b03a366Sstefano_zampini 
2743b03a366Sstefano_zampini   PetscFunctionBegin;
2753b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
27636e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
27736e030ebSStefano Zampini   pcbddc->DirichletBoundaries=DirichletBoundaries;
2783b03a366Sstefano_zampini   PetscFunctionReturn(0);
2793b03a366Sstefano_zampini }
2801e6b0712SBarry Smith 
2813b03a366Sstefano_zampini #undef __FUNCT__
2823b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
2833b03a366Sstefano_zampini /*@
284da1bb401SStefano Zampini  PCBDDCSetDirichletBoundaries - Set index set defining subdomain part (in local ordering)
285da1bb401SStefano Zampini                               of Dirichlet boundaries for the global problem.
2863b03a366Sstefano_zampini 
2873b03a366Sstefano_zampini    Not collective
2883b03a366Sstefano_zampini 
2893b03a366Sstefano_zampini    Input Parameters:
2903b03a366Sstefano_zampini +  pc - the preconditioning context
2910298fd71SBarry Smith -  DirichletBoundaries - sequential index set defining the subdomain part of Dirichlet boundaries (can be NULL)
2923b03a366Sstefano_zampini 
2933b03a366Sstefano_zampini    Level: intermediate
2943b03a366Sstefano_zampini 
2953b03a366Sstefano_zampini    Notes:
2963b03a366Sstefano_zampini 
2973b03a366Sstefano_zampini .seealso: PCBDDC
2983b03a366Sstefano_zampini @*/
2993b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3003b03a366Sstefano_zampini {
3013b03a366Sstefano_zampini   PetscErrorCode ierr;
3023b03a366Sstefano_zampini 
3033b03a366Sstefano_zampini   PetscFunctionBegin;
3043b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
305674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
3063b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
3073b03a366Sstefano_zampini   PetscFunctionReturn(0);
3083b03a366Sstefano_zampini }
3093b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
3101e6b0712SBarry Smith 
3113b03a366Sstefano_zampini #undef __FUNCT__
3120c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
31353cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
3140c7d97c5SJed Brown {
3150c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
31653cdbc3dSStefano Zampini   PetscErrorCode ierr;
3170c7d97c5SJed Brown 
3180c7d97c5SJed Brown   PetscFunctionBegin;
31953cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
32036e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
32136e030ebSStefano Zampini   pcbddc->NeumannBoundaries=NeumannBoundaries;
3220c7d97c5SJed Brown   PetscFunctionReturn(0);
3230c7d97c5SJed Brown }
3241e6b0712SBarry Smith 
3250c7d97c5SJed Brown #undef __FUNCT__
3260c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
32757527edcSJed Brown /*@
328da1bb401SStefano Zampini  PCBDDCSetNeumannBoundaries - Set index set defining subdomain part (in local ordering)
329da1bb401SStefano Zampini                               of Neumann boundaries for the global problem.
33057527edcSJed Brown 
3319c0446d6SStefano Zampini    Not collective
33257527edcSJed Brown 
33357527edcSJed Brown    Input Parameters:
33457527edcSJed Brown +  pc - the preconditioning context
3350298fd71SBarry Smith -  NeumannBoundaries - sequential index set defining the subdomain part of Neumann boundaries (can be NULL)
33657527edcSJed Brown 
33757527edcSJed Brown    Level: intermediate
33857527edcSJed Brown 
33957527edcSJed Brown    Notes:
34057527edcSJed Brown 
34157527edcSJed Brown .seealso: PCBDDC
34257527edcSJed Brown @*/
34353cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
3440c7d97c5SJed Brown {
3450c7d97c5SJed Brown   PetscErrorCode ierr;
3460c7d97c5SJed Brown 
3470c7d97c5SJed Brown   PetscFunctionBegin;
3480c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
349674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
35053cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
35153cdbc3dSStefano Zampini   PetscFunctionReturn(0);
35253cdbc3dSStefano Zampini }
35353cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
3541e6b0712SBarry Smith 
35553cdbc3dSStefano Zampini #undef __FUNCT__
356da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
357da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
358da1bb401SStefano Zampini {
359da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
360da1bb401SStefano Zampini 
361da1bb401SStefano Zampini   PetscFunctionBegin;
362da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
363da1bb401SStefano Zampini   PetscFunctionReturn(0);
364da1bb401SStefano Zampini }
3651e6b0712SBarry Smith 
366da1bb401SStefano Zampini #undef __FUNCT__
367da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
368da1bb401SStefano Zampini /*@
369da1bb401SStefano Zampini  PCBDDCGetDirichletBoundaries - Get index set defining subdomain part (in local ordering)
370da1bb401SStefano Zampini                                 of Dirichlet boundaries for the global problem.
371da1bb401SStefano Zampini 
372da1bb401SStefano Zampini    Not collective
373da1bb401SStefano Zampini 
374da1bb401SStefano Zampini    Input Parameters:
375da1bb401SStefano Zampini +  pc - the preconditioning context
376da1bb401SStefano Zampini 
377da1bb401SStefano Zampini    Output Parameters:
378da1bb401SStefano Zampini +  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
379da1bb401SStefano Zampini 
380da1bb401SStefano Zampini    Level: intermediate
381da1bb401SStefano Zampini 
382da1bb401SStefano Zampini    Notes:
383da1bb401SStefano Zampini 
384da1bb401SStefano Zampini .seealso: PCBDDC
385da1bb401SStefano Zampini @*/
386da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
387da1bb401SStefano Zampini {
388da1bb401SStefano Zampini   PetscErrorCode ierr;
389da1bb401SStefano Zampini 
390da1bb401SStefano Zampini   PetscFunctionBegin;
391da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
392da1bb401SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
393da1bb401SStefano Zampini   PetscFunctionReturn(0);
394da1bb401SStefano Zampini }
395da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
3961e6b0712SBarry Smith 
397da1bb401SStefano Zampini #undef __FUNCT__
39853cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
39953cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
40053cdbc3dSStefano Zampini {
40153cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
40253cdbc3dSStefano Zampini 
40353cdbc3dSStefano Zampini   PetscFunctionBegin;
40453cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
40553cdbc3dSStefano Zampini   PetscFunctionReturn(0);
40653cdbc3dSStefano Zampini }
4071e6b0712SBarry Smith 
40853cdbc3dSStefano Zampini #undef __FUNCT__
40953cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
41053cdbc3dSStefano Zampini /*@
411da1bb401SStefano Zampini  PCBDDCGetNeumannBoundaries - Get index set defining subdomain part (in local ordering)
412da1bb401SStefano Zampini                               of Neumann boundaries for the global problem.
41353cdbc3dSStefano Zampini 
4149c0446d6SStefano Zampini    Not collective
41553cdbc3dSStefano Zampini 
41653cdbc3dSStefano Zampini    Input Parameters:
41753cdbc3dSStefano Zampini +  pc - the preconditioning context
41853cdbc3dSStefano Zampini 
41953cdbc3dSStefano Zampini    Output Parameters:
42053cdbc3dSStefano Zampini +  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
42153cdbc3dSStefano Zampini 
42253cdbc3dSStefano Zampini    Level: intermediate
42353cdbc3dSStefano Zampini 
42453cdbc3dSStefano Zampini    Notes:
42553cdbc3dSStefano Zampini 
42653cdbc3dSStefano Zampini .seealso: PCBDDC
42753cdbc3dSStefano Zampini @*/
42853cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
42953cdbc3dSStefano Zampini {
43053cdbc3dSStefano Zampini   PetscErrorCode ierr;
43153cdbc3dSStefano Zampini 
43253cdbc3dSStefano Zampini   PetscFunctionBegin;
43353cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
43453cdbc3dSStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
4350c7d97c5SJed Brown   PetscFunctionReturn(0);
4360c7d97c5SJed Brown }
43736e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
4381e6b0712SBarry Smith 
43936e030ebSStefano Zampini #undef __FUNCT__
440da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
4411a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
44236e030ebSStefano Zampini {
44336e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
444da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
445da1bb401SStefano Zampini   PetscErrorCode ierr;
44636e030ebSStefano Zampini 
44736e030ebSStefano Zampini   PetscFunctionBegin;
448674ae819SStefano Zampini   /* free old CSR */
449674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
450fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
451674ae819SStefano Zampini   /* get CSR into graph structure */
452da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
453674ae819SStefano Zampini     ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr);
454674ae819SStefano Zampini     ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr);
455674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
456674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
457da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
4581a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
4591a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
460674ae819SStefano Zampini   }
461575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
46236e030ebSStefano Zampini   PetscFunctionReturn(0);
46336e030ebSStefano Zampini }
4641e6b0712SBarry Smith 
46536e030ebSStefano Zampini #undef __FUNCT__
466da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
46736e030ebSStefano Zampini /*@
468da1bb401SStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set CSR graph of local matrix for use of PCBDDC.
46936e030ebSStefano Zampini 
47036e030ebSStefano Zampini    Not collective
47136e030ebSStefano Zampini 
47236e030ebSStefano Zampini    Input Parameters:
47336e030ebSStefano Zampini +  pc - the preconditioning context
474da1bb401SStefano Zampini -  nvtxs - number of local vertices of the graph
475da1bb401SStefano Zampini -  xadj, adjncy - the CSR graph
476da1bb401SStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. In the former case the user must free the array passed in;
477da1bb401SStefano Zampini                                                              in the latter case, memory must be obtained with PetscMalloc.
47836e030ebSStefano Zampini 
47936e030ebSStefano Zampini    Level: intermediate
48036e030ebSStefano Zampini 
48136e030ebSStefano Zampini    Notes:
48236e030ebSStefano Zampini 
48336e030ebSStefano Zampini .seealso: PCBDDC
48436e030ebSStefano Zampini @*/
4851a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
48636e030ebSStefano Zampini {
487575ad6abSStefano Zampini   void (*f)(void) = 0;
48836e030ebSStefano Zampini   PetscErrorCode ierr;
48936e030ebSStefano Zampini 
49036e030ebSStefano Zampini   PetscFunctionBegin;
49136e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
492674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
493674ae819SStefano Zampini   PetscValidIntPointer(xadj,4);
494674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
495674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
496da1bb401SStefano Zampini   }
49736e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
498575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
499575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
500575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
501575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
502575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
50336e030ebSStefano Zampini   }
50436e030ebSStefano Zampini   PetscFunctionReturn(0);
50536e030ebSStefano Zampini }
5069c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
5071e6b0712SBarry Smith 
5089c0446d6SStefano Zampini #undef __FUNCT__
5099c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
5109c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
5119c0446d6SStefano Zampini {
5129c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5139c0446d6SStefano Zampini   PetscInt i;
5149c0446d6SStefano Zampini   PetscErrorCode ierr;
5159c0446d6SStefano Zampini 
5169c0446d6SStefano Zampini   PetscFunctionBegin;
517da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
5189c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
5199c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
5209c0446d6SStefano Zampini   }
521d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
522da1bb401SStefano Zampini   /* allocate space then set */
5239c0446d6SStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr);
5249c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
525da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
526da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
5279c0446d6SStefano Zampini   }
5289c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
5299c0446d6SStefano Zampini   PetscFunctionReturn(0);
5309c0446d6SStefano Zampini }
5311e6b0712SBarry Smith 
5329c0446d6SStefano Zampini #undef __FUNCT__
5339c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
5349c0446d6SStefano Zampini /*@
535da1bb401SStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of local mat.
5369c0446d6SStefano Zampini 
5379c0446d6SStefano Zampini    Not collective
5389c0446d6SStefano Zampini 
5399c0446d6SStefano Zampini    Input Parameters:
5409c0446d6SStefano Zampini +  pc - the preconditioning context
541da1bb401SStefano Zampini -  n - number of index sets defining the fields
542da1bb401SStefano Zampini -  IS[] - array of IS describing the fields
5439c0446d6SStefano Zampini 
5449c0446d6SStefano Zampini    Level: intermediate
5459c0446d6SStefano Zampini 
5469c0446d6SStefano Zampini    Notes:
5479c0446d6SStefano Zampini 
5489c0446d6SStefano Zampini .seealso: PCBDDC
5499c0446d6SStefano Zampini @*/
5509c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
5519c0446d6SStefano Zampini {
5522b510759SStefano Zampini   PetscInt       i;
5539c0446d6SStefano Zampini   PetscErrorCode ierr;
5549c0446d6SStefano Zampini 
5559c0446d6SStefano Zampini   PetscFunctionBegin;
5569c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5572b510759SStefano Zampini   for (i=0;i<n_is;i++) {
5582b510759SStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2);
5592b510759SStefano Zampini   }
5609c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
5619c0446d6SStefano Zampini   PetscFunctionReturn(0);
5629c0446d6SStefano Zampini }
563da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
564534831adSStefano Zampini #undef __FUNCT__
565534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
566534831adSStefano Zampini /* -------------------------------------------------------------------------- */
567534831adSStefano Zampini /*
568534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
569534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
5709c0446d6SStefano Zampini 
571534831adSStefano Zampini    Input Parameter:
572534831adSStefano Zampini +  pc - the preconditioner contex
573534831adSStefano Zampini 
574534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
575534831adSStefano Zampini 
576534831adSStefano Zampini    Notes:
577534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
578534831adSStefano Zampini    the user, but instead is called by KSPSolve().
579534831adSStefano Zampini */
580534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
581534831adSStefano Zampini {
582534831adSStefano Zampini   PetscErrorCode ierr;
583534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
584534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
585534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
586534831adSStefano Zampini   Mat            temp_mat;
5873972b0daSStefano Zampini   IS             dirIS;
5883972b0daSStefano Zampini   PetscInt       dirsize,i,*is_indices;
5893972b0daSStefano Zampini   PetscScalar    *array_x,*array_diagonal;
5903972b0daSStefano Zampini   Vec            used_vec;
5913972b0daSStefano Zampini   PetscBool      guess_nonzero;
592534831adSStefano Zampini 
593534831adSStefano Zampini   PetscFunctionBegin;
59462a6ff1dSStefano Zampini   /* Creates parallel work vectors used in presolve. */
59562a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
59662a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
59762a6ff1dSStefano Zampini   }
59862a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
59962a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
60062a6ff1dSStefano Zampini   }
6013972b0daSStefano Zampini   if (x) {
6023972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
6033972b0daSStefano Zampini     used_vec = x;
6043972b0daSStefano Zampini   } else {
6053972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
6063972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
6073972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6083972b0daSStefano Zampini   }
6093972b0daSStefano Zampini   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
6103972b0daSStefano Zampini   if (ksp) {
6113972b0daSStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
6123972b0daSStefano Zampini     if ( !guess_nonzero ) {
6133972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
6143972b0daSStefano Zampini     }
6153972b0daSStefano Zampini   }
6163308cffdSStefano Zampini 
61762a6ff1dSStefano Zampini   if (rhs) { /* TODO: wiser handling of rhs removal, which is only needed in case of zeroed rows */
6183972b0daSStefano Zampini     /* store the original rhs */
6193972b0daSStefano Zampini     ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
6203972b0daSStefano Zampini 
6213972b0daSStefano Zampini     /* Take into account zeroed rows -> change rhs and store solution removed */
6223972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
6233972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
6243972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6253972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6263972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6273972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6283972b0daSStefano Zampini     ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr);
6293972b0daSStefano Zampini     if (dirIS) {
6303972b0daSStefano Zampini       ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr);
6313972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6323972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6333972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6342fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
6353972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
6363972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
6373972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
6383972b0daSStefano Zampini     }
6393972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6403972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
641b76ba322SStefano Zampini 
6423972b0daSStefano Zampini     /* remove the computed solution from the rhs */
6433972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6443972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
6453972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
6463308cffdSStefano Zampini   }
647b76ba322SStefano Zampini 
648b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
6493972b0daSStefano Zampini   if (x) {
6503972b0daSStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
6513972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
65215aaf578SStefano Zampini     if (pcbddc->use_exact_dirichlet && !pcbddc->coarse_psi_B) {
653b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
654b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
655b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
656b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
657b76ba322SStefano Zampini       ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
658b76ba322SStefano Zampini       if (ksp) {
659b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
660b76ba322SStefano Zampini       }
661b76ba322SStefano Zampini     }
6623972b0daSStefano Zampini   }
663b76ba322SStefano Zampini 
664674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
665b76ba322SStefano Zampini     /* swap pointers for local matrices */
666b76ba322SStefano Zampini     temp_mat = matis->A;
667b76ba322SStefano Zampini     matis->A = pcbddc->local_mat;
668b76ba322SStefano Zampini     pcbddc->local_mat = temp_mat;
6693308cffdSStefano Zampini   }
6703308cffdSStefano Zampini   if (pcbddc->use_change_of_basis && rhs) {
671b76ba322SStefano Zampini     /* Get local rhs and apply transformation of basis */
672b76ba322SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
673b76ba322SStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
674b76ba322SStefano Zampini     /* from original basis to modified basis */
675b76ba322SStefano Zampini     ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
676b76ba322SStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
677b76ba322SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
678b76ba322SStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
679674ae819SStefano Zampini   }
6800bdf917eSStefano Zampini   if (ksp && pcbddc->NullSpace) {
681d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
682d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
683b76ba322SStefano Zampini   }
6840bdf917eSStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
685534831adSStefano Zampini   PetscFunctionReturn(0);
686534831adSStefano Zampini }
687534831adSStefano Zampini /* -------------------------------------------------------------------------- */
688534831adSStefano Zampini #undef __FUNCT__
689534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
690534831adSStefano Zampini /* -------------------------------------------------------------------------- */
691534831adSStefano Zampini /*
692534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
693534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
694534831adSStefano Zampini 
695534831adSStefano Zampini    Input Parameter:
696534831adSStefano Zampini +  pc - the preconditioner contex
697534831adSStefano Zampini 
698534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
699534831adSStefano Zampini 
700534831adSStefano Zampini    Notes:
701534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
702534831adSStefano Zampini    the user, but instead is called by KSPSolve().
703534831adSStefano Zampini */
704534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
705534831adSStefano Zampini {
706534831adSStefano Zampini   PetscErrorCode ierr;
707534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
708534831adSStefano Zampini   PC_IS          *pcis   = (PC_IS*)(pc->data);
709534831adSStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
710534831adSStefano Zampini   Mat            temp_mat;
711534831adSStefano Zampini 
712534831adSStefano Zampini   PetscFunctionBegin;
713674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
714534831adSStefano Zampini     /* swap pointers for local matrices */
715534831adSStefano Zampini     temp_mat = matis->A;
716534831adSStefano Zampini     matis->A = pcbddc->local_mat;
717534831adSStefano Zampini     pcbddc->local_mat = temp_mat;
7183425bc38SStefano Zampini   }
7193308cffdSStefano Zampini   if (pcbddc->use_change_of_basis && x) {
720534831adSStefano Zampini     /* Get Local boundary and apply transformation of basis to solution vector */
721534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
722534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
723534831adSStefano Zampini     /* from modified basis to original basis */
724534831adSStefano Zampini     ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr);
725534831adSStefano Zampini     /* put back modified values into the global vec using INSERT_VALUES copy mode */
726534831adSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
727534831adSStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
728534831adSStefano Zampini   }
7293972b0daSStefano Zampini   /* add solution removed in presolve */
7303425bc38SStefano Zampini   if (x) {
7313425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
7323425bc38SStefano Zampini   }
733fb223d50SStefano Zampini   /* restore rhs to its original state */
734fb223d50SStefano Zampini   if (rhs) {
735fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
736fb223d50SStefano Zampini   }
737534831adSStefano Zampini   PetscFunctionReturn(0);
738534831adSStefano Zampini }
739534831adSStefano Zampini /* -------------------------------------------------------------------------- */
74053cdbc3dSStefano Zampini #undef __FUNCT__
74153cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
7420c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
7430c7d97c5SJed Brown /*
7440c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
7450c7d97c5SJed Brown                   by setting data structures and options.
7460c7d97c5SJed Brown 
7470c7d97c5SJed Brown    Input Parameter:
74853cdbc3dSStefano Zampini +  pc - the preconditioner context
7490c7d97c5SJed Brown 
7500c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
7510c7d97c5SJed Brown 
7520c7d97c5SJed Brown    Notes:
7530c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
7540c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
7550c7d97c5SJed Brown */
75653cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
7570c7d97c5SJed Brown {
7580c7d97c5SJed Brown   PetscErrorCode ierr;
7590c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
760674ae819SStefano Zampini   MatStructure   flag;
761674ae819SStefano Zampini   PetscBool      computeis,computetopography,computesolvers;
7620c7d97c5SJed Brown 
7630c7d97c5SJed Brown   PetscFunctionBegin;
764674ae819SStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */
765*fcd409f5SStefano Zampini   /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */
7663b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
767*fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
7683b03a366Sstefano_zampini   /* Get stdout for dbg */
769674ae819SStefano Zampini   if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) {
770ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr);
771e269702eSStefano Zampini     ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
7722b510759SStefano Zampini     if (pcbddc->current_level) {
7732b510759SStefano Zampini       ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
7742b510759SStefano Zampini     }
775e269702eSStefano Zampini   }
776674ae819SStefano Zampini   /* first attempt to split work */
777674ae819SStefano Zampini   if (pc->setupcalled) {
778674ae819SStefano Zampini     computeis = PETSC_FALSE;
779674ae819SStefano Zampini     ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr);
780674ae819SStefano Zampini     if (flag == SAME_PRECONDITIONER) {
781674ae819SStefano Zampini       computetopography = PETSC_FALSE;
782674ae819SStefano Zampini       computesolvers = PETSC_FALSE;
783674ae819SStefano Zampini     } else if (flag == SAME_NONZERO_PATTERN) {
784674ae819SStefano Zampini       computetopography = PETSC_FALSE;
785674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
786674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
787674ae819SStefano Zampini       computetopography = PETSC_TRUE;
788674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
789674ae819SStefano Zampini     }
790674ae819SStefano Zampini   } else {
791674ae819SStefano Zampini     computeis = PETSC_TRUE;
792674ae819SStefano Zampini     computetopography = PETSC_TRUE;
793674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
794674ae819SStefano Zampini   }
795*fcd409f5SStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
796674ae819SStefano Zampini   if (computeis) {
797*fcd409f5SStefano Zampini     /* HACK INTO PCIS */
798*fcd409f5SStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
799*fcd409f5SStefano Zampini     pcis->computesolvers = PETSC_FALSE;
800674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
801674ae819SStefano Zampini   }
802674ae819SStefano Zampini   /* Analyze interface and set up local constraint and change of basis matrices */
803674ae819SStefano Zampini   if (computetopography) {
804674ae819SStefano Zampini     /* reset data */
805674ae819SStefano Zampini     ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
806674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
807674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
808674ae819SStefano Zampini   }
809674ae819SStefano Zampini   if (computesolvers) {
810674ae819SStefano Zampini     /* reset data */
811674ae819SStefano Zampini     ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
812674ae819SStefano Zampini     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
81399cc7994SStefano Zampini     /* Create coarse and local stuffs */
81499cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
815674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
8160c7d97c5SJed Brown   }
8172b510759SStefano Zampini   if (pcbddc->dbg_flag && pcbddc->current_level) {
8182b510759SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr);
8192b510759SStefano Zampini   }
8200c7d97c5SJed Brown   PetscFunctionReturn(0);
8210c7d97c5SJed Brown }
8220c7d97c5SJed Brown 
8230c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
8240c7d97c5SJed Brown /*
8250c7d97c5SJed Brown    PCApply_BDDC - Applies the BDDC preconditioner to a vector.
8260c7d97c5SJed Brown 
8270c7d97c5SJed Brown    Input Parameters:
8280c7d97c5SJed Brown .  pc - the preconditioner context
8290c7d97c5SJed Brown .  r - input vector (global)
8300c7d97c5SJed Brown 
8310c7d97c5SJed Brown    Output Parameter:
8320c7d97c5SJed Brown .  z - output vector (global)
8330c7d97c5SJed Brown 
8340c7d97c5SJed Brown    Application Interface Routine: PCApply()
8350c7d97c5SJed Brown  */
8360c7d97c5SJed Brown #undef __FUNCT__
8370c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
83853cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
8390c7d97c5SJed Brown {
8400c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
8410c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
8420c7d97c5SJed Brown   PetscErrorCode    ierr;
8433b03a366Sstefano_zampini   const PetscScalar one = 1.0;
8443b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
8452617d88aSStefano Zampini   const PetscScalar zero = 0.0;
8460c7d97c5SJed Brown 
8470c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
8480c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
8498eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
8500c7d97c5SJed Brown 
8510c7d97c5SJed Brown   PetscFunctionBegin;
85215aaf578SStefano Zampini   if (!pcbddc->use_exact_dirichlet || pcbddc->coarse_psi_B) {
8530c7d97c5SJed Brown     /* First Dirichlet solve */
8540c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8550c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
85653cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
8570c7d97c5SJed Brown     /*
8580c7d97c5SJed Brown       Assembling right hand side for BDDC operator
859674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
860674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
8610c7d97c5SJed Brown     */
8620c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
8630c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
8648eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
8650c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
8660c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
8670c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8680c7d97c5SJed Brown     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
869674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
870b76ba322SStefano Zampini   } else {
8710bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
872b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
873674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
874b76ba322SStefano Zampini   }
875b76ba322SStefano Zampini 
8762617d88aSStefano Zampini   /* Apply interface preconditioner
8772617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
8782617d88aSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr);
8792617d88aSStefano Zampini 
880674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
881674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
8820c7d97c5SJed Brown 
8833b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
8840c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8850c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8860c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
8878eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
88853cdbc3dSStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr);
8890c7d97c5SJed Brown   ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr);
8908eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
8910c7d97c5SJed Brown   ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr);
8920c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8930c7d97c5SJed Brown   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8940c7d97c5SJed Brown   PetscFunctionReturn(0);
8950c7d97c5SJed Brown }
896da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
897674ae819SStefano Zampini 
898da1bb401SStefano Zampini #undef __FUNCT__
899da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
900da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
901da1bb401SStefano Zampini {
902da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
903da1bb401SStefano Zampini   PetscErrorCode ierr;
904da1bb401SStefano Zampini 
905da1bb401SStefano Zampini   PetscFunctionBegin;
906da1bb401SStefano Zampini   /* free data created by PCIS */
907da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
908674ae819SStefano Zampini   /* free BDDC custom data  */
909674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
910674ae819SStefano Zampini   /* destroy objects related to topography */
911674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
912674ae819SStefano Zampini   /* free allocated graph structure */
913da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
914674ae819SStefano Zampini   /* free data for scaling operator */
915674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
916674ae819SStefano Zampini   /* free solvers stuff */
917674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
91833bc96a4SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
91933bc96a4SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
9209881197aSStefano Zampini   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
921ac78edfcSStefano Zampini   ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
92262a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
92362a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
92462a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
9253425bc38SStefano Zampini   /* remove functions */
926674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
927bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
9282b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
929b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
9302b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
931bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
932bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
933bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
934bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
935bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
936bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
937bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
938bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
939bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
940bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
941674ae819SStefano Zampini   /* Free the private data structure */
942674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
943da1bb401SStefano Zampini   PetscFunctionReturn(0);
944da1bb401SStefano Zampini }
9453425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
9461e6b0712SBarry Smith 
9473425bc38SStefano Zampini #undef __FUNCT__
9483425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
9493425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
9503425bc38SStefano Zampini {
951674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
9523425bc38SStefano Zampini   PC_IS*         pcis;
9533425bc38SStefano Zampini   PC_BDDC*       pcbddc;
9543425bc38SStefano Zampini   PetscErrorCode ierr;
9550c7d97c5SJed Brown 
9563425bc38SStefano Zampini   PetscFunctionBegin;
9573425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
9583425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
9593425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
9603425bc38SStefano Zampini 
9613425bc38SStefano Zampini   /* change of basis for physical rhs if needed
9623425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
9633308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
9643425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
9653425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9663425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
967fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
968fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
9693425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9703425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
971674ae819SStefano Zampini   /* Apply partition of unity */
9723425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
973674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
9748eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
9753425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
9763425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
9773425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
9783425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
9793425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
9803425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9813425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
982674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
9833425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9843425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9853425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
9863425bc38SStefano Zampini   }
9873425bc38SStefano Zampini   /* BDDC rhs */
9883425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
9898eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
9903425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
9913425bc38SStefano Zampini   }
9923425bc38SStefano Zampini   /* apply BDDC */
9933425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
9943425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
9953425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
9963425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
9973425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9983425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9993425bc38SStefano Zampini   /* restore original rhs */
10003425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
10013425bc38SStefano Zampini   PetscFunctionReturn(0);
10023425bc38SStefano Zampini }
10031e6b0712SBarry Smith 
10043425bc38SStefano Zampini #undef __FUNCT__
10053425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
10063425bc38SStefano Zampini /*@
10073425bc38SStefano Zampini  PCBDDCMatFETIDPGetRHS - Get rhs for FETIDP linear system.
10083425bc38SStefano Zampini 
10093425bc38SStefano Zampini    Collective
10103425bc38SStefano Zampini 
10113425bc38SStefano Zampini    Input Parameters:
10123425bc38SStefano Zampini +  fetidp_mat   - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators
10133425bc38SStefano Zampini +  standard_rhs - the rhs of your linear system
10143425bc38SStefano Zampini 
10153425bc38SStefano Zampini    Output Parameters:
10163425bc38SStefano Zampini +  fetidp_flux_rhs   - the rhs of the FETIDP linear system
10173425bc38SStefano Zampini 
10183425bc38SStefano Zampini    Level: developer
10193425bc38SStefano Zampini 
10203425bc38SStefano Zampini    Notes:
10213425bc38SStefano Zampini 
10223425bc38SStefano Zampini .seealso: PCBDDC
10233425bc38SStefano Zampini @*/
10243425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
10253425bc38SStefano Zampini {
1026674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10273425bc38SStefano Zampini   PetscErrorCode ierr;
10283425bc38SStefano Zampini 
10293425bc38SStefano Zampini   PetscFunctionBegin;
10303425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10313425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
10323425bc38SStefano Zampini   PetscFunctionReturn(0);
10333425bc38SStefano Zampini }
10343425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
10351e6b0712SBarry Smith 
10363425bc38SStefano Zampini #undef __FUNCT__
10373425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
10383425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
10393425bc38SStefano Zampini {
1040674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
10413425bc38SStefano Zampini   PC_IS*         pcis;
10423425bc38SStefano Zampini   PC_BDDC*       pcbddc;
10433425bc38SStefano Zampini   PetscErrorCode ierr;
10443425bc38SStefano Zampini 
10453425bc38SStefano Zampini   PetscFunctionBegin;
10463425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
10473425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
10483425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
10493425bc38SStefano Zampini 
10503425bc38SStefano Zampini   /* apply B_delta^T */
10513425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10523425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10533425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
10543425bc38SStefano Zampini   /* compute rhs for BDDC application */
10553425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
10568eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
10573425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10583425bc38SStefano Zampini   }
10593425bc38SStefano Zampini   /* apply BDDC */
10603425bc38SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr);
10613425bc38SStefano Zampini   /* put values into standard global vector */
10623425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10633425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10648eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
10653425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
10663425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
10673425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
10683425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
10693425bc38SStefano Zampini   }
10703425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10713425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10723425bc38SStefano Zampini   /* final change of basis if needed
10733425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
10743308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
10753425bc38SStefano Zampini   PetscFunctionReturn(0);
10763425bc38SStefano Zampini }
10771e6b0712SBarry Smith 
10783425bc38SStefano Zampini #undef __FUNCT__
10793425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
10803425bc38SStefano Zampini /*@
10813425bc38SStefano Zampini  PCBDDCMatFETIDPGetSolution - Get Solution for FETIDP linear system.
10823425bc38SStefano Zampini 
10833425bc38SStefano Zampini    Collective
10843425bc38SStefano Zampini 
10853425bc38SStefano Zampini    Input Parameters:
10863425bc38SStefano Zampini +  fetidp_mat        - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators
10873425bc38SStefano Zampini +  fetidp_flux_sol - the solution of the FETIDP linear system
10883425bc38SStefano Zampini 
10893425bc38SStefano Zampini    Output Parameters:
10903425bc38SStefano Zampini +  standard_sol      - the solution on the global domain
10913425bc38SStefano Zampini 
10923425bc38SStefano Zampini    Level: developer
10933425bc38SStefano Zampini 
10943425bc38SStefano Zampini    Notes:
10953425bc38SStefano Zampini 
10963425bc38SStefano Zampini .seealso: PCBDDC
10973425bc38SStefano Zampini @*/
10983425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
10993425bc38SStefano Zampini {
1100674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
11013425bc38SStefano Zampini   PetscErrorCode ierr;
11023425bc38SStefano Zampini 
11033425bc38SStefano Zampini   PetscFunctionBegin;
11043425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
11053425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
11063425bc38SStefano Zampini   PetscFunctionReturn(0);
11073425bc38SStefano Zampini }
11083425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
11091e6b0712SBarry Smith 
1110f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1111f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1112f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1113f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1114674ae819SStefano Zampini 
11153425bc38SStefano Zampini #undef __FUNCT__
11163425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
11173425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
11183425bc38SStefano Zampini {
1119674ae819SStefano Zampini 
1120674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
11213425bc38SStefano Zampini   Mat            newmat;
1122674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
11233425bc38SStefano Zampini   PC             newpc;
1124ce94432eSBarry Smith   MPI_Comm       comm;
11253425bc38SStefano Zampini   PetscErrorCode ierr;
11263425bc38SStefano Zampini 
11273425bc38SStefano Zampini   PetscFunctionBegin;
1128ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
11293425bc38SStefano Zampini   /* FETIDP linear matrix */
11303425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
11313425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
11323425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
11333425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
11343425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
11353425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
11363425bc38SStefano Zampini   /* FETIDP preconditioner */
11373425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
11383425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
11393425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
11403425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
11413425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
11423425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
11433425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
11443425bc38SStefano Zampini   ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr);
11453425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
11463425bc38SStefano Zampini   /* return pointers for objects created */
11473425bc38SStefano Zampini   *fetidp_mat=newmat;
11483425bc38SStefano Zampini   *fetidp_pc=newpc;
11493425bc38SStefano Zampini   PetscFunctionReturn(0);
11503425bc38SStefano Zampini }
11511e6b0712SBarry Smith 
11523425bc38SStefano Zampini #undef __FUNCT__
11533425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
11543425bc38SStefano Zampini /*@
11553425bc38SStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP.
11563425bc38SStefano Zampini 
11573425bc38SStefano Zampini    Collective
11583425bc38SStefano Zampini 
11593425bc38SStefano Zampini    Input Parameters:
11603425bc38SStefano Zampini +  pc - the BDDC preconditioning context (setup must be already called)
11613425bc38SStefano Zampini 
11623425bc38SStefano Zampini    Level: developer
11633425bc38SStefano Zampini 
11643425bc38SStefano Zampini    Notes:
11653425bc38SStefano Zampini 
11663425bc38SStefano Zampini .seealso: PCBDDC
11673425bc38SStefano Zampini @*/
11683425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
11693425bc38SStefano Zampini {
11703425bc38SStefano Zampini   PetscErrorCode ierr;
11713425bc38SStefano Zampini 
11723425bc38SStefano Zampini   PetscFunctionBegin;
11733425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11743425bc38SStefano Zampini   if (pc->setupcalled) {
1175516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1176f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
11773425bc38SStefano Zampini   PetscFunctionReturn(0);
11783425bc38SStefano Zampini }
11790c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1180da1bb401SStefano Zampini /*MC
1181da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
11820c7d97c5SJed Brown 
1183da1bb401SStefano Zampini    Options Database Keys:
1184da1bb401SStefano Zampini .    -pcbddc ??? -
1185da1bb401SStefano Zampini 
1186da1bb401SStefano Zampini    Level: intermediate
1187da1bb401SStefano Zampini 
1188da1bb401SStefano Zampini    Notes: The matrix used with this preconditioner must be of type MATIS
1189da1bb401SStefano Zampini 
1190da1bb401SStefano Zampini           Unlike more 'conventional' interface preconditioners, this iterates over ALL the
1191da1bb401SStefano Zampini           degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers
1192da1bb401SStefano Zampini           on the subdomains).
1193da1bb401SStefano Zampini 
1194da1bb401SStefano Zampini           Options for the coarse grid preconditioner can be set with -
1195da1bb401SStefano Zampini           Options for the Dirichlet subproblem can be set with -
1196da1bb401SStefano Zampini           Options for the Neumann subproblem can be set with -
1197da1bb401SStefano Zampini 
1198da1bb401SStefano Zampini    Contributed by Stefano Zampini
1199da1bb401SStefano Zampini 
1200da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1201da1bb401SStefano Zampini M*/
1202b2573a8aSBarry Smith 
1203da1bb401SStefano Zampini #undef __FUNCT__
1204da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
12058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1206da1bb401SStefano Zampini {
1207da1bb401SStefano Zampini   PetscErrorCode      ierr;
1208da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1209da1bb401SStefano Zampini 
1210da1bb401SStefano Zampini   PetscFunctionBegin;
1211da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1212da1bb401SStefano Zampini   ierr      = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr);
1213da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1214da1bb401SStefano Zampini 
1215da1bb401SStefano Zampini   /* create PCIS data structure */
1216da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1217da1bb401SStefano Zampini 
121847d04d0dSStefano Zampini   /* BDDC customization */
121947d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
122047d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
122147d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
122247d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
122347d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
122447d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
122547d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
122647d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
122747d04d0dSStefano Zampini 
1228674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
12290bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
12303972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1231534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1232534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1233534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1234da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1235da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1236da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1237da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1238da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
123915aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
124015aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1241da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1242da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1243da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1244da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1245da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1246da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1247da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1248da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1249da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1250da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
1251da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
1252da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
1253b76ba322SStefano Zampini   pcbddc->use_exact_dirichlet        = PETSC_TRUE;
125447d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
125547d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
12564fad6a16SStefano Zampini   pcbddc->current_level              = 0;
12572b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1258da1bb401SStefano Zampini 
1259674ae819SStefano Zampini   /* create local graph structure */
1260674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1261674ae819SStefano Zampini 
1262674ae819SStefano Zampini   /* scaling */
1263674ae819SStefano Zampini   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
1264674ae819SStefano Zampini   pcbddc->work_scaling               = 0;
1265da1bb401SStefano Zampini 
1266da1bb401SStefano Zampini   /* function pointers */
1267da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
1268da1bb401SStefano Zampini   pc->ops->applytranspose      = 0;
1269da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1270da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1271da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1272da1bb401SStefano Zampini   pc->ops->view                = 0;
1273da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1274da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1275da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1276534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1277534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1278da1bb401SStefano Zampini 
1279da1bb401SStefano Zampini   /* composing function */
1280674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1281bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
12822b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1283b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
12842b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1285bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1286bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1287bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1288bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
1289bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
1290bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
1291bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1292bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1293bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1294bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1295da1bb401SStefano Zampini   PetscFunctionReturn(0);
1296da1bb401SStefano Zampini }
12973425bc38SStefano Zampini 
1298