153cdbc3dSStefano Zampini /* TODOLIST 2af732b37SStefano Zampini Better management for BAIJ local mats. How to deal with SBAIJ? 385c4d303SStefano Zampini Provide PCApplyTranpose 4af732b37SStefano Zampini make runexe59 5be4bc593SStefano Zampini Man pages 6*28509bceSStefano Zampini Propagate nearnullspace info among levels 7be4bc593SStefano Zampini Move FETIDP code 8be4bc593SStefano Zampini Provide general case for subassembling 9be4bc593SStefano Zampini Preallocation routines in MatConvert_IS_AIJ 10*28509bceSStefano Zampini Allow different customizations among different linear solves (requires also reset/destroy of ksp_R and coarse_ksp) 11be4bc593SStefano Zampini Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver? 12be4bc593SStefano Zampini Better management in PCIS code 1385c4d303SStefano Zampini Is it possible working with PCBDDCGraph on boundary indices only? 14da1bb401SStefano Zampini DofSplitting and DM attached to pc? 15da1bb401SStefano Zampini Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet) 16b8ffe317SStefano Zampini BDDC with MG framework? 1753cdbc3dSStefano Zampini */ 180c7d97c5SJed Brown 1953cdbc3dSStefano Zampini /* ---------------------------------------------------------------------------------------------------------------------------------------------- 200c7d97c5SJed Brown Implementation of BDDC preconditioner based on: 210c7d97c5SJed Brown C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 2253cdbc3dSStefano Zampini ---------------------------------------------------------------------------------------------------------------------------------------------- */ 2353cdbc3dSStefano Zampini 24674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 25674ae819SStefano Zampini #include "bddcprivate.h" 263b03a366Sstefano_zampini #include <petscblaslapack.h> 27674ae819SStefano Zampini 280c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 290c7d97c5SJed Brown #undef __FUNCT__ 300c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC" 310c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc) 320c7d97c5SJed Brown { 330c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 340c7d97c5SJed Brown PetscErrorCode ierr; 350c7d97c5SJed Brown 360c7d97c5SJed Brown PetscFunctionBegin; 370c7d97c5SJed Brown ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 388eeda7d8SStefano Zampini /* Verbose debugging */ 398eeda7d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 408eeda7d8SStefano Zampini /* Primal space cumstomization */ 418eeda7d8SStefano 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); 428eeda7d8SStefano 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); 438eeda7d8SStefano 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); 448eeda7d8SStefano Zampini /* Change of basis */ 458eeda7d8SStefano 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); 468eeda7d8SStefano 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); 47674ae819SStefano Zampini if (!pcbddc->use_change_of_basis) { 48674ae819SStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 49674ae819SStefano Zampini } 508eeda7d8SStefano Zampini /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 518eeda7d8SStefano 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); 520298fd71SBarry Smith ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 532b510759SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 54674ae819SStefano 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); 550c7d97c5SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 560c7d97c5SJed Brown PetscFunctionReturn(0); 570c7d97c5SJed Brown } 580c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 59674ae819SStefano Zampini #undef __FUNCT__ 60674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 61674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 62674ae819SStefano Zampini { 63674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 64674ae819SStefano Zampini PetscErrorCode ierr; 651e6b0712SBarry Smith 66674ae819SStefano Zampini PetscFunctionBegin; 67674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 68674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 69674ae819SStefano Zampini pcbddc->user_primal_vertices = PrimalVertices; 70674ae819SStefano Zampini PetscFunctionReturn(0); 71674ae819SStefano Zampini } 72674ae819SStefano Zampini #undef __FUNCT__ 73674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 74674ae819SStefano Zampini /*@ 75*28509bceSStefano Zampini PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 76674ae819SStefano Zampini 77674ae819SStefano Zampini Not collective 78674ae819SStefano Zampini 79674ae819SStefano Zampini Input Parameters: 80674ae819SStefano Zampini + pc - the preconditioning context 81*28509bceSStefano Zampini - PrimalVertices - index set of primal vertices in local numbering 82674ae819SStefano Zampini 83674ae819SStefano Zampini Level: intermediate 84674ae819SStefano Zampini 85674ae819SStefano Zampini Notes: 86674ae819SStefano Zampini 87674ae819SStefano Zampini .seealso: PCBDDC 88674ae819SStefano Zampini @*/ 89674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 90674ae819SStefano Zampini { 91674ae819SStefano Zampini PetscErrorCode ierr; 92674ae819SStefano Zampini 93674ae819SStefano Zampini PetscFunctionBegin; 94674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 95674ae819SStefano Zampini PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 96674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 97674ae819SStefano Zampini PetscFunctionReturn(0); 98674ae819SStefano Zampini } 99674ae819SStefano Zampini /* -------------------------------------------------------------------------- */ 1000c7d97c5SJed Brown #undef __FUNCT__ 1014fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 1024fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 1034fad6a16SStefano Zampini { 1044fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1054fad6a16SStefano Zampini 1064fad6a16SStefano Zampini PetscFunctionBegin; 1074fad6a16SStefano Zampini pcbddc->coarsening_ratio = k; 1084fad6a16SStefano Zampini PetscFunctionReturn(0); 1094fad6a16SStefano Zampini } 1101e6b0712SBarry Smith 1114fad6a16SStefano Zampini #undef __FUNCT__ 1124fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio" 1134fad6a16SStefano Zampini /*@ 114*28509bceSStefano Zampini PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 1154fad6a16SStefano Zampini 1164fad6a16SStefano Zampini Logically collective on PC 1174fad6a16SStefano Zampini 1184fad6a16SStefano Zampini Input Parameters: 1194fad6a16SStefano Zampini + pc - the preconditioning context 120*28509bceSStefano Zampini - k - coarsening ratio (H/h at the coarser level) 1214fad6a16SStefano Zampini 122*28509bceSStefano Zampini Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 1234fad6a16SStefano Zampini 1244fad6a16SStefano Zampini Level: intermediate 1254fad6a16SStefano Zampini 1264fad6a16SStefano Zampini Notes: 1274fad6a16SStefano Zampini 1284fad6a16SStefano Zampini .seealso: PCBDDC 1294fad6a16SStefano Zampini @*/ 1304fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 1314fad6a16SStefano Zampini { 1324fad6a16SStefano Zampini PetscErrorCode ierr; 1334fad6a16SStefano Zampini 1344fad6a16SStefano Zampini PetscFunctionBegin; 1354fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1362b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,k,2); 1374fad6a16SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 1384fad6a16SStefano Zampini PetscFunctionReturn(0); 1394fad6a16SStefano Zampini } 1402b510759SStefano Zampini 141b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 1422b510759SStefano Zampini #undef __FUNCT__ 143b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 144b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 145b8ffe317SStefano Zampini { 146b8ffe317SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 147b8ffe317SStefano Zampini 148b8ffe317SStefano Zampini PetscFunctionBegin; 14985c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = flg; 150b8ffe317SStefano Zampini PetscFunctionReturn(0); 151b8ffe317SStefano Zampini } 152b8ffe317SStefano Zampini 153b8ffe317SStefano Zampini #undef __FUNCT__ 154b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 155b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 1562b510759SStefano Zampini { 1572b510759SStefano Zampini PetscErrorCode ierr; 1582b510759SStefano Zampini 1592b510759SStefano Zampini PetscFunctionBegin; 1602b510759SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 161b8ffe317SStefano Zampini PetscValidLogicalCollectiveBool(pc,flg,2); 162b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1632b510759SStefano Zampini PetscFunctionReturn(0); 1642b510759SStefano Zampini } 1651e6b0712SBarry Smith 1664fad6a16SStefano Zampini #undef __FUNCT__ 1672b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC" 1682b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 1694fad6a16SStefano Zampini { 1704fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1714fad6a16SStefano Zampini 1724fad6a16SStefano Zampini PetscFunctionBegin; 1732b510759SStefano Zampini pcbddc->current_level = level; 1744fad6a16SStefano Zampini PetscFunctionReturn(0); 1754fad6a16SStefano Zampini } 1761e6b0712SBarry Smith 1774fad6a16SStefano Zampini #undef __FUNCT__ 178b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel" 179b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 180b8ffe317SStefano Zampini { 181b8ffe317SStefano Zampini PetscErrorCode ierr; 182b8ffe317SStefano Zampini 183b8ffe317SStefano Zampini PetscFunctionBegin; 184b8ffe317SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 185b8ffe317SStefano Zampini PetscValidLogicalCollectiveInt(pc,level,2); 186b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 187b8ffe317SStefano Zampini PetscFunctionReturn(0); 188b8ffe317SStefano Zampini } 189b8ffe317SStefano Zampini 190b8ffe317SStefano Zampini #undef __FUNCT__ 1912b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC" 1922b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 1932b510759SStefano Zampini { 1942b510759SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1952b510759SStefano Zampini 1962b510759SStefano Zampini PetscFunctionBegin; 1972b510759SStefano Zampini pcbddc->max_levels = levels; 1982b510759SStefano Zampini PetscFunctionReturn(0); 1992b510759SStefano Zampini } 2002b510759SStefano Zampini 2012b510759SStefano Zampini #undef __FUNCT__ 2022b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels" 2034fad6a16SStefano Zampini /*@ 204*28509bceSStefano Zampini PCBDDCSetLevels - Sets the maximum number of levels for multilevel 2054fad6a16SStefano Zampini 2064fad6a16SStefano Zampini Logically collective on PC 2074fad6a16SStefano Zampini 2084fad6a16SStefano Zampini Input Parameters: 2094fad6a16SStefano Zampini + pc - the preconditioning context 210*28509bceSStefano Zampini - levels - the maximum number of levels (max 9) 2114fad6a16SStefano Zampini 212*28509bceSStefano Zampini Default value is 0, i.e. traditional one-level BDDC 2134fad6a16SStefano Zampini 2144fad6a16SStefano Zampini Level: intermediate 2154fad6a16SStefano Zampini 2164fad6a16SStefano Zampini Notes: 2174fad6a16SStefano Zampini 2184fad6a16SStefano Zampini .seealso: PCBDDC 2194fad6a16SStefano Zampini @*/ 2202b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 2214fad6a16SStefano Zampini { 2224fad6a16SStefano Zampini PetscErrorCode ierr; 2234fad6a16SStefano Zampini 2244fad6a16SStefano Zampini PetscFunctionBegin; 2254fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2262b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,levels,2); 2272b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 2284fad6a16SStefano Zampini PetscFunctionReturn(0); 2294fad6a16SStefano Zampini } 2304fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */ 2311e6b0712SBarry Smith 2324fad6a16SStefano Zampini #undef __FUNCT__ 2330bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 2340bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 2350bdf917eSStefano Zampini { 2360bdf917eSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2370bdf917eSStefano Zampini PetscErrorCode ierr; 2380bdf917eSStefano Zampini 2390bdf917eSStefano Zampini PetscFunctionBegin; 2400bdf917eSStefano Zampini ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 2410bdf917eSStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 2420bdf917eSStefano Zampini pcbddc->NullSpace = NullSpace; 2430bdf917eSStefano Zampini PetscFunctionReturn(0); 2440bdf917eSStefano Zampini } 2451e6b0712SBarry Smith 2460bdf917eSStefano Zampini #undef __FUNCT__ 2470bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace" 2480bdf917eSStefano Zampini /*@ 249*28509bceSStefano Zampini PCBDDCSetNullSpace - Set nullspace for BDDC operator 2500bdf917eSStefano Zampini 2510bdf917eSStefano Zampini Logically collective on PC and MatNullSpace 2520bdf917eSStefano Zampini 2530bdf917eSStefano Zampini Input Parameters: 2540bdf917eSStefano Zampini + pc - the preconditioning context 255*28509bceSStefano Zampini - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 2560bdf917eSStefano Zampini 2570bdf917eSStefano Zampini Level: intermediate 2580bdf917eSStefano Zampini 2590bdf917eSStefano Zampini Notes: 2600bdf917eSStefano Zampini 2610bdf917eSStefano Zampini .seealso: PCBDDC 2620bdf917eSStefano Zampini @*/ 2630bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 2640bdf917eSStefano Zampini { 2650bdf917eSStefano Zampini PetscErrorCode ierr; 2660bdf917eSStefano Zampini 2670bdf917eSStefano Zampini PetscFunctionBegin; 2680bdf917eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 269674ae819SStefano Zampini PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 2702b510759SStefano Zampini PetscCheckSameComm(pc,1,NullSpace,2); 2710bdf917eSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 2720bdf917eSStefano Zampini PetscFunctionReturn(0); 2730bdf917eSStefano Zampini } 2740bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */ 2751e6b0712SBarry Smith 2760bdf917eSStefano Zampini #undef __FUNCT__ 2773b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 2783b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 2793b03a366Sstefano_zampini { 2803b03a366Sstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2813b03a366Sstefano_zampini PetscErrorCode ierr; 2823b03a366Sstefano_zampini 2833b03a366Sstefano_zampini PetscFunctionBegin; 2843b03a366Sstefano_zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 28536e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 28636e030ebSStefano Zampini pcbddc->DirichletBoundaries=DirichletBoundaries; 2873b03a366Sstefano_zampini PetscFunctionReturn(0); 2883b03a366Sstefano_zampini } 2891e6b0712SBarry Smith 2903b03a366Sstefano_zampini #undef __FUNCT__ 2913b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 2923b03a366Sstefano_zampini /*@ 293*28509bceSStefano Zampini PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 2943b03a366Sstefano_zampini 2953b03a366Sstefano_zampini Not collective 2963b03a366Sstefano_zampini 2973b03a366Sstefano_zampini Input Parameters: 2983b03a366Sstefano_zampini + pc - the preconditioning context 299*28509bceSStefano Zampini - DirichletBoundaries - sequential IS defining the subdomain part of Dirichlet boundaries (in local ordering) 3003b03a366Sstefano_zampini 3013b03a366Sstefano_zampini Level: intermediate 3023b03a366Sstefano_zampini 3033b03a366Sstefano_zampini Notes: 3043b03a366Sstefano_zampini 3053b03a366Sstefano_zampini .seealso: PCBDDC 3063b03a366Sstefano_zampini @*/ 3073b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 3083b03a366Sstefano_zampini { 3093b03a366Sstefano_zampini PetscErrorCode ierr; 3103b03a366Sstefano_zampini 3113b03a366Sstefano_zampini PetscFunctionBegin; 3123b03a366Sstefano_zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 313674ae819SStefano Zampini PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 3143b03a366Sstefano_zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 3153b03a366Sstefano_zampini PetscFunctionReturn(0); 3163b03a366Sstefano_zampini } 3173b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */ 3181e6b0712SBarry Smith 3193b03a366Sstefano_zampini #undef __FUNCT__ 3200c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 32153cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 3220c7d97c5SJed Brown { 3230c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 32453cdbc3dSStefano Zampini PetscErrorCode ierr; 3250c7d97c5SJed Brown 3260c7d97c5SJed Brown PetscFunctionBegin; 32753cdbc3dSStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 32836e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 32936e030ebSStefano Zampini pcbddc->NeumannBoundaries=NeumannBoundaries; 3300c7d97c5SJed Brown PetscFunctionReturn(0); 3310c7d97c5SJed Brown } 3321e6b0712SBarry Smith 3330c7d97c5SJed Brown #undef __FUNCT__ 3340c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 33557527edcSJed Brown /*@ 336*28509bceSStefano Zampini PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 33757527edcSJed Brown 3389c0446d6SStefano Zampini Not collective 33957527edcSJed Brown 34057527edcSJed Brown Input Parameters: 34157527edcSJed Brown + pc - the preconditioning context 342*28509bceSStefano Zampini - NeumannBoundaries - sequential IS defining the subdomain part of Neumann boundaries (in local ordering) 34357527edcSJed Brown 34457527edcSJed Brown Level: intermediate 34557527edcSJed Brown 34657527edcSJed Brown Notes: 34757527edcSJed Brown 34857527edcSJed Brown .seealso: PCBDDC 34957527edcSJed Brown @*/ 35053cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 3510c7d97c5SJed Brown { 3520c7d97c5SJed Brown PetscErrorCode ierr; 3530c7d97c5SJed Brown 3540c7d97c5SJed Brown PetscFunctionBegin; 3550c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 356674ae819SStefano Zampini PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 35753cdbc3dSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 35853cdbc3dSStefano Zampini PetscFunctionReturn(0); 35953cdbc3dSStefano Zampini } 36053cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 3611e6b0712SBarry Smith 36253cdbc3dSStefano Zampini #undef __FUNCT__ 363da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 364da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 365da1bb401SStefano Zampini { 366da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 367da1bb401SStefano Zampini 368da1bb401SStefano Zampini PetscFunctionBegin; 369da1bb401SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundaries; 370da1bb401SStefano Zampini PetscFunctionReturn(0); 371da1bb401SStefano Zampini } 3721e6b0712SBarry Smith 373da1bb401SStefano Zampini #undef __FUNCT__ 374da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 375da1bb401SStefano Zampini /*@ 376*28509bceSStefano Zampini PCBDDCGetDirichletBoundaries - Get IS for local Dirichlet boundaries 377da1bb401SStefano Zampini 378da1bb401SStefano Zampini Not collective 379da1bb401SStefano Zampini 380da1bb401SStefano Zampini Input Parameters: 381*28509bceSStefano Zampini . pc - the preconditioning context 382da1bb401SStefano Zampini 383da1bb401SStefano Zampini Output Parameters: 384*28509bceSStefano Zampini . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 385da1bb401SStefano Zampini 386da1bb401SStefano Zampini Level: intermediate 387da1bb401SStefano Zampini 388da1bb401SStefano Zampini Notes: 389da1bb401SStefano Zampini 390da1bb401SStefano Zampini .seealso: PCBDDC 391da1bb401SStefano Zampini @*/ 392da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 393da1bb401SStefano Zampini { 394da1bb401SStefano Zampini PetscErrorCode ierr; 395da1bb401SStefano Zampini 396da1bb401SStefano Zampini PetscFunctionBegin; 397da1bb401SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 398da1bb401SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 399da1bb401SStefano Zampini PetscFunctionReturn(0); 400da1bb401SStefano Zampini } 401da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 4021e6b0712SBarry Smith 403da1bb401SStefano Zampini #undef __FUNCT__ 40453cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 40553cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 40653cdbc3dSStefano Zampini { 40753cdbc3dSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 40853cdbc3dSStefano Zampini 40953cdbc3dSStefano Zampini PetscFunctionBegin; 41053cdbc3dSStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundaries; 41153cdbc3dSStefano Zampini PetscFunctionReturn(0); 41253cdbc3dSStefano Zampini } 4131e6b0712SBarry Smith 41453cdbc3dSStefano Zampini #undef __FUNCT__ 41553cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 41653cdbc3dSStefano Zampini /*@ 417*28509bceSStefano Zampini PCBDDCGetNeumannBoundaries - Get IS for local Neumann boundaries 41853cdbc3dSStefano Zampini 4199c0446d6SStefano Zampini Not collective 42053cdbc3dSStefano Zampini 42153cdbc3dSStefano Zampini Input Parameters: 422*28509bceSStefano Zampini . pc - the preconditioning context 42353cdbc3dSStefano Zampini 42453cdbc3dSStefano Zampini Output Parameters: 425*28509bceSStefano Zampini . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 42653cdbc3dSStefano Zampini 42753cdbc3dSStefano Zampini Level: intermediate 42853cdbc3dSStefano Zampini 42953cdbc3dSStefano Zampini Notes: 43053cdbc3dSStefano Zampini 43153cdbc3dSStefano Zampini .seealso: PCBDDC 43253cdbc3dSStefano Zampini @*/ 43353cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 43453cdbc3dSStefano Zampini { 43553cdbc3dSStefano Zampini PetscErrorCode ierr; 43653cdbc3dSStefano Zampini 43753cdbc3dSStefano Zampini PetscFunctionBegin; 43853cdbc3dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 43953cdbc3dSStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 4400c7d97c5SJed Brown PetscFunctionReturn(0); 4410c7d97c5SJed Brown } 44236e030ebSStefano Zampini /* -------------------------------------------------------------------------- */ 4431e6b0712SBarry Smith 44436e030ebSStefano Zampini #undef __FUNCT__ 445da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 4461a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 44736e030ebSStefano Zampini { 44836e030ebSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 449da1bb401SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 450da1bb401SStefano Zampini PetscErrorCode ierr; 45136e030ebSStefano Zampini 45236e030ebSStefano Zampini PetscFunctionBegin; 453674ae819SStefano Zampini /* free old CSR */ 454674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 455fb223d50SStefano Zampini /* TODO: PCBDDCGraphSetAdjacency */ 456674ae819SStefano Zampini /* get CSR into graph structure */ 457da1bb401SStefano Zampini if (copymode == PETSC_COPY_VALUES) { 458674ae819SStefano Zampini ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr); 459674ae819SStefano Zampini ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr); 460674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 461674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 462da1bb401SStefano Zampini } else if (copymode == PETSC_OWN_POINTER) { 4631a83f524SJed Brown mat_graph->xadj = (PetscInt*)xadj; 4641a83f524SJed Brown mat_graph->adjncy = (PetscInt*)adjncy; 465674ae819SStefano Zampini } 466575ad6abSStefano Zampini mat_graph->nvtxs_csr = nvtxs; 46736e030ebSStefano Zampini PetscFunctionReturn(0); 46836e030ebSStefano Zampini } 4691e6b0712SBarry Smith 47036e030ebSStefano Zampini #undef __FUNCT__ 471da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 47236e030ebSStefano Zampini /*@ 473*28509bceSStefano Zampini PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix 47436e030ebSStefano Zampini 47536e030ebSStefano Zampini Not collective 47636e030ebSStefano Zampini 47736e030ebSStefano Zampini Input Parameters: 47836e030ebSStefano Zampini + pc - the preconditioning context 479*28509bceSStefano Zampini . nvtxs - number of local vertices of the graph (i.e., the local size of your problem) 480*28509bceSStefano Zampini . xadj, adjncy - the CSR graph 481*28509bceSStefano Zampini - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 48236e030ebSStefano Zampini 48336e030ebSStefano Zampini Level: intermediate 48436e030ebSStefano Zampini 48536e030ebSStefano Zampini Notes: 48636e030ebSStefano Zampini 487*28509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode 48836e030ebSStefano Zampini @*/ 4891a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 49036e030ebSStefano Zampini { 491575ad6abSStefano Zampini void (*f)(void) = 0; 49236e030ebSStefano Zampini PetscErrorCode ierr; 49336e030ebSStefano Zampini 49436e030ebSStefano Zampini PetscFunctionBegin; 49536e030ebSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 496674ae819SStefano Zampini PetscValidIntPointer(xadj,3); 497674ae819SStefano Zampini PetscValidIntPointer(xadj,4); 498674ae819SStefano Zampini if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 499674ae819SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 500da1bb401SStefano Zampini } 50136e030ebSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 502575ad6abSStefano Zampini /* free arrays if PCBDDC is not the PC type */ 503575ad6abSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 504575ad6abSStefano Zampini if (!f && copymode == PETSC_OWN_POINTER) { 505575ad6abSStefano Zampini ierr = PetscFree(xadj);CHKERRQ(ierr); 506575ad6abSStefano Zampini ierr = PetscFree(adjncy);CHKERRQ(ierr); 50736e030ebSStefano Zampini } 50836e030ebSStefano Zampini PetscFunctionReturn(0); 50936e030ebSStefano Zampini } 5109c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */ 5111e6b0712SBarry Smith 5129c0446d6SStefano Zampini #undef __FUNCT__ 5139c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 5149c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 5159c0446d6SStefano Zampini { 5169c0446d6SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 5179c0446d6SStefano Zampini PetscInt i; 5189c0446d6SStefano Zampini PetscErrorCode ierr; 5199c0446d6SStefano Zampini 5209c0446d6SStefano Zampini PetscFunctionBegin; 521da1bb401SStefano Zampini /* Destroy ISes if they were already set */ 5229c0446d6SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 5239c0446d6SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 5249c0446d6SStefano Zampini } 525d11ae9bbSstefano_zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 526da1bb401SStefano Zampini /* allocate space then set */ 5279c0446d6SStefano Zampini ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr); 5289c0446d6SStefano Zampini for (i=0;i<n_is;i++) { 529da1bb401SStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 530da1bb401SStefano Zampini pcbddc->ISForDofs[i]=ISForDofs[i]; 5319c0446d6SStefano Zampini } 5329c0446d6SStefano Zampini pcbddc->n_ISForDofs=n_is; 5339c0446d6SStefano Zampini PetscFunctionReturn(0); 5349c0446d6SStefano Zampini } 5351e6b0712SBarry Smith 5369c0446d6SStefano Zampini #undef __FUNCT__ 5379c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting" 5389c0446d6SStefano Zampini /*@ 539*28509bceSStefano Zampini PCBDDCSetDofsSplitting - Set index sets defining fields of the local Neumann matrix 5409c0446d6SStefano Zampini 5419c0446d6SStefano Zampini Not collective 5429c0446d6SStefano Zampini 5439c0446d6SStefano Zampini Input Parameters: 5449c0446d6SStefano Zampini + pc - the preconditioning context 545*28509bceSStefano Zampini - n_is - number of index sets defining the fields 546*28509bceSStefano Zampini . ISForDofs - array of IS describing the fields 5479c0446d6SStefano Zampini 5489c0446d6SStefano Zampini Level: intermediate 5499c0446d6SStefano Zampini 5509c0446d6SStefano Zampini Notes: 5519c0446d6SStefano Zampini 5529c0446d6SStefano Zampini .seealso: PCBDDC 5539c0446d6SStefano Zampini @*/ 5549c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 5559c0446d6SStefano Zampini { 5562b510759SStefano Zampini PetscInt i; 5579c0446d6SStefano Zampini PetscErrorCode ierr; 5589c0446d6SStefano Zampini 5599c0446d6SStefano Zampini PetscFunctionBegin; 5609c0446d6SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5612b510759SStefano Zampini for (i=0;i<n_is;i++) { 5622b510759SStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2); 5632b510759SStefano Zampini } 5649c0446d6SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 5659c0446d6SStefano Zampini PetscFunctionReturn(0); 5669c0446d6SStefano Zampini } 567da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 568534831adSStefano Zampini #undef __FUNCT__ 569534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC" 570534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 571534831adSStefano Zampini /* 572534831adSStefano Zampini PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 573534831adSStefano Zampini guess if a transformation of basis approach has been selected. 5749c0446d6SStefano Zampini 575534831adSStefano Zampini Input Parameter: 576534831adSStefano Zampini + pc - the preconditioner contex 577534831adSStefano Zampini 578534831adSStefano Zampini Application Interface Routine: PCPreSolve() 579534831adSStefano Zampini 580534831adSStefano Zampini Notes: 581534831adSStefano Zampini The interface routine PCPreSolve() is not usually called directly by 582534831adSStefano Zampini the user, but instead is called by KSPSolve(). 583534831adSStefano Zampini */ 584534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 585534831adSStefano Zampini { 586534831adSStefano Zampini PetscErrorCode ierr; 587534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 588534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 589534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 590534831adSStefano Zampini Mat temp_mat; 5913972b0daSStefano Zampini IS dirIS; 5923972b0daSStefano Zampini PetscInt dirsize,i,*is_indices; 5933972b0daSStefano Zampini PetscScalar *array_x,*array_diagonal; 5943972b0daSStefano Zampini Vec used_vec; 59592e3dcfbSStefano Zampini PetscBool guess_nonzero,flg,bddc_has_dirichlet_boundaries; 596534831adSStefano Zampini 597534831adSStefano Zampini PetscFunctionBegin; 59885c4d303SStefano Zampini /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 59985c4d303SStefano Zampini if (ksp) { 60085c4d303SStefano Zampini PetscBool iscg; 60185c4d303SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 60285c4d303SStefano Zampini if (!iscg) { 60385c4d303SStefano Zampini ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 60485c4d303SStefano Zampini } 60585c4d303SStefano Zampini } 60685c4d303SStefano Zampini /* Creates parallel work vectors used in presolve */ 60762a6ff1dSStefano Zampini if (!pcbddc->original_rhs) { 60862a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 60962a6ff1dSStefano Zampini } 61062a6ff1dSStefano Zampini if (!pcbddc->temp_solution) { 61162a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 61262a6ff1dSStefano Zampini } 6133972b0daSStefano Zampini if (x) { 6143972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 6153972b0daSStefano Zampini used_vec = x; 6163972b0daSStefano Zampini } else { 6173972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 6183972b0daSStefano Zampini used_vec = pcbddc->temp_solution; 6193972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 6203972b0daSStefano Zampini } 6213972b0daSStefano Zampini /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 6223972b0daSStefano Zampini if (ksp) { 6233972b0daSStefano Zampini ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 6243972b0daSStefano Zampini if (!guess_nonzero) { 6253972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 6263972b0daSStefano Zampini } 6273972b0daSStefano Zampini } 6283308cffdSStefano Zampini 62992e3dcfbSStefano Zampini /* TODO: remove when Dirichlet boundaries will be shared */ 63092e3dcfbSStefano Zampini ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr); 63192e3dcfbSStefano Zampini flg = PETSC_FALSE; 63292e3dcfbSStefano Zampini if (dirIS) flg = PETSC_TRUE; 63392e3dcfbSStefano Zampini ierr = MPI_Allreduce(&flg,&bddc_has_dirichlet_boundaries,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 63492e3dcfbSStefano Zampini 6353972b0daSStefano Zampini /* store the original rhs */ 6363972b0daSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 6373972b0daSStefano Zampini 6383972b0daSStefano Zampini /* Take into account zeroed rows -> change rhs and store solution removed */ 63992e3dcfbSStefano Zampini if (rhs && bddc_has_dirichlet_boundaries) { 6403972b0daSStefano Zampini ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 6413972b0daSStefano Zampini ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 6423972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6433972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6443972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6453972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6463972b0daSStefano Zampini if (dirIS) { 6473972b0daSStefano Zampini ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr); 6483972b0daSStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 6493972b0daSStefano Zampini ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 6503972b0daSStefano Zampini ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 6512fa5cd67SKarl Rupp for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 6523972b0daSStefano Zampini ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 6533972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 6543972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 6553972b0daSStefano Zampini } 6563972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6573972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 658b76ba322SStefano Zampini 6593972b0daSStefano Zampini /* remove the computed solution from the rhs */ 6603972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 6613972b0daSStefano Zampini ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 6623972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 6633308cffdSStefano Zampini } 664b76ba322SStefano Zampini 665b76ba322SStefano Zampini /* store partially computed solution and set initial guess */ 6663972b0daSStefano Zampini if (x) { 6673972b0daSStefano Zampini ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 6683972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 66985c4d303SStefano Zampini if (pcbddc->use_exact_dirichlet_trick) { 670b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 671b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 672b76ba322SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 673b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 674b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 675b76ba322SStefano Zampini if (ksp) { 676b76ba322SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 677b76ba322SStefano Zampini } 678b76ba322SStefano Zampini } 6793972b0daSStefano Zampini } 680b76ba322SStefano Zampini 68192e3dcfbSStefano Zampini /* prepare MatMult and rhs for solver */ 682674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 683b76ba322SStefano Zampini /* swap pointers for local matrices */ 684b76ba322SStefano Zampini temp_mat = matis->A; 685b76ba322SStefano Zampini matis->A = pcbddc->local_mat; 686b76ba322SStefano Zampini pcbddc->local_mat = temp_mat; 68792e3dcfbSStefano Zampini if (rhs) { 688b76ba322SStefano Zampini /* Get local rhs and apply transformation of basis */ 689b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 690b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 691b76ba322SStefano Zampini /* from original basis to modified basis */ 692b76ba322SStefano Zampini ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 693b76ba322SStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 694b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 695b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 696674ae819SStefano Zampini } 69792e3dcfbSStefano Zampini } 69892e3dcfbSStefano Zampini 69992e3dcfbSStefano Zampini /* remove nullspace if present */ 7000bdf917eSStefano Zampini if (ksp && pcbddc->NullSpace) { 701d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 702d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 703b76ba322SStefano Zampini } 7040bdf917eSStefano Zampini ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 705534831adSStefano Zampini PetscFunctionReturn(0); 706534831adSStefano Zampini } 707534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 708534831adSStefano Zampini #undef __FUNCT__ 709534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC" 710534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 711534831adSStefano Zampini /* 712534831adSStefano Zampini PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 713534831adSStefano Zampini approach has been selected. Also, restores rhs to its original state. 714534831adSStefano Zampini 715534831adSStefano Zampini Input Parameter: 716534831adSStefano Zampini + pc - the preconditioner contex 717534831adSStefano Zampini 718534831adSStefano Zampini Application Interface Routine: PCPostSolve() 719534831adSStefano Zampini 720534831adSStefano Zampini Notes: 721534831adSStefano Zampini The interface routine PCPostSolve() is not usually called directly by 722534831adSStefano Zampini the user, but instead is called by KSPSolve(). 723534831adSStefano Zampini */ 724534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 725534831adSStefano Zampini { 726534831adSStefano Zampini PetscErrorCode ierr; 727534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 728534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 729534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 730534831adSStefano Zampini Mat temp_mat; 731534831adSStefano Zampini 732534831adSStefano Zampini PetscFunctionBegin; 733674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 734534831adSStefano Zampini /* swap pointers for local matrices */ 735534831adSStefano Zampini temp_mat = matis->A; 736534831adSStefano Zampini matis->A = pcbddc->local_mat; 737534831adSStefano Zampini pcbddc->local_mat = temp_mat; 7383425bc38SStefano Zampini } 7393308cffdSStefano Zampini if (pcbddc->use_change_of_basis && x) { 740534831adSStefano Zampini /* Get Local boundary and apply transformation of basis to solution vector */ 741534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 742534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 743534831adSStefano Zampini /* from modified basis to original basis */ 744534831adSStefano Zampini ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 745534831adSStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 746534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 747534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 748534831adSStefano Zampini } 7493972b0daSStefano Zampini /* add solution removed in presolve */ 7503425bc38SStefano Zampini if (x) { 7513425bc38SStefano Zampini ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 7523425bc38SStefano Zampini } 753fb223d50SStefano Zampini /* restore rhs to its original state */ 754fb223d50SStefano Zampini if (rhs) { 755fb223d50SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 756fb223d50SStefano Zampini } 757534831adSStefano Zampini PetscFunctionReturn(0); 758534831adSStefano Zampini } 759534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 76053cdbc3dSStefano Zampini #undef __FUNCT__ 76153cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC" 7620c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 7630c7d97c5SJed Brown /* 7640c7d97c5SJed Brown PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 7650c7d97c5SJed Brown by setting data structures and options. 7660c7d97c5SJed Brown 7670c7d97c5SJed Brown Input Parameter: 76853cdbc3dSStefano Zampini + pc - the preconditioner context 7690c7d97c5SJed Brown 7700c7d97c5SJed Brown Application Interface Routine: PCSetUp() 7710c7d97c5SJed Brown 7720c7d97c5SJed Brown Notes: 7730c7d97c5SJed Brown The interface routine PCSetUp() is not usually called directly by 7740c7d97c5SJed Brown the user, but instead is called by PCApply() if necessary. 7750c7d97c5SJed Brown */ 77653cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc) 7770c7d97c5SJed Brown { 7780c7d97c5SJed Brown PetscErrorCode ierr; 7790c7d97c5SJed Brown PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 780674ae819SStefano Zampini MatStructure flag; 781674ae819SStefano Zampini PetscBool computeis,computetopography,computesolvers; 7820c7d97c5SJed Brown 7830c7d97c5SJed Brown PetscFunctionBegin; 784674ae819SStefano Zampini /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */ 785fcd409f5SStefano Zampini /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */ 7863b03a366Sstefano_zampini /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 787fcd409f5SStefano Zampini Also, BDDC directly build the Dirichlet problem */ 7883b03a366Sstefano_zampini /* Get stdout for dbg */ 789674ae819SStefano Zampini if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) { 790ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr); 791e269702eSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 7922b510759SStefano Zampini if (pcbddc->current_level) { 7932b510759SStefano Zampini ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 7942b510759SStefano Zampini } 795e269702eSStefano Zampini } 796674ae819SStefano Zampini /* first attempt to split work */ 797674ae819SStefano Zampini if (pc->setupcalled) { 798674ae819SStefano Zampini computeis = PETSC_FALSE; 799674ae819SStefano Zampini ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr); 800674ae819SStefano Zampini if (flag == SAME_PRECONDITIONER) { 801674ae819SStefano Zampini computetopography = PETSC_FALSE; 802674ae819SStefano Zampini computesolvers = PETSC_FALSE; 803674ae819SStefano Zampini } else if (flag == SAME_NONZERO_PATTERN) { 804674ae819SStefano Zampini computetopography = PETSC_FALSE; 805674ae819SStefano Zampini computesolvers = PETSC_TRUE; 806674ae819SStefano Zampini } else { /* DIFFERENT_NONZERO_PATTERN */ 807674ae819SStefano Zampini computetopography = PETSC_TRUE; 808674ae819SStefano Zampini computesolvers = PETSC_TRUE; 809674ae819SStefano Zampini } 810674ae819SStefano Zampini } else { 811674ae819SStefano Zampini computeis = PETSC_TRUE; 812674ae819SStefano Zampini computetopography = PETSC_TRUE; 813674ae819SStefano Zampini computesolvers = PETSC_TRUE; 814674ae819SStefano Zampini } 815fcd409f5SStefano Zampini /* Set up all the "iterative substructuring" common block without computing solvers */ 816674ae819SStefano Zampini if (computeis) { 817fcd409f5SStefano Zampini /* HACK INTO PCIS */ 818fcd409f5SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 819fcd409f5SStefano Zampini pcis->computesolvers = PETSC_FALSE; 820674ae819SStefano Zampini ierr = PCISSetUp(pc);CHKERRQ(ierr); 821674ae819SStefano Zampini } 822674ae819SStefano Zampini /* Analyze interface and set up local constraint and change of basis matrices */ 823674ae819SStefano Zampini if (computetopography) { 824674ae819SStefano Zampini /* reset data */ 825674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 826674ae819SStefano Zampini ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 827674ae819SStefano Zampini ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 828674ae819SStefano Zampini } 829674ae819SStefano Zampini if (computesolvers) { 830674ae819SStefano Zampini /* reset data */ 831674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 832674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 83399cc7994SStefano Zampini /* Create coarse and local stuffs */ 83499cc7994SStefano Zampini ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 835674ae819SStefano Zampini ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 8360c7d97c5SJed Brown } 8372b510759SStefano Zampini if (pcbddc->dbg_flag && pcbddc->current_level) { 8382b510759SStefano Zampini ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 8392b510759SStefano Zampini } 8400c7d97c5SJed Brown PetscFunctionReturn(0); 8410c7d97c5SJed Brown } 8420c7d97c5SJed Brown 8430c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 8440c7d97c5SJed Brown /* 8450c7d97c5SJed Brown PCApply_BDDC - Applies the BDDC preconditioner to a vector. 8460c7d97c5SJed Brown 8470c7d97c5SJed Brown Input Parameters: 8480c7d97c5SJed Brown . pc - the preconditioner context 8490c7d97c5SJed Brown . r - input vector (global) 8500c7d97c5SJed Brown 8510c7d97c5SJed Brown Output Parameter: 8520c7d97c5SJed Brown . z - output vector (global) 8530c7d97c5SJed Brown 8540c7d97c5SJed Brown Application Interface Routine: PCApply() 8550c7d97c5SJed Brown */ 8560c7d97c5SJed Brown #undef __FUNCT__ 8570c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC" 85853cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 8590c7d97c5SJed Brown { 8600c7d97c5SJed Brown PC_IS *pcis = (PC_IS*)(pc->data); 8610c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 8620c7d97c5SJed Brown PetscErrorCode ierr; 8633b03a366Sstefano_zampini const PetscScalar one = 1.0; 8643b03a366Sstefano_zampini const PetscScalar m_one = -1.0; 8652617d88aSStefano Zampini const PetscScalar zero = 0.0; 8660c7d97c5SJed Brown 8670c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN 8680c7d97c5SJed Brown NN interface preconditioner changed to BDDC 8698eeda7d8SStefano Zampini Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 8700c7d97c5SJed Brown 8710c7d97c5SJed Brown PetscFunctionBegin; 87285c4d303SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 8730c7d97c5SJed Brown /* First Dirichlet solve */ 8740c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8750c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 87653cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 8770c7d97c5SJed Brown /* 8780c7d97c5SJed Brown Assembling right hand side for BDDC operator 879674ae819SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 880674ae819SStefano Zampini - pcis->vec1_B the interface part of the global vector z 8810c7d97c5SJed Brown */ 8820c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 8830c7d97c5SJed Brown ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 8848eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 8850c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 8860c7d97c5SJed Brown ierr = VecCopy(r,z);CHKERRQ(ierr); 8870c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8880c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 889674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 890b76ba322SStefano Zampini } else { 8910bdf917eSStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 892b76ba322SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 893674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 894b76ba322SStefano Zampini } 895b76ba322SStefano Zampini 8962617d88aSStefano Zampini /* Apply interface preconditioner 8972617d88aSStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 8982617d88aSStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr); 8992617d88aSStefano Zampini 900674ae819SStefano Zampini /* Apply transpose of partition of unity operator */ 901674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 9020c7d97c5SJed Brown 9033b03a366Sstefano_zampini /* Second Dirichlet solve and assembling of output */ 9040c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9050c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9060c7d97c5SJed Brown ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 9078eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 90853cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr); 9090c7d97c5SJed Brown ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr); 9108eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 9110c7d97c5SJed Brown ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr); 9120c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9130c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9140c7d97c5SJed Brown PetscFunctionReturn(0); 9150c7d97c5SJed Brown } 916da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 917674ae819SStefano Zampini 918da1bb401SStefano Zampini #undef __FUNCT__ 919da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC" 920da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc) 921da1bb401SStefano Zampini { 922da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 923da1bb401SStefano Zampini PetscErrorCode ierr; 924da1bb401SStefano Zampini 925da1bb401SStefano Zampini PetscFunctionBegin; 926da1bb401SStefano Zampini /* free data created by PCIS */ 927da1bb401SStefano Zampini ierr = PCISDestroy(pc);CHKERRQ(ierr); 928674ae819SStefano Zampini /* free BDDC custom data */ 929674ae819SStefano Zampini ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 930674ae819SStefano Zampini /* destroy objects related to topography */ 931674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 932674ae819SStefano Zampini /* free allocated graph structure */ 933da1bb401SStefano Zampini ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 934674ae819SStefano Zampini /* free data for scaling operator */ 935674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 936674ae819SStefano Zampini /* free solvers stuff */ 937674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 93833bc96a4SStefano Zampini ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr); 93933bc96a4SStefano Zampini ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr); 9409881197aSStefano Zampini ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr); 941ac78edfcSStefano Zampini ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 94262a6ff1dSStefano Zampini /* free global vectors needed in presolve */ 94362a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 94462a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 9453425bc38SStefano Zampini /* remove functions */ 946674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 947bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 9482b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 949b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 9502b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 951bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 952bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 953bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 954bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 955bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 956bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 957bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 958bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 959bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 960bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 961674ae819SStefano Zampini /* Free the private data structure */ 962674ae819SStefano Zampini ierr = PetscFree(pc->data);CHKERRQ(ierr); 963da1bb401SStefano Zampini PetscFunctionReturn(0); 964da1bb401SStefano Zampini } 9653425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 9661e6b0712SBarry Smith 9673425bc38SStefano Zampini #undef __FUNCT__ 9683425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 9693425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 9703425bc38SStefano Zampini { 971674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 9723425bc38SStefano Zampini PC_IS* pcis; 9733425bc38SStefano Zampini PC_BDDC* pcbddc; 9743425bc38SStefano Zampini PetscErrorCode ierr; 9750c7d97c5SJed Brown 9763425bc38SStefano Zampini PetscFunctionBegin; 9773425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 9783425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 9793425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 9803425bc38SStefano Zampini 9813425bc38SStefano Zampini /* change of basis for physical rhs if needed 9823425bc38SStefano Zampini It also changes the rhs in case of dirichlet boundaries */ 9833308cffdSStefano Zampini ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 9843425bc38SStefano Zampini /* store vectors for computation of fetidp final solution */ 9853425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9863425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 987fb223d50SStefano Zampini /* scale rhs since it should be unassembled */ 988fb223d50SStefano Zampini /* TODO use counter scaling? (also below) */ 9893425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9903425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 991674ae819SStefano Zampini /* Apply partition of unity */ 9923425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 993674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 9948eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 9953425bc38SStefano Zampini /* compute partially subassembled Schur complement right-hand side */ 9963425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 9973425bc38SStefano Zampini ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 9983425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 9993425bc38SStefano Zampini ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 10003425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10013425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1002674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 10033425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10043425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10053425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 10063425bc38SStefano Zampini } 10073425bc38SStefano Zampini /* BDDC rhs */ 10083425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 10098eeda7d8SStefano Zampini if (pcbddc->switch_static) { 10103425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10113425bc38SStefano Zampini } 10123425bc38SStefano Zampini /* apply BDDC */ 10133425bc38SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 10143425bc38SStefano Zampini /* Application of B_delta and assembling of rhs for fetidp fluxes */ 10153425bc38SStefano Zampini ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 10163425bc38SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 10173425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10183425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10193425bc38SStefano Zampini /* restore original rhs */ 10203425bc38SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 10213425bc38SStefano Zampini PetscFunctionReturn(0); 10223425bc38SStefano Zampini } 10231e6b0712SBarry Smith 10243425bc38SStefano Zampini #undef __FUNCT__ 10253425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 10263425bc38SStefano Zampini /*@ 1027*28509bceSStefano Zampini PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system 10283425bc38SStefano Zampini 10293425bc38SStefano Zampini Collective 10303425bc38SStefano Zampini 10313425bc38SStefano Zampini Input Parameters: 1032*28509bceSStefano Zampini + fetidp_mat - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators 1033*28509bceSStefano Zampini . standard_rhs - the right-hand side for your linear system 10343425bc38SStefano Zampini 10353425bc38SStefano Zampini Output Parameters: 1036*28509bceSStefano Zampini - fetidp_flux_rhs - the right-hand side for the FETIDP linear system 10373425bc38SStefano Zampini 10383425bc38SStefano Zampini Level: developer 10393425bc38SStefano Zampini 10403425bc38SStefano Zampini Notes: 10413425bc38SStefano Zampini 1042*28509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 10433425bc38SStefano Zampini @*/ 10443425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 10453425bc38SStefano Zampini { 1046674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10473425bc38SStefano Zampini PetscErrorCode ierr; 10483425bc38SStefano Zampini 10493425bc38SStefano Zampini PetscFunctionBegin; 10503425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10513425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 10523425bc38SStefano Zampini PetscFunctionReturn(0); 10533425bc38SStefano Zampini } 10543425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 10551e6b0712SBarry Smith 10563425bc38SStefano Zampini #undef __FUNCT__ 10573425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 10583425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 10593425bc38SStefano Zampini { 1060674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10613425bc38SStefano Zampini PC_IS* pcis; 10623425bc38SStefano Zampini PC_BDDC* pcbddc; 10633425bc38SStefano Zampini PetscErrorCode ierr; 10643425bc38SStefano Zampini 10653425bc38SStefano Zampini PetscFunctionBegin; 10663425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10673425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 10683425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 10693425bc38SStefano Zampini 10703425bc38SStefano Zampini /* apply B_delta^T */ 10713425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10723425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10733425bc38SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 10743425bc38SStefano Zampini /* compute rhs for BDDC application */ 10753425bc38SStefano Zampini ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 10768eeda7d8SStefano Zampini if (pcbddc->switch_static) { 10773425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10783425bc38SStefano Zampini } 10793425bc38SStefano Zampini /* apply BDDC */ 10803425bc38SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 10813425bc38SStefano Zampini /* put values into standard global vector */ 10823425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10833425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10848eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 10853425bc38SStefano Zampini /* compute values into the interior if solved for the partially subassembled Schur complement */ 10863425bc38SStefano Zampini ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 10873425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 10883425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10893425bc38SStefano Zampini } 10903425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10913425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10923425bc38SStefano Zampini /* final change of basis if needed 10933425bc38SStefano Zampini Is also sums the dirichlet part removed during RHS assembling */ 10943308cffdSStefano Zampini ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 10953425bc38SStefano Zampini PetscFunctionReturn(0); 10963425bc38SStefano Zampini } 10971e6b0712SBarry Smith 10983425bc38SStefano Zampini #undef __FUNCT__ 10993425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 11003425bc38SStefano Zampini /*@ 1101*28509bceSStefano Zampini PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system 11023425bc38SStefano Zampini 11033425bc38SStefano Zampini Collective 11043425bc38SStefano Zampini 11053425bc38SStefano Zampini Input Parameters: 1106*28509bceSStefano Zampini + fetidp_mat - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators 1107*28509bceSStefano Zampini . fetidp_flux_sol - the solution of the FETIDP linear system 11083425bc38SStefano Zampini 11093425bc38SStefano Zampini Output Parameters: 1110*28509bceSStefano Zampini - standard_sol - the solution defined on the physical domain 11113425bc38SStefano Zampini 11123425bc38SStefano Zampini Level: developer 11133425bc38SStefano Zampini 11143425bc38SStefano Zampini Notes: 11153425bc38SStefano Zampini 1116*28509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 11173425bc38SStefano Zampini @*/ 11183425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 11193425bc38SStefano Zampini { 1120674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 11213425bc38SStefano Zampini PetscErrorCode ierr; 11223425bc38SStefano Zampini 11233425bc38SStefano Zampini PetscFunctionBegin; 11243425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 11253425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 11263425bc38SStefano Zampini PetscFunctionReturn(0); 11273425bc38SStefano Zampini } 11283425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 11291e6b0712SBarry Smith 1130f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1131f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1132f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1133f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1134674ae819SStefano Zampini 11353425bc38SStefano Zampini #undef __FUNCT__ 11363425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 11373425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 11383425bc38SStefano Zampini { 1139674ae819SStefano Zampini 1140674ae819SStefano Zampini FETIDPMat_ctx fetidpmat_ctx; 11413425bc38SStefano Zampini Mat newmat; 1142674ae819SStefano Zampini FETIDPPC_ctx fetidppc_ctx; 11433425bc38SStefano Zampini PC newpc; 1144ce94432eSBarry Smith MPI_Comm comm; 11453425bc38SStefano Zampini PetscErrorCode ierr; 11463425bc38SStefano Zampini 11473425bc38SStefano Zampini PetscFunctionBegin; 1148ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 11493425bc38SStefano Zampini /* FETIDP linear matrix */ 11503425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 11513425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 11523425bc38SStefano Zampini ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 11533425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 11543425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 11553425bc38SStefano Zampini ierr = MatSetUp(newmat);CHKERRQ(ierr); 11563425bc38SStefano Zampini /* FETIDP preconditioner */ 11573425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 11583425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 11593425bc38SStefano Zampini ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 11603425bc38SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 11613425bc38SStefano Zampini ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 11623425bc38SStefano Zampini ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 11633425bc38SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 11643425bc38SStefano Zampini ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr); 11653425bc38SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 11663425bc38SStefano Zampini /* return pointers for objects created */ 11673425bc38SStefano Zampini *fetidp_mat=newmat; 11683425bc38SStefano Zampini *fetidp_pc=newpc; 11693425bc38SStefano Zampini PetscFunctionReturn(0); 11703425bc38SStefano Zampini } 11711e6b0712SBarry Smith 11723425bc38SStefano Zampini #undef __FUNCT__ 11733425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 11743425bc38SStefano Zampini /*@ 1175*28509bceSStefano Zampini PCBDDCCreateFETIDPOperators - Create operators for FETIDP 11763425bc38SStefano Zampini 11773425bc38SStefano Zampini Collective 11783425bc38SStefano Zampini 11793425bc38SStefano Zampini Input Parameters: 1180*28509bceSStefano Zampini + pc - the BDDC preconditioning context already setup 1181*28509bceSStefano Zampini 1182*28509bceSStefano Zampini Output Parameters: 1183*28509bceSStefano Zampini . fetidp_mat - shell FETIDP matrix object 1184*28509bceSStefano Zampini . fetidp_pc - shell Dirichlet preconditioner for FETIDP matrix 1185*28509bceSStefano Zampini 1186*28509bceSStefano Zampini Options Database Keys: 1187*28509bceSStefano Zampini - -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers 11883425bc38SStefano Zampini 11893425bc38SStefano Zampini Level: developer 11903425bc38SStefano Zampini 11913425bc38SStefano Zampini Notes: 1192*28509bceSStefano Zampini Currently the only operation provided for FETIDP matrix is MatMult 11933425bc38SStefano Zampini 11943425bc38SStefano Zampini .seealso: PCBDDC 11953425bc38SStefano Zampini @*/ 11963425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 11973425bc38SStefano Zampini { 11983425bc38SStefano Zampini PetscErrorCode ierr; 11993425bc38SStefano Zampini 12003425bc38SStefano Zampini PetscFunctionBegin; 12013425bc38SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12023425bc38SStefano Zampini if (pc->setupcalled) { 1203516d51a7SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1204f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 12053425bc38SStefano Zampini PetscFunctionReturn(0); 12063425bc38SStefano Zampini } 12070c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 1208da1bb401SStefano Zampini /*MC 1209da1bb401SStefano Zampini PCBDDC - Balancing Domain Decomposition by Constraints. 12100c7d97c5SJed Brown 1211*28509bceSStefano Zampini An implementation of the BDDC preconditioner based on 1212*28509bceSStefano Zampini 1213*28509bceSStefano Zampini .vb 1214*28509bceSStefano Zampini [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 1215*28509bceSStefano Zampini [2] A. Klawonn and O. B. Widlund. "Dual-Primal FETI Methods for Linear Elasticity", http://cs.nyu.edu/csweb/Research/TechReports/TR2004-855/TR2004-855.pdf 1216*28509bceSStefano Zampini [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 1217*28509bceSStefano Zampini .ve 1218*28509bceSStefano Zampini 1219*28509bceSStefano Zampini The matrix to be preconditioned (Pmat) must be of type MATIS. 1220*28509bceSStefano Zampini 1221*28509bceSStefano Zampini Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ or MATSEQBAIJ, either with real or complex numbers. 1222*28509bceSStefano Zampini 1223*28509bceSStefano Zampini Unlike 'conventional' interface preconditioners, PCBDDC iterates over all degrees of freedom, not just those on the interface. This allows the use of approximate solvers on the subdomains. 1224*28509bceSStefano Zampini 1225*28509bceSStefano Zampini It also works with unsymmetric and indefinite problems. 1226*28509bceSStefano Zampini 1227*28509bceSStefano Zampini Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace 1228*28509bceSStefano Zampini 1229*28509bceSStefano Zampini Boundary nodes are split in vertices, edges and faces using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph 1230*28509bceSStefano Zampini 1231*28509bceSStefano Zampini Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace. 1232*28509bceSStefano Zampini 1233*28509bceSStefano Zampini Change of basis is performed similarly to [2]. When more the one constraint is present on a single connected component (i.e. an edge or a face), a robust method based on local QR factorizations is used. 1234*28509bceSStefano Zampini 1235*28509bceSStefano Zampini The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object. 1236*28509bceSStefano Zampini 1237da1bb401SStefano Zampini Options Database Keys: 1238*28509bceSStefano Zampini 1239*28509bceSStefano Zampini . -pc_bddc_use_vertices <1> - use or not vertices in primal space 1240*28509bceSStefano Zampini . -pc_bddc_use_edges <1> - use or not edges in primal space 1241*28509bceSStefano Zampini . -pc_bddc_use_faces <1> - use or not faces in primal space 1242*28509bceSStefano Zampini . -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only) 1243*28509bceSStefano Zampini . -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested 1244*28509bceSStefano Zampini . -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1]) 1245*28509bceSStefano Zampini . -pc_bddc_levels <0> - maximum number of levels for multilevel 1246*28509bceSStefano Zampini . -pc_bddc_coarsening_ratio - H/h ratio at the coarser level 1247*28509bceSStefano Zampini - -pc_bddc_check_level <0> - set verbosity level of debugging output 1248*28509bceSStefano Zampini 1249*28509bceSStefano Zampini Options for Dirichlet, Neumann or coarse solver can be set with 1250*28509bceSStefano Zampini .vb 1251*28509bceSStefano Zampini -pc_bddc_dirichlet_ 1252*28509bceSStefano Zampini -pc_bddc_neumann_ 1253*28509bceSStefano Zampini -pc_bddc_coarse_ 1254*28509bceSStefano Zampini .ve 1255*28509bceSStefano Zampini e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg 1256*28509bceSStefano Zampini 1257*28509bceSStefano Zampini When using a multilevel approach, solvers' options at the N-th level can be specified as 1258*28509bceSStefano Zampini .vb 1259*28509bceSStefano Zampini -pc_bddc_dirichlet_N_ 1260*28509bceSStefano Zampini -pc_bddc_neumann_N_ 1261*28509bceSStefano Zampini -pc_bddc_coarse_N_ 1262*28509bceSStefano Zampini .ve 1263*28509bceSStefano Zampini Note that level number ranges from the finest 0 to the coarsest N 1264da1bb401SStefano Zampini 1265da1bb401SStefano Zampini Level: intermediate 1266da1bb401SStefano Zampini 1267*28509bceSStefano Zampini Notes: 1268*28509bceSStefano Zampini Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose 1269da1bb401SStefano Zampini 1270*28509bceSStefano Zampini New deluxe scaling operator will be available soon. 1271da1bb401SStefano Zampini 1272da1bb401SStefano Zampini Contributed by Stefano Zampini 1273da1bb401SStefano Zampini 1274da1bb401SStefano Zampini .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1275da1bb401SStefano Zampini M*/ 1276b2573a8aSBarry Smith 1277da1bb401SStefano Zampini #undef __FUNCT__ 1278da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC" 12798cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1280da1bb401SStefano Zampini { 1281da1bb401SStefano Zampini PetscErrorCode ierr; 1282da1bb401SStefano Zampini PC_BDDC *pcbddc; 1283da1bb401SStefano Zampini 1284da1bb401SStefano Zampini PetscFunctionBegin; 1285da1bb401SStefano Zampini /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1286da1bb401SStefano Zampini ierr = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr); 1287da1bb401SStefano Zampini pc->data = (void*)pcbddc; 1288da1bb401SStefano Zampini 1289da1bb401SStefano Zampini /* create PCIS data structure */ 1290da1bb401SStefano Zampini ierr = PCISCreate(pc);CHKERRQ(ierr); 1291da1bb401SStefano Zampini 129247d04d0dSStefano Zampini /* BDDC customization */ 129347d04d0dSStefano Zampini pcbddc->use_vertices = PETSC_TRUE; 129447d04d0dSStefano Zampini pcbddc->use_edges = PETSC_TRUE; 129547d04d0dSStefano Zampini pcbddc->use_faces = PETSC_FALSE; 129647d04d0dSStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 129747d04d0dSStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 129847d04d0dSStefano Zampini pcbddc->switch_static = PETSC_FALSE; 129947d04d0dSStefano Zampini pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 130047d04d0dSStefano Zampini pcbddc->dbg_flag = 0; 130147d04d0dSStefano Zampini 1302674ae819SStefano Zampini pcbddc->user_primal_vertices = 0; 13030bdf917eSStefano Zampini pcbddc->NullSpace = 0; 13043972b0daSStefano Zampini pcbddc->temp_solution = 0; 1305534831adSStefano Zampini pcbddc->original_rhs = 0; 1306534831adSStefano Zampini pcbddc->local_mat = 0; 1307534831adSStefano Zampini pcbddc->ChangeOfBasisMatrix = 0; 1308da1bb401SStefano Zampini pcbddc->coarse_vec = 0; 1309da1bb401SStefano Zampini pcbddc->coarse_rhs = 0; 1310da1bb401SStefano Zampini pcbddc->coarse_ksp = 0; 1311da1bb401SStefano Zampini pcbddc->coarse_phi_B = 0; 1312da1bb401SStefano Zampini pcbddc->coarse_phi_D = 0; 131315aaf578SStefano Zampini pcbddc->coarse_psi_B = 0; 131415aaf578SStefano Zampini pcbddc->coarse_psi_D = 0; 1315da1bb401SStefano Zampini pcbddc->vec1_P = 0; 1316da1bb401SStefano Zampini pcbddc->vec1_R = 0; 1317da1bb401SStefano Zampini pcbddc->vec2_R = 0; 1318da1bb401SStefano Zampini pcbddc->local_auxmat1 = 0; 1319da1bb401SStefano Zampini pcbddc->local_auxmat2 = 0; 1320da1bb401SStefano Zampini pcbddc->R_to_B = 0; 1321da1bb401SStefano Zampini pcbddc->R_to_D = 0; 1322da1bb401SStefano Zampini pcbddc->ksp_D = 0; 1323da1bb401SStefano Zampini pcbddc->ksp_R = 0; 1324da1bb401SStefano Zampini pcbddc->NeumannBoundaries = 0; 1325da1bb401SStefano Zampini pcbddc->ISForDofs = 0; 1326da1bb401SStefano Zampini pcbddc->ConstraintMatrix = 0; 132785c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 132847d04d0dSStefano Zampini pcbddc->coarse_loc_to_glob = 0; 132947d04d0dSStefano Zampini pcbddc->coarsening_ratio = 8; 13304fad6a16SStefano Zampini pcbddc->current_level = 0; 13312b510759SStefano Zampini pcbddc->max_levels = 0; 1332da1bb401SStefano Zampini 1333674ae819SStefano Zampini /* create local graph structure */ 1334674ae819SStefano Zampini ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1335674ae819SStefano Zampini 1336674ae819SStefano Zampini /* scaling */ 1337674ae819SStefano Zampini pcbddc->use_deluxe_scaling = PETSC_FALSE; 1338674ae819SStefano Zampini pcbddc->work_scaling = 0; 1339da1bb401SStefano Zampini 1340da1bb401SStefano Zampini /* function pointers */ 1341da1bb401SStefano Zampini pc->ops->apply = PCApply_BDDC; 1342da1bb401SStefano Zampini pc->ops->applytranspose = 0; 1343da1bb401SStefano Zampini pc->ops->setup = PCSetUp_BDDC; 1344da1bb401SStefano Zampini pc->ops->destroy = PCDestroy_BDDC; 1345da1bb401SStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1346da1bb401SStefano Zampini pc->ops->view = 0; 1347da1bb401SStefano Zampini pc->ops->applyrichardson = 0; 1348da1bb401SStefano Zampini pc->ops->applysymmetricleft = 0; 1349da1bb401SStefano Zampini pc->ops->applysymmetricright = 0; 1350534831adSStefano Zampini pc->ops->presolve = PCPreSolve_BDDC; 1351534831adSStefano Zampini pc->ops->postsolve = PCPostSolve_BDDC; 1352da1bb401SStefano Zampini 1353da1bb401SStefano Zampini /* composing function */ 1354674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1355bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 13562b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1357b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 13582b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1359bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1360bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1361bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1362bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1363bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1364bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 1365bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1366bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1367bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1368bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1369da1bb401SStefano Zampini PetscFunctionReturn(0); 1370da1bb401SStefano Zampini } 13713425bc38SStefano Zampini 1372