153cdbc3dSStefano Zampini /* TODOLIST 2*85c4d303SStefano Zampini Dont use qr when number of primal dof per cc is 1 3*85c4d303SStefano Zampini Provide PCApplyTranpose 4*85c4d303SStefano Zampini why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver? 5*85c4d303SStefano Zampini Is it possible working with PCBDDCGraph on boundary indices only? 6da1bb401SStefano Zampini DofSplitting and DM attached to pc? 7da1bb401SStefano Zampini Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet) 8b8ffe317SStefano Zampini BDDC with MG framework? 9a0ba757dSStefano Zampini provide other ops? Ask to developers 10a0ba757dSStefano Zampini man pages 1153cdbc3dSStefano Zampini */ 120c7d97c5SJed Brown 1353cdbc3dSStefano Zampini /* ---------------------------------------------------------------------------------------------------------------------------------------------- 140c7d97c5SJed Brown Implementation of BDDC preconditioner based on: 150c7d97c5SJed Brown C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 1653cdbc3dSStefano Zampini ---------------------------------------------------------------------------------------------------------------------------------------------- */ 1753cdbc3dSStefano Zampini 18674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 19674ae819SStefano Zampini #include "bddcprivate.h" 203b03a366Sstefano_zampini #include <petscblaslapack.h> 21674ae819SStefano Zampini 220c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 230c7d97c5SJed Brown #undef __FUNCT__ 240c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC" 250c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc) 260c7d97c5SJed Brown { 270c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 280c7d97c5SJed Brown PetscErrorCode ierr; 290c7d97c5SJed Brown 300c7d97c5SJed Brown PetscFunctionBegin; 310c7d97c5SJed Brown ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 328eeda7d8SStefano Zampini /* Verbose debugging */ 338eeda7d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 348eeda7d8SStefano Zampini /* Primal space cumstomization */ 358eeda7d8SStefano 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); 368eeda7d8SStefano 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); 378eeda7d8SStefano 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); 388eeda7d8SStefano Zampini /* Change of basis */ 398eeda7d8SStefano 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); 408eeda7d8SStefano 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); 41674ae819SStefano Zampini if (!pcbddc->use_change_of_basis) { 42674ae819SStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 43674ae819SStefano Zampini } 448eeda7d8SStefano Zampini /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 458eeda7d8SStefano 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); 460298fd71SBarry Smith ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 472b510759SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 48674ae819SStefano 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); 490c7d97c5SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 500c7d97c5SJed Brown PetscFunctionReturn(0); 510c7d97c5SJed Brown } 520c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 53674ae819SStefano Zampini #undef __FUNCT__ 54674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 55674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 56674ae819SStefano Zampini { 57674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 58674ae819SStefano Zampini PetscErrorCode ierr; 591e6b0712SBarry Smith 60674ae819SStefano Zampini PetscFunctionBegin; 61674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 62674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 63674ae819SStefano Zampini pcbddc->user_primal_vertices = PrimalVertices; 64674ae819SStefano Zampini PetscFunctionReturn(0); 65674ae819SStefano Zampini } 66674ae819SStefano Zampini #undef __FUNCT__ 67674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 68674ae819SStefano Zampini /*@ 69674ae819SStefano Zampini PCBDDCSetPrimalVerticesLocalIS - Set user defined primal vertices in PCBDDC. 70674ae819SStefano Zampini 71674ae819SStefano Zampini Not collective 72674ae819SStefano Zampini 73674ae819SStefano Zampini Input Parameters: 74674ae819SStefano Zampini + pc - the preconditioning context 75674ae819SStefano Zampini - PrimalVertices - index sets of primal vertices in local numbering 76674ae819SStefano Zampini 77674ae819SStefano Zampini Level: intermediate 78674ae819SStefano Zampini 79674ae819SStefano Zampini Notes: 80674ae819SStefano Zampini 81674ae819SStefano Zampini .seealso: PCBDDC 82674ae819SStefano Zampini @*/ 83674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 84674ae819SStefano Zampini { 85674ae819SStefano Zampini PetscErrorCode ierr; 86674ae819SStefano Zampini 87674ae819SStefano Zampini PetscFunctionBegin; 88674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 89674ae819SStefano Zampini PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 90674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 91674ae819SStefano Zampini PetscFunctionReturn(0); 92674ae819SStefano Zampini } 93674ae819SStefano Zampini /* -------------------------------------------------------------------------- */ 940c7d97c5SJed Brown #undef __FUNCT__ 954fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 964fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 974fad6a16SStefano Zampini { 984fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 994fad6a16SStefano Zampini 1004fad6a16SStefano Zampini PetscFunctionBegin; 1014fad6a16SStefano Zampini pcbddc->coarsening_ratio = k; 1024fad6a16SStefano Zampini PetscFunctionReturn(0); 1034fad6a16SStefano Zampini } 1041e6b0712SBarry Smith 1054fad6a16SStefano Zampini #undef __FUNCT__ 1064fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio" 1074fad6a16SStefano Zampini /*@ 1084fad6a16SStefano Zampini PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel coarsening 1094fad6a16SStefano Zampini 1104fad6a16SStefano Zampini Logically collective on PC 1114fad6a16SStefano Zampini 1124fad6a16SStefano Zampini Input Parameters: 1134fad6a16SStefano Zampini + pc - the preconditioning context 1144fad6a16SStefano Zampini - k - coarsening ratio 1154fad6a16SStefano Zampini 1164fad6a16SStefano Zampini Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level. 1174fad6a16SStefano Zampini 1184fad6a16SStefano Zampini Level: intermediate 1194fad6a16SStefano Zampini 1204fad6a16SStefano Zampini Notes: 1214fad6a16SStefano Zampini 1224fad6a16SStefano Zampini .seealso: PCBDDC 1234fad6a16SStefano Zampini @*/ 1244fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 1254fad6a16SStefano Zampini { 1264fad6a16SStefano Zampini PetscErrorCode ierr; 1274fad6a16SStefano Zampini 1284fad6a16SStefano Zampini PetscFunctionBegin; 1294fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1302b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,k,2); 1314fad6a16SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 1324fad6a16SStefano Zampini PetscFunctionReturn(0); 1334fad6a16SStefano Zampini } 1342b510759SStefano Zampini 135b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 1362b510759SStefano Zampini #undef __FUNCT__ 137b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 138b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 139b8ffe317SStefano Zampini { 140b8ffe317SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 141b8ffe317SStefano Zampini 142b8ffe317SStefano Zampini PetscFunctionBegin; 143*85c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = flg; 144b8ffe317SStefano Zampini PetscFunctionReturn(0); 145b8ffe317SStefano Zampini } 146b8ffe317SStefano Zampini 147b8ffe317SStefano Zampini #undef __FUNCT__ 148b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 149b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 1502b510759SStefano Zampini { 1512b510759SStefano Zampini PetscErrorCode ierr; 1522b510759SStefano Zampini 1532b510759SStefano Zampini PetscFunctionBegin; 1542b510759SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 155b8ffe317SStefano Zampini PetscValidLogicalCollectiveBool(pc,flg,2); 156b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1572b510759SStefano Zampini PetscFunctionReturn(0); 1582b510759SStefano Zampini } 1591e6b0712SBarry Smith 1604fad6a16SStefano Zampini #undef __FUNCT__ 1612b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC" 1622b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 1634fad6a16SStefano Zampini { 1644fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1654fad6a16SStefano Zampini 1664fad6a16SStefano Zampini PetscFunctionBegin; 1672b510759SStefano Zampini pcbddc->current_level = level; 1684fad6a16SStefano Zampini PetscFunctionReturn(0); 1694fad6a16SStefano Zampini } 1701e6b0712SBarry Smith 1714fad6a16SStefano Zampini #undef __FUNCT__ 172b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel" 173b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 174b8ffe317SStefano Zampini { 175b8ffe317SStefano Zampini PetscErrorCode ierr; 176b8ffe317SStefano Zampini 177b8ffe317SStefano Zampini PetscFunctionBegin; 178b8ffe317SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 179b8ffe317SStefano Zampini PetscValidLogicalCollectiveInt(pc,level,2); 180b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 181b8ffe317SStefano Zampini PetscFunctionReturn(0); 182b8ffe317SStefano Zampini } 183b8ffe317SStefano Zampini 184b8ffe317SStefano Zampini #undef __FUNCT__ 1852b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC" 1862b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 1872b510759SStefano Zampini { 1882b510759SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1892b510759SStefano Zampini 1902b510759SStefano Zampini PetscFunctionBegin; 1912b510759SStefano Zampini pcbddc->max_levels = levels; 1922b510759SStefano Zampini PetscFunctionReturn(0); 1932b510759SStefano Zampini } 1942b510759SStefano Zampini 1952b510759SStefano Zampini #undef __FUNCT__ 1962b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels" 1974fad6a16SStefano Zampini /*@ 1982b510759SStefano Zampini PCBDDCSetLevels - Sets the maximum number of levels within the multilevel approach. 1994fad6a16SStefano Zampini 2004fad6a16SStefano Zampini Logically collective on PC 2014fad6a16SStefano Zampini 2024fad6a16SStefano Zampini Input Parameters: 2034fad6a16SStefano Zampini + pc - the preconditioning context 2042b510759SStefano Zampini - levels - the maximum number of levels 2054fad6a16SStefano Zampini 2062b510759SStefano Zampini Default value is 0, i.e. coarse problem will be solved exactly 2074fad6a16SStefano Zampini 2084fad6a16SStefano Zampini Level: intermediate 2094fad6a16SStefano Zampini 2104fad6a16SStefano Zampini Notes: 2114fad6a16SStefano Zampini 2124fad6a16SStefano Zampini .seealso: PCBDDC 2134fad6a16SStefano Zampini @*/ 2142b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 2154fad6a16SStefano Zampini { 2164fad6a16SStefano Zampini PetscErrorCode ierr; 2174fad6a16SStefano Zampini 2184fad6a16SStefano Zampini PetscFunctionBegin; 2194fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2202b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,levels,2); 2212b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 2224fad6a16SStefano Zampini PetscFunctionReturn(0); 2234fad6a16SStefano Zampini } 2244fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */ 2251e6b0712SBarry Smith 2264fad6a16SStefano Zampini #undef __FUNCT__ 2270bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 2280bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 2290bdf917eSStefano Zampini { 2300bdf917eSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2310bdf917eSStefano Zampini PetscErrorCode ierr; 2320bdf917eSStefano Zampini 2330bdf917eSStefano Zampini PetscFunctionBegin; 2340bdf917eSStefano Zampini ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 2350bdf917eSStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 2360bdf917eSStefano Zampini pcbddc->NullSpace = NullSpace; 2370bdf917eSStefano Zampini PetscFunctionReturn(0); 2380bdf917eSStefano Zampini } 2391e6b0712SBarry Smith 2400bdf917eSStefano Zampini #undef __FUNCT__ 2410bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace" 2420bdf917eSStefano Zampini /*@ 2430bdf917eSStefano Zampini PCBDDCSetNullSpace - Set NullSpace of global operator of BDDC preconditioned mat. 2440bdf917eSStefano Zampini 2450bdf917eSStefano Zampini Logically collective on PC and MatNullSpace 2460bdf917eSStefano Zampini 2470bdf917eSStefano Zampini Input Parameters: 2480bdf917eSStefano Zampini + pc - the preconditioning context 2490bdf917eSStefano Zampini - NullSpace - Null space of the linear operator to be preconditioned. 2500bdf917eSStefano Zampini 2510bdf917eSStefano Zampini Level: intermediate 2520bdf917eSStefano Zampini 2530bdf917eSStefano Zampini Notes: 2540bdf917eSStefano Zampini 2550bdf917eSStefano Zampini .seealso: PCBDDC 2560bdf917eSStefano Zampini @*/ 2570bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 2580bdf917eSStefano Zampini { 2590bdf917eSStefano Zampini PetscErrorCode ierr; 2600bdf917eSStefano Zampini 2610bdf917eSStefano Zampini PetscFunctionBegin; 2620bdf917eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 263674ae819SStefano Zampini PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 2642b510759SStefano Zampini PetscCheckSameComm(pc,1,NullSpace,2); 2650bdf917eSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 2660bdf917eSStefano Zampini PetscFunctionReturn(0); 2670bdf917eSStefano Zampini } 2680bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */ 2691e6b0712SBarry Smith 2700bdf917eSStefano Zampini #undef __FUNCT__ 2713b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 2723b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 2733b03a366Sstefano_zampini { 2743b03a366Sstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2753b03a366Sstefano_zampini PetscErrorCode ierr; 2763b03a366Sstefano_zampini 2773b03a366Sstefano_zampini PetscFunctionBegin; 2783b03a366Sstefano_zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 27936e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 28036e030ebSStefano Zampini pcbddc->DirichletBoundaries=DirichletBoundaries; 2813b03a366Sstefano_zampini PetscFunctionReturn(0); 2823b03a366Sstefano_zampini } 2831e6b0712SBarry Smith 2843b03a366Sstefano_zampini #undef __FUNCT__ 2853b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 2863b03a366Sstefano_zampini /*@ 287da1bb401SStefano Zampini PCBDDCSetDirichletBoundaries - Set index set defining subdomain part (in local ordering) 288da1bb401SStefano Zampini of Dirichlet boundaries for the global problem. 2893b03a366Sstefano_zampini 2903b03a366Sstefano_zampini Not collective 2913b03a366Sstefano_zampini 2923b03a366Sstefano_zampini Input Parameters: 2933b03a366Sstefano_zampini + pc - the preconditioning context 2940298fd71SBarry Smith - DirichletBoundaries - sequential index set defining the subdomain part of Dirichlet boundaries (can be NULL) 2953b03a366Sstefano_zampini 2963b03a366Sstefano_zampini Level: intermediate 2973b03a366Sstefano_zampini 2983b03a366Sstefano_zampini Notes: 2993b03a366Sstefano_zampini 3003b03a366Sstefano_zampini .seealso: PCBDDC 3013b03a366Sstefano_zampini @*/ 3023b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 3033b03a366Sstefano_zampini { 3043b03a366Sstefano_zampini PetscErrorCode ierr; 3053b03a366Sstefano_zampini 3063b03a366Sstefano_zampini PetscFunctionBegin; 3073b03a366Sstefano_zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 308674ae819SStefano Zampini PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 3093b03a366Sstefano_zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 3103b03a366Sstefano_zampini PetscFunctionReturn(0); 3113b03a366Sstefano_zampini } 3123b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */ 3131e6b0712SBarry Smith 3143b03a366Sstefano_zampini #undef __FUNCT__ 3150c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 31653cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 3170c7d97c5SJed Brown { 3180c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 31953cdbc3dSStefano Zampini PetscErrorCode ierr; 3200c7d97c5SJed Brown 3210c7d97c5SJed Brown PetscFunctionBegin; 32253cdbc3dSStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 32336e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 32436e030ebSStefano Zampini pcbddc->NeumannBoundaries=NeumannBoundaries; 3250c7d97c5SJed Brown PetscFunctionReturn(0); 3260c7d97c5SJed Brown } 3271e6b0712SBarry Smith 3280c7d97c5SJed Brown #undef __FUNCT__ 3290c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 33057527edcSJed Brown /*@ 331da1bb401SStefano Zampini PCBDDCSetNeumannBoundaries - Set index set defining subdomain part (in local ordering) 332da1bb401SStefano Zampini of Neumann boundaries for the global problem. 33357527edcSJed Brown 3349c0446d6SStefano Zampini Not collective 33557527edcSJed Brown 33657527edcSJed Brown Input Parameters: 33757527edcSJed Brown + pc - the preconditioning context 3380298fd71SBarry Smith - NeumannBoundaries - sequential index set defining the subdomain part of Neumann boundaries (can be NULL) 33957527edcSJed Brown 34057527edcSJed Brown Level: intermediate 34157527edcSJed Brown 34257527edcSJed Brown Notes: 34357527edcSJed Brown 34457527edcSJed Brown .seealso: PCBDDC 34557527edcSJed Brown @*/ 34653cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 3470c7d97c5SJed Brown { 3480c7d97c5SJed Brown PetscErrorCode ierr; 3490c7d97c5SJed Brown 3500c7d97c5SJed Brown PetscFunctionBegin; 3510c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 352674ae819SStefano Zampini PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 35353cdbc3dSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 35453cdbc3dSStefano Zampini PetscFunctionReturn(0); 35553cdbc3dSStefano Zampini } 35653cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 3571e6b0712SBarry Smith 35853cdbc3dSStefano Zampini #undef __FUNCT__ 359da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 360da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 361da1bb401SStefano Zampini { 362da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 363da1bb401SStefano Zampini 364da1bb401SStefano Zampini PetscFunctionBegin; 365da1bb401SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundaries; 366da1bb401SStefano Zampini PetscFunctionReturn(0); 367da1bb401SStefano Zampini } 3681e6b0712SBarry Smith 369da1bb401SStefano Zampini #undef __FUNCT__ 370da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 371da1bb401SStefano Zampini /*@ 372da1bb401SStefano Zampini PCBDDCGetDirichletBoundaries - Get index set defining subdomain part (in local ordering) 373da1bb401SStefano Zampini of Dirichlet boundaries for the global problem. 374da1bb401SStefano Zampini 375da1bb401SStefano Zampini Not collective 376da1bb401SStefano Zampini 377da1bb401SStefano Zampini Input Parameters: 378da1bb401SStefano Zampini + pc - the preconditioning context 379da1bb401SStefano Zampini 380da1bb401SStefano Zampini Output Parameters: 381da1bb401SStefano Zampini + DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 382da1bb401SStefano Zampini 383da1bb401SStefano Zampini Level: intermediate 384da1bb401SStefano Zampini 385da1bb401SStefano Zampini Notes: 386da1bb401SStefano Zampini 387da1bb401SStefano Zampini .seealso: PCBDDC 388da1bb401SStefano Zampini @*/ 389da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 390da1bb401SStefano Zampini { 391da1bb401SStefano Zampini PetscErrorCode ierr; 392da1bb401SStefano Zampini 393da1bb401SStefano Zampini PetscFunctionBegin; 394da1bb401SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 395da1bb401SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 396da1bb401SStefano Zampini PetscFunctionReturn(0); 397da1bb401SStefano Zampini } 398da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 3991e6b0712SBarry Smith 400da1bb401SStefano Zampini #undef __FUNCT__ 40153cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 40253cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 40353cdbc3dSStefano Zampini { 40453cdbc3dSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 40553cdbc3dSStefano Zampini 40653cdbc3dSStefano Zampini PetscFunctionBegin; 40753cdbc3dSStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundaries; 40853cdbc3dSStefano Zampini PetscFunctionReturn(0); 40953cdbc3dSStefano Zampini } 4101e6b0712SBarry Smith 41153cdbc3dSStefano Zampini #undef __FUNCT__ 41253cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 41353cdbc3dSStefano Zampini /*@ 414da1bb401SStefano Zampini PCBDDCGetNeumannBoundaries - Get index set defining subdomain part (in local ordering) 415da1bb401SStefano Zampini of Neumann boundaries for the global problem. 41653cdbc3dSStefano Zampini 4179c0446d6SStefano Zampini Not collective 41853cdbc3dSStefano Zampini 41953cdbc3dSStefano Zampini Input Parameters: 42053cdbc3dSStefano Zampini + pc - the preconditioning context 42153cdbc3dSStefano Zampini 42253cdbc3dSStefano Zampini Output Parameters: 42353cdbc3dSStefano Zampini + NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 42453cdbc3dSStefano Zampini 42553cdbc3dSStefano Zampini Level: intermediate 42653cdbc3dSStefano Zampini 42753cdbc3dSStefano Zampini Notes: 42853cdbc3dSStefano Zampini 42953cdbc3dSStefano Zampini .seealso: PCBDDC 43053cdbc3dSStefano Zampini @*/ 43153cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 43253cdbc3dSStefano Zampini { 43353cdbc3dSStefano Zampini PetscErrorCode ierr; 43453cdbc3dSStefano Zampini 43553cdbc3dSStefano Zampini PetscFunctionBegin; 43653cdbc3dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 43753cdbc3dSStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 4380c7d97c5SJed Brown PetscFunctionReturn(0); 4390c7d97c5SJed Brown } 44036e030ebSStefano Zampini /* -------------------------------------------------------------------------- */ 4411e6b0712SBarry Smith 44236e030ebSStefano Zampini #undef __FUNCT__ 443da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 4441a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 44536e030ebSStefano Zampini { 44636e030ebSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 447da1bb401SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 448da1bb401SStefano Zampini PetscErrorCode ierr; 44936e030ebSStefano Zampini 45036e030ebSStefano Zampini PetscFunctionBegin; 451674ae819SStefano Zampini /* free old CSR */ 452674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 453fb223d50SStefano Zampini /* TODO: PCBDDCGraphSetAdjacency */ 454674ae819SStefano Zampini /* get CSR into graph structure */ 455da1bb401SStefano Zampini if (copymode == PETSC_COPY_VALUES) { 456674ae819SStefano Zampini ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr); 457674ae819SStefano Zampini ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr); 458674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 459674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 460da1bb401SStefano Zampini } else if (copymode == PETSC_OWN_POINTER) { 4611a83f524SJed Brown mat_graph->xadj = (PetscInt*)xadj; 4621a83f524SJed Brown mat_graph->adjncy = (PetscInt*)adjncy; 463674ae819SStefano Zampini } 464575ad6abSStefano Zampini mat_graph->nvtxs_csr = nvtxs; 46536e030ebSStefano Zampini PetscFunctionReturn(0); 46636e030ebSStefano Zampini } 4671e6b0712SBarry Smith 46836e030ebSStefano Zampini #undef __FUNCT__ 469da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 47036e030ebSStefano Zampini /*@ 471da1bb401SStefano Zampini PCBDDCSetLocalAdjacencyGraph - Set CSR graph of local matrix for use of PCBDDC. 47236e030ebSStefano Zampini 47336e030ebSStefano Zampini Not collective 47436e030ebSStefano Zampini 47536e030ebSStefano Zampini Input Parameters: 47636e030ebSStefano Zampini + pc - the preconditioning context 477da1bb401SStefano Zampini - nvtxs - number of local vertices of the graph 478da1bb401SStefano Zampini - xadj, adjncy - the CSR graph 479da1bb401SStefano Zampini - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. In the former case the user must free the array passed in; 480da1bb401SStefano Zampini in the latter case, memory must be obtained with PetscMalloc. 48136e030ebSStefano Zampini 48236e030ebSStefano Zampini Level: intermediate 48336e030ebSStefano Zampini 48436e030ebSStefano Zampini Notes: 48536e030ebSStefano Zampini 48636e030ebSStefano Zampini .seealso: PCBDDC 48736e030ebSStefano Zampini @*/ 4881a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 48936e030ebSStefano Zampini { 490575ad6abSStefano Zampini void (*f)(void) = 0; 49136e030ebSStefano Zampini PetscErrorCode ierr; 49236e030ebSStefano Zampini 49336e030ebSStefano Zampini PetscFunctionBegin; 49436e030ebSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 495674ae819SStefano Zampini PetscValidIntPointer(xadj,3); 496674ae819SStefano Zampini PetscValidIntPointer(xadj,4); 497674ae819SStefano Zampini if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 498674ae819SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 499da1bb401SStefano Zampini } 50036e030ebSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 501575ad6abSStefano Zampini /* free arrays if PCBDDC is not the PC type */ 502575ad6abSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 503575ad6abSStefano Zampini if (!f && copymode == PETSC_OWN_POINTER) { 504575ad6abSStefano Zampini ierr = PetscFree(xadj);CHKERRQ(ierr); 505575ad6abSStefano Zampini ierr = PetscFree(adjncy);CHKERRQ(ierr); 50636e030ebSStefano Zampini } 50736e030ebSStefano Zampini PetscFunctionReturn(0); 50836e030ebSStefano Zampini } 5099c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */ 5101e6b0712SBarry Smith 5119c0446d6SStefano Zampini #undef __FUNCT__ 5129c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 5139c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 5149c0446d6SStefano Zampini { 5159c0446d6SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 5169c0446d6SStefano Zampini PetscInt i; 5179c0446d6SStefano Zampini PetscErrorCode ierr; 5189c0446d6SStefano Zampini 5199c0446d6SStefano Zampini PetscFunctionBegin; 520da1bb401SStefano Zampini /* Destroy ISes if they were already set */ 5219c0446d6SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 5229c0446d6SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 5239c0446d6SStefano Zampini } 524d11ae9bbSstefano_zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 525da1bb401SStefano Zampini /* allocate space then set */ 5269c0446d6SStefano Zampini ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr); 5279c0446d6SStefano Zampini for (i=0;i<n_is;i++) { 528da1bb401SStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 529da1bb401SStefano Zampini pcbddc->ISForDofs[i]=ISForDofs[i]; 5309c0446d6SStefano Zampini } 5319c0446d6SStefano Zampini pcbddc->n_ISForDofs=n_is; 5329c0446d6SStefano Zampini PetscFunctionReturn(0); 5339c0446d6SStefano Zampini } 5341e6b0712SBarry Smith 5359c0446d6SStefano Zampini #undef __FUNCT__ 5369c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting" 5379c0446d6SStefano Zampini /*@ 538da1bb401SStefano Zampini PCBDDCSetDofsSplitting - Set index sets defining fields of local mat. 5399c0446d6SStefano Zampini 5409c0446d6SStefano Zampini Not collective 5419c0446d6SStefano Zampini 5429c0446d6SStefano Zampini Input Parameters: 5439c0446d6SStefano Zampini + pc - the preconditioning context 544da1bb401SStefano Zampini - n - number of index sets defining the fields 545da1bb401SStefano Zampini - IS[] - array of IS describing the fields 5469c0446d6SStefano Zampini 5479c0446d6SStefano Zampini Level: intermediate 5489c0446d6SStefano Zampini 5499c0446d6SStefano Zampini Notes: 5509c0446d6SStefano Zampini 5519c0446d6SStefano Zampini .seealso: PCBDDC 5529c0446d6SStefano Zampini @*/ 5539c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 5549c0446d6SStefano Zampini { 5552b510759SStefano Zampini PetscInt i; 5569c0446d6SStefano Zampini PetscErrorCode ierr; 5579c0446d6SStefano Zampini 5589c0446d6SStefano Zampini PetscFunctionBegin; 5599c0446d6SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5602b510759SStefano Zampini for (i=0;i<n_is;i++) { 5612b510759SStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2); 5622b510759SStefano Zampini } 5639c0446d6SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 5649c0446d6SStefano Zampini PetscFunctionReturn(0); 5659c0446d6SStefano Zampini } 566da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 567534831adSStefano Zampini #undef __FUNCT__ 568534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC" 569534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 570534831adSStefano Zampini /* 571534831adSStefano Zampini PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 572534831adSStefano Zampini guess if a transformation of basis approach has been selected. 5739c0446d6SStefano Zampini 574534831adSStefano Zampini Input Parameter: 575534831adSStefano Zampini + pc - the preconditioner contex 576534831adSStefano Zampini 577534831adSStefano Zampini Application Interface Routine: PCPreSolve() 578534831adSStefano Zampini 579534831adSStefano Zampini Notes: 580534831adSStefano Zampini The interface routine PCPreSolve() is not usually called directly by 581534831adSStefano Zampini the user, but instead is called by KSPSolve(). 582534831adSStefano Zampini */ 583534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 584534831adSStefano Zampini { 585534831adSStefano Zampini PetscErrorCode ierr; 586534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 587534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 588534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 589534831adSStefano Zampini Mat temp_mat; 5903972b0daSStefano Zampini IS dirIS; 5913972b0daSStefano Zampini PetscInt dirsize,i,*is_indices; 5923972b0daSStefano Zampini PetscScalar *array_x,*array_diagonal; 5933972b0daSStefano Zampini Vec used_vec; 59492e3dcfbSStefano Zampini PetscBool guess_nonzero,flg,bddc_has_dirichlet_boundaries; 595534831adSStefano Zampini 596534831adSStefano Zampini PetscFunctionBegin; 597*85c4d303SStefano Zampini /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 598*85c4d303SStefano Zampini if (ksp) { 599*85c4d303SStefano Zampini PetscBool iscg; 600*85c4d303SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 601*85c4d303SStefano Zampini if (!iscg) { 602*85c4d303SStefano Zampini ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 603*85c4d303SStefano Zampini } 604*85c4d303SStefano Zampini } 605*85c4d303SStefano Zampini /* Creates parallel work vectors used in presolve */ 60662a6ff1dSStefano Zampini if (!pcbddc->original_rhs) { 60762a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 60862a6ff1dSStefano Zampini } 60962a6ff1dSStefano Zampini if (!pcbddc->temp_solution) { 61062a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 61162a6ff1dSStefano Zampini } 6123972b0daSStefano Zampini if (x) { 6133972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 6143972b0daSStefano Zampini used_vec = x; 6153972b0daSStefano Zampini } else { 6163972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 6173972b0daSStefano Zampini used_vec = pcbddc->temp_solution; 6183972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 6193972b0daSStefano Zampini } 6203972b0daSStefano Zampini /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 6213972b0daSStefano Zampini if (ksp) { 6223972b0daSStefano Zampini ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 6233972b0daSStefano Zampini if (!guess_nonzero) { 6243972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 6253972b0daSStefano Zampini } 6263972b0daSStefano Zampini } 6273308cffdSStefano Zampini 62892e3dcfbSStefano Zampini /* TODO: remove when Dirichlet boundaries will be shared */ 62992e3dcfbSStefano Zampini ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr); 63092e3dcfbSStefano Zampini flg = PETSC_FALSE; 63192e3dcfbSStefano Zampini if (dirIS) flg = PETSC_TRUE; 63292e3dcfbSStefano Zampini ierr = MPI_Allreduce(&flg,&bddc_has_dirichlet_boundaries,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 63392e3dcfbSStefano Zampini 6343972b0daSStefano Zampini /* store the original rhs */ 6353972b0daSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 6363972b0daSStefano Zampini 6373972b0daSStefano Zampini /* Take into account zeroed rows -> change rhs and store solution removed */ 63892e3dcfbSStefano Zampini if (rhs && bddc_has_dirichlet_boundaries) { 6393972b0daSStefano Zampini ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 6403972b0daSStefano Zampini ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 6413972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6423972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6433972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6443972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6453972b0daSStefano Zampini if (dirIS) { 6463972b0daSStefano Zampini ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr); 6473972b0daSStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 6483972b0daSStefano Zampini ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 6493972b0daSStefano Zampini ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 6502fa5cd67SKarl Rupp for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 6513972b0daSStefano Zampini ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 6523972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 6533972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 6543972b0daSStefano Zampini } 6553972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6563972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 657b76ba322SStefano Zampini 6583972b0daSStefano Zampini /* remove the computed solution from the rhs */ 6593972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 6603972b0daSStefano Zampini ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 6613972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 6623308cffdSStefano Zampini } 663b76ba322SStefano Zampini 664b76ba322SStefano Zampini /* store partially computed solution and set initial guess */ 6653972b0daSStefano Zampini if (x) { 6663972b0daSStefano Zampini ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 6673972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 668*85c4d303SStefano Zampini if (pcbddc->use_exact_dirichlet_trick) { 669b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 670b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 671b76ba322SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 672b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 673b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 674b76ba322SStefano Zampini if (ksp) { 675b76ba322SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 676b76ba322SStefano Zampini } 677b76ba322SStefano Zampini } 6783972b0daSStefano Zampini } 679b76ba322SStefano Zampini 68092e3dcfbSStefano Zampini /* prepare MatMult and rhs for solver */ 681674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 682b76ba322SStefano Zampini /* swap pointers for local matrices */ 683b76ba322SStefano Zampini temp_mat = matis->A; 684b76ba322SStefano Zampini matis->A = pcbddc->local_mat; 685b76ba322SStefano Zampini pcbddc->local_mat = temp_mat; 68692e3dcfbSStefano Zampini if (rhs) { 687b76ba322SStefano Zampini /* Get local rhs and apply transformation of basis */ 688b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 689b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 690b76ba322SStefano Zampini /* from original basis to modified basis */ 691b76ba322SStefano Zampini ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 692b76ba322SStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 693b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 694b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 695674ae819SStefano Zampini } 69692e3dcfbSStefano Zampini } 69792e3dcfbSStefano Zampini 69892e3dcfbSStefano Zampini /* remove nullspace if present */ 6990bdf917eSStefano Zampini if (ksp && pcbddc->NullSpace) { 700d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 701d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 702b76ba322SStefano Zampini } 7030bdf917eSStefano Zampini ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 704534831adSStefano Zampini PetscFunctionReturn(0); 705534831adSStefano Zampini } 706534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 707534831adSStefano Zampini #undef __FUNCT__ 708534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC" 709534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 710534831adSStefano Zampini /* 711534831adSStefano Zampini PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 712534831adSStefano Zampini approach has been selected. Also, restores rhs to its original state. 713534831adSStefano Zampini 714534831adSStefano Zampini Input Parameter: 715534831adSStefano Zampini + pc - the preconditioner contex 716534831adSStefano Zampini 717534831adSStefano Zampini Application Interface Routine: PCPostSolve() 718534831adSStefano Zampini 719534831adSStefano Zampini Notes: 720534831adSStefano Zampini The interface routine PCPostSolve() is not usually called directly by 721534831adSStefano Zampini the user, but instead is called by KSPSolve(). 722534831adSStefano Zampini */ 723534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 724534831adSStefano Zampini { 725534831adSStefano Zampini PetscErrorCode ierr; 726534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 727534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 728534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 729534831adSStefano Zampini Mat temp_mat; 730534831adSStefano Zampini 731534831adSStefano Zampini PetscFunctionBegin; 732674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 733534831adSStefano Zampini /* swap pointers for local matrices */ 734534831adSStefano Zampini temp_mat = matis->A; 735534831adSStefano Zampini matis->A = pcbddc->local_mat; 736534831adSStefano Zampini pcbddc->local_mat = temp_mat; 7373425bc38SStefano Zampini } 7383308cffdSStefano Zampini if (pcbddc->use_change_of_basis && x) { 739534831adSStefano Zampini /* Get Local boundary and apply transformation of basis to solution vector */ 740534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 741534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 742534831adSStefano Zampini /* from modified basis to original basis */ 743534831adSStefano Zampini ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 744534831adSStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 745534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 746534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 747534831adSStefano Zampini } 7483972b0daSStefano Zampini /* add solution removed in presolve */ 7493425bc38SStefano Zampini if (x) { 7503425bc38SStefano Zampini ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 7513425bc38SStefano Zampini } 752fb223d50SStefano Zampini /* restore rhs to its original state */ 753fb223d50SStefano Zampini if (rhs) { 754fb223d50SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 755fb223d50SStefano Zampini } 756534831adSStefano Zampini PetscFunctionReturn(0); 757534831adSStefano Zampini } 758534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 75953cdbc3dSStefano Zampini #undef __FUNCT__ 76053cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC" 7610c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 7620c7d97c5SJed Brown /* 7630c7d97c5SJed Brown PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 7640c7d97c5SJed Brown by setting data structures and options. 7650c7d97c5SJed Brown 7660c7d97c5SJed Brown Input Parameter: 76753cdbc3dSStefano Zampini + pc - the preconditioner context 7680c7d97c5SJed Brown 7690c7d97c5SJed Brown Application Interface Routine: PCSetUp() 7700c7d97c5SJed Brown 7710c7d97c5SJed Brown Notes: 7720c7d97c5SJed Brown The interface routine PCSetUp() is not usually called directly by 7730c7d97c5SJed Brown the user, but instead is called by PCApply() if necessary. 7740c7d97c5SJed Brown */ 77553cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc) 7760c7d97c5SJed Brown { 7770c7d97c5SJed Brown PetscErrorCode ierr; 7780c7d97c5SJed Brown PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 779674ae819SStefano Zampini MatStructure flag; 780674ae819SStefano Zampini PetscBool computeis,computetopography,computesolvers; 7810c7d97c5SJed Brown 7820c7d97c5SJed Brown PetscFunctionBegin; 783674ae819SStefano Zampini /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */ 784fcd409f5SStefano Zampini /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */ 7853b03a366Sstefano_zampini /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 786fcd409f5SStefano Zampini Also, BDDC directly build the Dirichlet problem */ 7873b03a366Sstefano_zampini /* Get stdout for dbg */ 788674ae819SStefano Zampini if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) { 789ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr); 790e269702eSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 7912b510759SStefano Zampini if (pcbddc->current_level) { 7922b510759SStefano Zampini ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 7932b510759SStefano Zampini } 794e269702eSStefano Zampini } 795674ae819SStefano Zampini /* first attempt to split work */ 796674ae819SStefano Zampini if (pc->setupcalled) { 797674ae819SStefano Zampini computeis = PETSC_FALSE; 798674ae819SStefano Zampini ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr); 799674ae819SStefano Zampini if (flag == SAME_PRECONDITIONER) { 800674ae819SStefano Zampini computetopography = PETSC_FALSE; 801674ae819SStefano Zampini computesolvers = PETSC_FALSE; 802674ae819SStefano Zampini } else if (flag == SAME_NONZERO_PATTERN) { 803674ae819SStefano Zampini computetopography = PETSC_FALSE; 804674ae819SStefano Zampini computesolvers = PETSC_TRUE; 805674ae819SStefano Zampini } else { /* DIFFERENT_NONZERO_PATTERN */ 806674ae819SStefano Zampini computetopography = PETSC_TRUE; 807674ae819SStefano Zampini computesolvers = PETSC_TRUE; 808674ae819SStefano Zampini } 809674ae819SStefano Zampini } else { 810674ae819SStefano Zampini computeis = PETSC_TRUE; 811674ae819SStefano Zampini computetopography = PETSC_TRUE; 812674ae819SStefano Zampini computesolvers = PETSC_TRUE; 813674ae819SStefano Zampini } 814fcd409f5SStefano Zampini /* Set up all the "iterative substructuring" common block without computing solvers */ 815674ae819SStefano Zampini if (computeis) { 816fcd409f5SStefano Zampini /* HACK INTO PCIS */ 817fcd409f5SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 818fcd409f5SStefano Zampini pcis->computesolvers = PETSC_FALSE; 819674ae819SStefano Zampini ierr = PCISSetUp(pc);CHKERRQ(ierr); 820674ae819SStefano Zampini } 821674ae819SStefano Zampini /* Analyze interface and set up local constraint and change of basis matrices */ 822674ae819SStefano Zampini if (computetopography) { 823674ae819SStefano Zampini /* reset data */ 824674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 825674ae819SStefano Zampini ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 826674ae819SStefano Zampini ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 827674ae819SStefano Zampini } 828674ae819SStefano Zampini if (computesolvers) { 829674ae819SStefano Zampini /* reset data */ 830674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 831674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 83299cc7994SStefano Zampini /* Create coarse and local stuffs */ 83399cc7994SStefano Zampini ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 834674ae819SStefano Zampini ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 8350c7d97c5SJed Brown } 8362b510759SStefano Zampini if (pcbddc->dbg_flag && pcbddc->current_level) { 8372b510759SStefano Zampini ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 8382b510759SStefano Zampini } 8390c7d97c5SJed Brown PetscFunctionReturn(0); 8400c7d97c5SJed Brown } 8410c7d97c5SJed Brown 8420c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 8430c7d97c5SJed Brown /* 8440c7d97c5SJed Brown PCApply_BDDC - Applies the BDDC preconditioner to a vector. 8450c7d97c5SJed Brown 8460c7d97c5SJed Brown Input Parameters: 8470c7d97c5SJed Brown . pc - the preconditioner context 8480c7d97c5SJed Brown . r - input vector (global) 8490c7d97c5SJed Brown 8500c7d97c5SJed Brown Output Parameter: 8510c7d97c5SJed Brown . z - output vector (global) 8520c7d97c5SJed Brown 8530c7d97c5SJed Brown Application Interface Routine: PCApply() 8540c7d97c5SJed Brown */ 8550c7d97c5SJed Brown #undef __FUNCT__ 8560c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC" 85753cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 8580c7d97c5SJed Brown { 8590c7d97c5SJed Brown PC_IS *pcis = (PC_IS*)(pc->data); 8600c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 8610c7d97c5SJed Brown PetscErrorCode ierr; 8623b03a366Sstefano_zampini const PetscScalar one = 1.0; 8633b03a366Sstefano_zampini const PetscScalar m_one = -1.0; 8642617d88aSStefano Zampini const PetscScalar zero = 0.0; 8650c7d97c5SJed Brown 8660c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN 8670c7d97c5SJed Brown NN interface preconditioner changed to BDDC 8688eeda7d8SStefano Zampini Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 8690c7d97c5SJed Brown 8700c7d97c5SJed Brown PetscFunctionBegin; 871*85c4d303SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 8720c7d97c5SJed Brown /* First Dirichlet solve */ 8730c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8740c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 87553cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 8760c7d97c5SJed Brown /* 8770c7d97c5SJed Brown Assembling right hand side for BDDC operator 878674ae819SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 879674ae819SStefano Zampini - pcis->vec1_B the interface part of the global vector z 8800c7d97c5SJed Brown */ 8810c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 8820c7d97c5SJed Brown ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 8838eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 8840c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 8850c7d97c5SJed Brown ierr = VecCopy(r,z);CHKERRQ(ierr); 8860c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8870c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 888674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 889b76ba322SStefano Zampini } else { 8900bdf917eSStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 891b76ba322SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 892674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 893b76ba322SStefano Zampini } 894b76ba322SStefano Zampini 8952617d88aSStefano Zampini /* Apply interface preconditioner 8962617d88aSStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 8972617d88aSStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr); 8982617d88aSStefano Zampini 899674ae819SStefano Zampini /* Apply transpose of partition of unity operator */ 900674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 9010c7d97c5SJed Brown 9023b03a366Sstefano_zampini /* Second Dirichlet solve and assembling of output */ 9030c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9040c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9050c7d97c5SJed Brown ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 9068eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 90753cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr); 9080c7d97c5SJed Brown ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr); 9098eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 9100c7d97c5SJed Brown ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr); 9110c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9120c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9130c7d97c5SJed Brown PetscFunctionReturn(0); 9140c7d97c5SJed Brown } 915da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 916674ae819SStefano Zampini 917da1bb401SStefano Zampini #undef __FUNCT__ 918da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC" 919da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc) 920da1bb401SStefano Zampini { 921da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 922da1bb401SStefano Zampini PetscErrorCode ierr; 923da1bb401SStefano Zampini 924da1bb401SStefano Zampini PetscFunctionBegin; 925da1bb401SStefano Zampini /* free data created by PCIS */ 926da1bb401SStefano Zampini ierr = PCISDestroy(pc);CHKERRQ(ierr); 927674ae819SStefano Zampini /* free BDDC custom data */ 928674ae819SStefano Zampini ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 929674ae819SStefano Zampini /* destroy objects related to topography */ 930674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 931674ae819SStefano Zampini /* free allocated graph structure */ 932da1bb401SStefano Zampini ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 933674ae819SStefano Zampini /* free data for scaling operator */ 934674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 935674ae819SStefano Zampini /* free solvers stuff */ 936674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 93733bc96a4SStefano Zampini ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr); 93833bc96a4SStefano Zampini ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr); 9399881197aSStefano Zampini ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr); 940ac78edfcSStefano Zampini ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 94162a6ff1dSStefano Zampini /* free global vectors needed in presolve */ 94262a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 94362a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 9443425bc38SStefano Zampini /* remove functions */ 945674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 946bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 9472b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 948b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 9492b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 950bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 951bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 952bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 953bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 954bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 955bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 956bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 957bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 958bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 959bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 960674ae819SStefano Zampini /* Free the private data structure */ 961674ae819SStefano Zampini ierr = PetscFree(pc->data);CHKERRQ(ierr); 962da1bb401SStefano Zampini PetscFunctionReturn(0); 963da1bb401SStefano Zampini } 9643425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 9651e6b0712SBarry Smith 9663425bc38SStefano Zampini #undef __FUNCT__ 9673425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 9683425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 9693425bc38SStefano Zampini { 970674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 9713425bc38SStefano Zampini PC_IS* pcis; 9723425bc38SStefano Zampini PC_BDDC* pcbddc; 9733425bc38SStefano Zampini PetscErrorCode ierr; 9740c7d97c5SJed Brown 9753425bc38SStefano Zampini PetscFunctionBegin; 9763425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 9773425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 9783425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 9793425bc38SStefano Zampini 9803425bc38SStefano Zampini /* change of basis for physical rhs if needed 9813425bc38SStefano Zampini It also changes the rhs in case of dirichlet boundaries */ 9823308cffdSStefano Zampini ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 9833425bc38SStefano Zampini /* store vectors for computation of fetidp final solution */ 9843425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9853425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 986fb223d50SStefano Zampini /* scale rhs since it should be unassembled */ 987fb223d50SStefano Zampini /* TODO use counter scaling? (also below) */ 9883425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9893425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 990674ae819SStefano Zampini /* Apply partition of unity */ 9913425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 992674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 9938eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 9943425bc38SStefano Zampini /* compute partially subassembled Schur complement right-hand side */ 9953425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 9963425bc38SStefano Zampini ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 9973425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 9983425bc38SStefano Zampini ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 9993425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10003425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1001674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 10023425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10033425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10043425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 10053425bc38SStefano Zampini } 10063425bc38SStefano Zampini /* BDDC rhs */ 10073425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 10088eeda7d8SStefano Zampini if (pcbddc->switch_static) { 10093425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10103425bc38SStefano Zampini } 10113425bc38SStefano Zampini /* apply BDDC */ 10123425bc38SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 10133425bc38SStefano Zampini /* Application of B_delta and assembling of rhs for fetidp fluxes */ 10143425bc38SStefano Zampini ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 10153425bc38SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 10163425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10173425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10183425bc38SStefano Zampini /* restore original rhs */ 10193425bc38SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 10203425bc38SStefano Zampini PetscFunctionReturn(0); 10213425bc38SStefano Zampini } 10221e6b0712SBarry Smith 10233425bc38SStefano Zampini #undef __FUNCT__ 10243425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 10253425bc38SStefano Zampini /*@ 10263425bc38SStefano Zampini PCBDDCMatFETIDPGetRHS - Get rhs for FETIDP linear system. 10273425bc38SStefano Zampini 10283425bc38SStefano Zampini Collective 10293425bc38SStefano Zampini 10303425bc38SStefano Zampini Input Parameters: 10313425bc38SStefano Zampini + fetidp_mat - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators 10323425bc38SStefano Zampini + standard_rhs - the rhs of your linear system 10333425bc38SStefano Zampini 10343425bc38SStefano Zampini Output Parameters: 10353425bc38SStefano Zampini + fetidp_flux_rhs - the rhs of the FETIDP linear system 10363425bc38SStefano Zampini 10373425bc38SStefano Zampini Level: developer 10383425bc38SStefano Zampini 10393425bc38SStefano Zampini Notes: 10403425bc38SStefano Zampini 10413425bc38SStefano Zampini .seealso: PCBDDC 10423425bc38SStefano Zampini @*/ 10433425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 10443425bc38SStefano Zampini { 1045674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10463425bc38SStefano Zampini PetscErrorCode ierr; 10473425bc38SStefano Zampini 10483425bc38SStefano Zampini PetscFunctionBegin; 10493425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10503425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 10513425bc38SStefano Zampini PetscFunctionReturn(0); 10523425bc38SStefano Zampini } 10533425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 10541e6b0712SBarry Smith 10553425bc38SStefano Zampini #undef __FUNCT__ 10563425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 10573425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 10583425bc38SStefano Zampini { 1059674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10603425bc38SStefano Zampini PC_IS* pcis; 10613425bc38SStefano Zampini PC_BDDC* pcbddc; 10623425bc38SStefano Zampini PetscErrorCode ierr; 10633425bc38SStefano Zampini 10643425bc38SStefano Zampini PetscFunctionBegin; 10653425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10663425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 10673425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 10683425bc38SStefano Zampini 10693425bc38SStefano Zampini /* apply B_delta^T */ 10703425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10713425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10723425bc38SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 10733425bc38SStefano Zampini /* compute rhs for BDDC application */ 10743425bc38SStefano Zampini ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 10758eeda7d8SStefano Zampini if (pcbddc->switch_static) { 10763425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10773425bc38SStefano Zampini } 10783425bc38SStefano Zampini /* apply BDDC */ 10793425bc38SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 10803425bc38SStefano Zampini /* put values into standard global vector */ 10813425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10823425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10838eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 10843425bc38SStefano Zampini /* compute values into the interior if solved for the partially subassembled Schur complement */ 10853425bc38SStefano Zampini ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 10863425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 10873425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10883425bc38SStefano Zampini } 10893425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10903425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10913425bc38SStefano Zampini /* final change of basis if needed 10923425bc38SStefano Zampini Is also sums the dirichlet part removed during RHS assembling */ 10933308cffdSStefano Zampini ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 10943425bc38SStefano Zampini PetscFunctionReturn(0); 10953425bc38SStefano Zampini } 10961e6b0712SBarry Smith 10973425bc38SStefano Zampini #undef __FUNCT__ 10983425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 10993425bc38SStefano Zampini /*@ 11003425bc38SStefano Zampini PCBDDCMatFETIDPGetSolution - Get Solution for FETIDP linear system. 11013425bc38SStefano Zampini 11023425bc38SStefano Zampini Collective 11033425bc38SStefano Zampini 11043425bc38SStefano Zampini Input Parameters: 11053425bc38SStefano Zampini + fetidp_mat - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators 11063425bc38SStefano Zampini + fetidp_flux_sol - the solution of the FETIDP linear system 11073425bc38SStefano Zampini 11083425bc38SStefano Zampini Output Parameters: 11093425bc38SStefano Zampini + standard_sol - the solution on the global domain 11103425bc38SStefano Zampini 11113425bc38SStefano Zampini Level: developer 11123425bc38SStefano Zampini 11133425bc38SStefano Zampini Notes: 11143425bc38SStefano Zampini 11153425bc38SStefano Zampini .seealso: PCBDDC 11163425bc38SStefano Zampini @*/ 11173425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 11183425bc38SStefano Zampini { 1119674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 11203425bc38SStefano Zampini PetscErrorCode ierr; 11213425bc38SStefano Zampini 11223425bc38SStefano Zampini PetscFunctionBegin; 11233425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 11243425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 11253425bc38SStefano Zampini PetscFunctionReturn(0); 11263425bc38SStefano Zampini } 11273425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 11281e6b0712SBarry Smith 1129f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1130f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1131f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1132f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1133674ae819SStefano Zampini 11343425bc38SStefano Zampini #undef __FUNCT__ 11353425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 11363425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 11373425bc38SStefano Zampini { 1138674ae819SStefano Zampini 1139674ae819SStefano Zampini FETIDPMat_ctx fetidpmat_ctx; 11403425bc38SStefano Zampini Mat newmat; 1141674ae819SStefano Zampini FETIDPPC_ctx fetidppc_ctx; 11423425bc38SStefano Zampini PC newpc; 1143ce94432eSBarry Smith MPI_Comm comm; 11443425bc38SStefano Zampini PetscErrorCode ierr; 11453425bc38SStefano Zampini 11463425bc38SStefano Zampini PetscFunctionBegin; 1147ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 11483425bc38SStefano Zampini /* FETIDP linear matrix */ 11493425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 11503425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 11513425bc38SStefano Zampini ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 11523425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 11533425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 11543425bc38SStefano Zampini ierr = MatSetUp(newmat);CHKERRQ(ierr); 11553425bc38SStefano Zampini /* FETIDP preconditioner */ 11563425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 11573425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 11583425bc38SStefano Zampini ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 11593425bc38SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 11603425bc38SStefano Zampini ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 11613425bc38SStefano Zampini ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 11623425bc38SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 11633425bc38SStefano Zampini ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr); 11643425bc38SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 11653425bc38SStefano Zampini /* return pointers for objects created */ 11663425bc38SStefano Zampini *fetidp_mat=newmat; 11673425bc38SStefano Zampini *fetidp_pc=newpc; 11683425bc38SStefano Zampini PetscFunctionReturn(0); 11693425bc38SStefano Zampini } 11701e6b0712SBarry Smith 11713425bc38SStefano Zampini #undef __FUNCT__ 11723425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 11733425bc38SStefano Zampini /*@ 11743425bc38SStefano Zampini PCBDDCCreateFETIDPOperators - Create operators for FETIDP. 11753425bc38SStefano Zampini 11763425bc38SStefano Zampini Collective 11773425bc38SStefano Zampini 11783425bc38SStefano Zampini Input Parameters: 11793425bc38SStefano Zampini + pc - the BDDC preconditioning context (setup must be already called) 11803425bc38SStefano Zampini 11813425bc38SStefano Zampini Level: developer 11823425bc38SStefano Zampini 11833425bc38SStefano Zampini Notes: 11843425bc38SStefano Zampini 11853425bc38SStefano Zampini .seealso: PCBDDC 11863425bc38SStefano Zampini @*/ 11873425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 11883425bc38SStefano Zampini { 11893425bc38SStefano Zampini PetscErrorCode ierr; 11903425bc38SStefano Zampini 11913425bc38SStefano Zampini PetscFunctionBegin; 11923425bc38SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11933425bc38SStefano Zampini if (pc->setupcalled) { 1194516d51a7SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1195f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 11963425bc38SStefano Zampini PetscFunctionReturn(0); 11973425bc38SStefano Zampini } 11980c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 1199da1bb401SStefano Zampini /*MC 1200da1bb401SStefano Zampini PCBDDC - Balancing Domain Decomposition by Constraints. 12010c7d97c5SJed Brown 1202da1bb401SStefano Zampini Options Database Keys: 1203da1bb401SStefano Zampini . -pcbddc ??? - 1204da1bb401SStefano Zampini 1205da1bb401SStefano Zampini Level: intermediate 1206da1bb401SStefano Zampini 1207da1bb401SStefano Zampini Notes: The matrix used with this preconditioner must be of type MATIS 1208da1bb401SStefano Zampini 1209da1bb401SStefano Zampini Unlike more 'conventional' interface preconditioners, this iterates over ALL the 1210da1bb401SStefano Zampini degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers 1211da1bb401SStefano Zampini on the subdomains). 1212da1bb401SStefano Zampini 1213da1bb401SStefano Zampini Options for the coarse grid preconditioner can be set with - 1214da1bb401SStefano Zampini Options for the Dirichlet subproblem can be set with - 1215da1bb401SStefano Zampini Options for the Neumann subproblem can be set with - 1216da1bb401SStefano Zampini 1217da1bb401SStefano Zampini Contributed by Stefano Zampini 1218da1bb401SStefano Zampini 1219da1bb401SStefano Zampini .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1220da1bb401SStefano Zampini M*/ 1221b2573a8aSBarry Smith 1222da1bb401SStefano Zampini #undef __FUNCT__ 1223da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC" 12248cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1225da1bb401SStefano Zampini { 1226da1bb401SStefano Zampini PetscErrorCode ierr; 1227da1bb401SStefano Zampini PC_BDDC *pcbddc; 1228da1bb401SStefano Zampini 1229da1bb401SStefano Zampini PetscFunctionBegin; 1230da1bb401SStefano Zampini /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1231da1bb401SStefano Zampini ierr = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr); 1232da1bb401SStefano Zampini pc->data = (void*)pcbddc; 1233da1bb401SStefano Zampini 1234da1bb401SStefano Zampini /* create PCIS data structure */ 1235da1bb401SStefano Zampini ierr = PCISCreate(pc);CHKERRQ(ierr); 1236da1bb401SStefano Zampini 123747d04d0dSStefano Zampini /* BDDC customization */ 123847d04d0dSStefano Zampini pcbddc->use_vertices = PETSC_TRUE; 123947d04d0dSStefano Zampini pcbddc->use_edges = PETSC_TRUE; 124047d04d0dSStefano Zampini pcbddc->use_faces = PETSC_FALSE; 124147d04d0dSStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 124247d04d0dSStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 124347d04d0dSStefano Zampini pcbddc->switch_static = PETSC_FALSE; 124447d04d0dSStefano Zampini pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 124547d04d0dSStefano Zampini pcbddc->dbg_flag = 0; 124647d04d0dSStefano Zampini 1247674ae819SStefano Zampini pcbddc->user_primal_vertices = 0; 12480bdf917eSStefano Zampini pcbddc->NullSpace = 0; 12493972b0daSStefano Zampini pcbddc->temp_solution = 0; 1250534831adSStefano Zampini pcbddc->original_rhs = 0; 1251534831adSStefano Zampini pcbddc->local_mat = 0; 1252534831adSStefano Zampini pcbddc->ChangeOfBasisMatrix = 0; 1253da1bb401SStefano Zampini pcbddc->coarse_vec = 0; 1254da1bb401SStefano Zampini pcbddc->coarse_rhs = 0; 1255da1bb401SStefano Zampini pcbddc->coarse_ksp = 0; 1256da1bb401SStefano Zampini pcbddc->coarse_phi_B = 0; 1257da1bb401SStefano Zampini pcbddc->coarse_phi_D = 0; 125815aaf578SStefano Zampini pcbddc->coarse_psi_B = 0; 125915aaf578SStefano Zampini pcbddc->coarse_psi_D = 0; 1260da1bb401SStefano Zampini pcbddc->vec1_P = 0; 1261da1bb401SStefano Zampini pcbddc->vec1_R = 0; 1262da1bb401SStefano Zampini pcbddc->vec2_R = 0; 1263da1bb401SStefano Zampini pcbddc->local_auxmat1 = 0; 1264da1bb401SStefano Zampini pcbddc->local_auxmat2 = 0; 1265da1bb401SStefano Zampini pcbddc->R_to_B = 0; 1266da1bb401SStefano Zampini pcbddc->R_to_D = 0; 1267da1bb401SStefano Zampini pcbddc->ksp_D = 0; 1268da1bb401SStefano Zampini pcbddc->ksp_R = 0; 1269da1bb401SStefano Zampini pcbddc->NeumannBoundaries = 0; 1270da1bb401SStefano Zampini pcbddc->ISForDofs = 0; 1271da1bb401SStefano Zampini pcbddc->ConstraintMatrix = 0; 1272*85c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 127347d04d0dSStefano Zampini pcbddc->coarse_loc_to_glob = 0; 127447d04d0dSStefano Zampini pcbddc->coarsening_ratio = 8; 12754fad6a16SStefano Zampini pcbddc->current_level = 0; 12762b510759SStefano Zampini pcbddc->max_levels = 0; 1277da1bb401SStefano Zampini 1278674ae819SStefano Zampini /* create local graph structure */ 1279674ae819SStefano Zampini ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1280674ae819SStefano Zampini 1281674ae819SStefano Zampini /* scaling */ 1282674ae819SStefano Zampini pcbddc->use_deluxe_scaling = PETSC_FALSE; 1283674ae819SStefano Zampini pcbddc->work_scaling = 0; 1284da1bb401SStefano Zampini 1285da1bb401SStefano Zampini /* function pointers */ 1286da1bb401SStefano Zampini pc->ops->apply = PCApply_BDDC; 1287da1bb401SStefano Zampini pc->ops->applytranspose = 0; 1288da1bb401SStefano Zampini pc->ops->setup = PCSetUp_BDDC; 1289da1bb401SStefano Zampini pc->ops->destroy = PCDestroy_BDDC; 1290da1bb401SStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1291da1bb401SStefano Zampini pc->ops->view = 0; 1292da1bb401SStefano Zampini pc->ops->applyrichardson = 0; 1293da1bb401SStefano Zampini pc->ops->applysymmetricleft = 0; 1294da1bb401SStefano Zampini pc->ops->applysymmetricright = 0; 1295534831adSStefano Zampini pc->ops->presolve = PCPreSolve_BDDC; 1296534831adSStefano Zampini pc->ops->postsolve = PCPostSolve_BDDC; 1297da1bb401SStefano Zampini 1298da1bb401SStefano Zampini /* composing function */ 1299674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1300bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 13012b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1302b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 13032b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1304bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1305bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1306bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1307bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1308bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1309bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 1310bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1311bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1312bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1313bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1314da1bb401SStefano Zampini PetscFunctionReturn(0); 1315da1bb401SStefano Zampini } 13163425bc38SStefano Zampini 1317