153cdbc3dSStefano Zampini /* TODOLIST 2*eb97c9d2SStefano Zampini 3*eb97c9d2SStefano Zampini ConstraintsSetup 4*eb97c9d2SStefano Zampini - assure same constraints between neighbours by sorting vals by global index before SVD! 5*eb97c9d2SStefano Zampini - tolerances for constraints as an option (take care of single precision!) 6*eb97c9d2SStefano Zampini - Allow different constraints customizations among different linear solves (requires also reset/destroy of ksp_R and coarse_ksp) 7*eb97c9d2SStefano Zampini - MAT_IGNORE_ZERO_ENTRIES for Constraints Matrix 8*eb97c9d2SStefano Zampini 9*eb97c9d2SStefano Zampini Solvers 10*eb97c9d2SStefano Zampini - Better management for block size > 1 for idx_R_local and others 11*eb97c9d2SStefano Zampini - Try to reduce the work when reusing the solvers 12*eb97c9d2SStefano Zampini - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers) 13*eb97c9d2SStefano Zampini - reuse already allocated coarse matrix if possible 14*eb97c9d2SStefano Zampini - Propagate ksp prefixes for solvers to mat objects? 15*eb97c9d2SStefano Zampini - Propagate nearnullspace info among levels 16*eb97c9d2SStefano Zampini 17*eb97c9d2SStefano Zampini User interface 18*eb97c9d2SStefano Zampini - Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet) 19*eb97c9d2SStefano Zampini - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access) 20*eb97c9d2SStefano Zampini - Provide PCApplyTranpose_BDDC 21*eb97c9d2SStefano Zampini - DofSplitting and DM attached to pc? 22*eb97c9d2SStefano Zampini 23*eb97c9d2SStefano Zampini Debugging output 24*eb97c9d2SStefano Zampini - Better management of verbosity levels of debugging output 25*eb97c9d2SStefano Zampini - Crashes on some architecture -> call SynchronizedAllow before every SynchronizedPrintf 26*eb97c9d2SStefano Zampini 27*eb97c9d2SStefano Zampini Build 28*eb97c9d2SStefano Zampini - make runexe59 29*eb97c9d2SStefano Zampini 30*eb97c9d2SStefano Zampini Extra 31*eb97c9d2SStefano Zampini - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)? 32*eb97c9d2SStefano Zampini - Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver? 33*eb97c9d2SStefano Zampini - add support for computing h,H and related using coordinates? 34*eb97c9d2SStefano Zampini - Change of basis approach does not work with my nonlinear mechanics example. why? maybe an issue with l2gmap? 35*eb97c9d2SStefano Zampini - Better management in PCIS code 36*eb97c9d2SStefano Zampini - BDDC with MG framework? 37*eb97c9d2SStefano Zampini 38*eb97c9d2SStefano Zampini FETIDP 39*eb97c9d2SStefano Zampini - Move FETIDP code to its own classes 40*eb97c9d2SStefano Zampini 41*eb97c9d2SStefano Zampini MATIS related operations contained in BDDC code 42*eb97c9d2SStefano Zampini - Provide general case for subassembling 43*eb97c9d2SStefano Zampini - Preallocation routines in MatConvert_IS_AIJ 44*eb97c9d2SStefano Zampini 4553cdbc3dSStefano Zampini */ 460c7d97c5SJed Brown 4753cdbc3dSStefano Zampini /* ---------------------------------------------------------------------------------------------------------------------------------------------- 480c7d97c5SJed Brown Implementation of BDDC preconditioner based on: 490c7d97c5SJed Brown C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 5053cdbc3dSStefano Zampini ---------------------------------------------------------------------------------------------------------------------------------------------- */ 5153cdbc3dSStefano Zampini 52674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 53674ae819SStefano Zampini #include "bddcprivate.h" 543b03a366Sstefano_zampini #include <petscblaslapack.h> 55674ae819SStefano Zampini 560c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 570c7d97c5SJed Brown #undef __FUNCT__ 580c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC" 590c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc) 600c7d97c5SJed Brown { 610c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 620c7d97c5SJed Brown PetscErrorCode ierr; 630c7d97c5SJed Brown 640c7d97c5SJed Brown PetscFunctionBegin; 650c7d97c5SJed Brown ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 668eeda7d8SStefano Zampini /* Verbose debugging */ 678eeda7d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 688eeda7d8SStefano Zampini /* Primal space cumstomization */ 698eeda7d8SStefano 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); 708eeda7d8SStefano 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); 718eeda7d8SStefano 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); 728eeda7d8SStefano Zampini /* Change of basis */ 738eeda7d8SStefano 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); 748eeda7d8SStefano 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); 75674ae819SStefano Zampini if (!pcbddc->use_change_of_basis) { 76674ae819SStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 77674ae819SStefano Zampini } 788eeda7d8SStefano Zampini /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 798eeda7d8SStefano 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); 800298fd71SBarry Smith ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 812b510759SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 82674ae819SStefano 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); 830c7d97c5SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 840c7d97c5SJed Brown PetscFunctionReturn(0); 850c7d97c5SJed Brown } 860c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 87674ae819SStefano Zampini #undef __FUNCT__ 88674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 89674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 90674ae819SStefano Zampini { 91674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 92674ae819SStefano Zampini PetscErrorCode ierr; 931e6b0712SBarry Smith 94674ae819SStefano Zampini PetscFunctionBegin; 95674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 96674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 97674ae819SStefano Zampini pcbddc->user_primal_vertices = PrimalVertices; 98674ae819SStefano Zampini PetscFunctionReturn(0); 99674ae819SStefano Zampini } 100674ae819SStefano Zampini #undef __FUNCT__ 101674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 102674ae819SStefano Zampini /*@ 10328509bceSStefano Zampini PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 104674ae819SStefano Zampini 105674ae819SStefano Zampini Not collective 106674ae819SStefano Zampini 107674ae819SStefano Zampini Input Parameters: 108674ae819SStefano Zampini + pc - the preconditioning context 10928509bceSStefano Zampini - PrimalVertices - index set of primal vertices in local numbering 110674ae819SStefano Zampini 111674ae819SStefano Zampini Level: intermediate 112674ae819SStefano Zampini 113674ae819SStefano Zampini Notes: 114674ae819SStefano Zampini 115674ae819SStefano Zampini .seealso: PCBDDC 116674ae819SStefano Zampini @*/ 117674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 118674ae819SStefano Zampini { 119674ae819SStefano Zampini PetscErrorCode ierr; 120674ae819SStefano Zampini 121674ae819SStefano Zampini PetscFunctionBegin; 122674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 123674ae819SStefano Zampini PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 124674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 125674ae819SStefano Zampini PetscFunctionReturn(0); 126674ae819SStefano Zampini } 127674ae819SStefano Zampini /* -------------------------------------------------------------------------- */ 1280c7d97c5SJed Brown #undef __FUNCT__ 1294fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 1304fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 1314fad6a16SStefano Zampini { 1324fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1334fad6a16SStefano Zampini 1344fad6a16SStefano Zampini PetscFunctionBegin; 1354fad6a16SStefano Zampini pcbddc->coarsening_ratio = k; 1364fad6a16SStefano Zampini PetscFunctionReturn(0); 1374fad6a16SStefano Zampini } 1381e6b0712SBarry Smith 1394fad6a16SStefano Zampini #undef __FUNCT__ 1404fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio" 1414fad6a16SStefano Zampini /*@ 14228509bceSStefano Zampini PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 1434fad6a16SStefano Zampini 1444fad6a16SStefano Zampini Logically collective on PC 1454fad6a16SStefano Zampini 1464fad6a16SStefano Zampini Input Parameters: 1474fad6a16SStefano Zampini + pc - the preconditioning context 14828509bceSStefano Zampini - k - coarsening ratio (H/h at the coarser level) 1494fad6a16SStefano Zampini 15028509bceSStefano Zampini Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 1514fad6a16SStefano Zampini 1524fad6a16SStefano Zampini Level: intermediate 1534fad6a16SStefano Zampini 1544fad6a16SStefano Zampini Notes: 1554fad6a16SStefano Zampini 1564fad6a16SStefano Zampini .seealso: PCBDDC 1574fad6a16SStefano Zampini @*/ 1584fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 1594fad6a16SStefano Zampini { 1604fad6a16SStefano Zampini PetscErrorCode ierr; 1614fad6a16SStefano Zampini 1624fad6a16SStefano Zampini PetscFunctionBegin; 1634fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1642b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,k,2); 1654fad6a16SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 1664fad6a16SStefano Zampini PetscFunctionReturn(0); 1674fad6a16SStefano Zampini } 1682b510759SStefano Zampini 169b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 1702b510759SStefano Zampini #undef __FUNCT__ 171b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 172b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 173b8ffe317SStefano Zampini { 174b8ffe317SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 175b8ffe317SStefano Zampini 176b8ffe317SStefano Zampini PetscFunctionBegin; 17785c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = flg; 178b8ffe317SStefano Zampini PetscFunctionReturn(0); 179b8ffe317SStefano Zampini } 180b8ffe317SStefano Zampini 181b8ffe317SStefano Zampini #undef __FUNCT__ 182b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 183b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 1842b510759SStefano Zampini { 1852b510759SStefano Zampini PetscErrorCode ierr; 1862b510759SStefano Zampini 1872b510759SStefano Zampini PetscFunctionBegin; 1882b510759SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 189b8ffe317SStefano Zampini PetscValidLogicalCollectiveBool(pc,flg,2); 190b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1912b510759SStefano Zampini PetscFunctionReturn(0); 1922b510759SStefano Zampini } 1931e6b0712SBarry Smith 1944fad6a16SStefano Zampini #undef __FUNCT__ 1952b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC" 1962b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 1974fad6a16SStefano Zampini { 1984fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1994fad6a16SStefano Zampini 2004fad6a16SStefano Zampini PetscFunctionBegin; 2012b510759SStefano Zampini pcbddc->current_level = level; 2024fad6a16SStefano Zampini PetscFunctionReturn(0); 2034fad6a16SStefano Zampini } 2041e6b0712SBarry Smith 2054fad6a16SStefano Zampini #undef __FUNCT__ 206b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel" 207b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 208b8ffe317SStefano Zampini { 209b8ffe317SStefano Zampini PetscErrorCode ierr; 210b8ffe317SStefano Zampini 211b8ffe317SStefano Zampini PetscFunctionBegin; 212b8ffe317SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 213b8ffe317SStefano Zampini PetscValidLogicalCollectiveInt(pc,level,2); 214b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 215b8ffe317SStefano Zampini PetscFunctionReturn(0); 216b8ffe317SStefano Zampini } 217b8ffe317SStefano Zampini 218b8ffe317SStefano Zampini #undef __FUNCT__ 2192b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC" 2202b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 2212b510759SStefano Zampini { 2222b510759SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2232b510759SStefano Zampini 2242b510759SStefano Zampini PetscFunctionBegin; 2252b510759SStefano Zampini pcbddc->max_levels = levels; 2262b510759SStefano Zampini PetscFunctionReturn(0); 2272b510759SStefano Zampini } 2282b510759SStefano Zampini 2292b510759SStefano Zampini #undef __FUNCT__ 2302b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels" 2314fad6a16SStefano Zampini /*@ 23228509bceSStefano Zampini PCBDDCSetLevels - Sets the maximum number of levels for multilevel 2334fad6a16SStefano Zampini 2344fad6a16SStefano Zampini Logically collective on PC 2354fad6a16SStefano Zampini 2364fad6a16SStefano Zampini Input Parameters: 2374fad6a16SStefano Zampini + pc - the preconditioning context 23828509bceSStefano Zampini - levels - the maximum number of levels (max 9) 2394fad6a16SStefano Zampini 24028509bceSStefano Zampini Default value is 0, i.e. traditional one-level BDDC 2414fad6a16SStefano Zampini 2424fad6a16SStefano Zampini Level: intermediate 2434fad6a16SStefano Zampini 2444fad6a16SStefano Zampini Notes: 2454fad6a16SStefano Zampini 2464fad6a16SStefano Zampini .seealso: PCBDDC 2474fad6a16SStefano Zampini @*/ 2482b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 2494fad6a16SStefano Zampini { 2504fad6a16SStefano Zampini PetscErrorCode ierr; 2514fad6a16SStefano Zampini 2524fad6a16SStefano Zampini PetscFunctionBegin; 2534fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2542b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,levels,2); 2552b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 2564fad6a16SStefano Zampini PetscFunctionReturn(0); 2574fad6a16SStefano Zampini } 2584fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */ 2591e6b0712SBarry Smith 2604fad6a16SStefano Zampini #undef __FUNCT__ 2610bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 2620bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 2630bdf917eSStefano Zampini { 2640bdf917eSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2650bdf917eSStefano Zampini PetscErrorCode ierr; 2660bdf917eSStefano Zampini 2670bdf917eSStefano Zampini PetscFunctionBegin; 2680bdf917eSStefano Zampini ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 2690bdf917eSStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 2700bdf917eSStefano Zampini pcbddc->NullSpace = NullSpace; 2710bdf917eSStefano Zampini PetscFunctionReturn(0); 2720bdf917eSStefano Zampini } 2731e6b0712SBarry Smith 2740bdf917eSStefano Zampini #undef __FUNCT__ 2750bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace" 2760bdf917eSStefano Zampini /*@ 27728509bceSStefano Zampini PCBDDCSetNullSpace - Set nullspace for BDDC operator 2780bdf917eSStefano Zampini 2790bdf917eSStefano Zampini Logically collective on PC and MatNullSpace 2800bdf917eSStefano Zampini 2810bdf917eSStefano Zampini Input Parameters: 2820bdf917eSStefano Zampini + pc - the preconditioning context 28328509bceSStefano Zampini - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 2840bdf917eSStefano Zampini 2850bdf917eSStefano Zampini Level: intermediate 2860bdf917eSStefano Zampini 2870bdf917eSStefano Zampini Notes: 2880bdf917eSStefano Zampini 2890bdf917eSStefano Zampini .seealso: PCBDDC 2900bdf917eSStefano Zampini @*/ 2910bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 2920bdf917eSStefano Zampini { 2930bdf917eSStefano Zampini PetscErrorCode ierr; 2940bdf917eSStefano Zampini 2950bdf917eSStefano Zampini PetscFunctionBegin; 2960bdf917eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 297674ae819SStefano Zampini PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 2982b510759SStefano Zampini PetscCheckSameComm(pc,1,NullSpace,2); 2990bdf917eSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 3000bdf917eSStefano Zampini PetscFunctionReturn(0); 3010bdf917eSStefano Zampini } 3020bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */ 3031e6b0712SBarry Smith 3040bdf917eSStefano Zampini #undef __FUNCT__ 3053b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 3063b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 3073b03a366Sstefano_zampini { 3083b03a366Sstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3093b03a366Sstefano_zampini PetscErrorCode ierr; 3103b03a366Sstefano_zampini 3113b03a366Sstefano_zampini PetscFunctionBegin; 3123b03a366Sstefano_zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 31336e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 31436e030ebSStefano Zampini pcbddc->DirichletBoundaries=DirichletBoundaries; 3153b03a366Sstefano_zampini PetscFunctionReturn(0); 3163b03a366Sstefano_zampini } 3171e6b0712SBarry Smith 3183b03a366Sstefano_zampini #undef __FUNCT__ 3193b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 3203b03a366Sstefano_zampini /*@ 32128509bceSStefano Zampini PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 3223b03a366Sstefano_zampini 3233b03a366Sstefano_zampini Not collective 3243b03a366Sstefano_zampini 3253b03a366Sstefano_zampini Input Parameters: 3263b03a366Sstefano_zampini + pc - the preconditioning context 32728509bceSStefano Zampini - DirichletBoundaries - sequential IS defining the subdomain part of Dirichlet boundaries (in local ordering) 3283b03a366Sstefano_zampini 3293b03a366Sstefano_zampini Level: intermediate 3303b03a366Sstefano_zampini 3313b03a366Sstefano_zampini Notes: 3323b03a366Sstefano_zampini 3333b03a366Sstefano_zampini .seealso: PCBDDC 3343b03a366Sstefano_zampini @*/ 3353b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 3363b03a366Sstefano_zampini { 3373b03a366Sstefano_zampini PetscErrorCode ierr; 3383b03a366Sstefano_zampini 3393b03a366Sstefano_zampini PetscFunctionBegin; 3403b03a366Sstefano_zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 341674ae819SStefano Zampini PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 3423b03a366Sstefano_zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 3433b03a366Sstefano_zampini PetscFunctionReturn(0); 3443b03a366Sstefano_zampini } 3453b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */ 3461e6b0712SBarry Smith 3473b03a366Sstefano_zampini #undef __FUNCT__ 3480c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 34953cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 3500c7d97c5SJed Brown { 3510c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 35253cdbc3dSStefano Zampini PetscErrorCode ierr; 3530c7d97c5SJed Brown 3540c7d97c5SJed Brown PetscFunctionBegin; 35553cdbc3dSStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 35636e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 35736e030ebSStefano Zampini pcbddc->NeumannBoundaries=NeumannBoundaries; 3580c7d97c5SJed Brown PetscFunctionReturn(0); 3590c7d97c5SJed Brown } 3601e6b0712SBarry Smith 3610c7d97c5SJed Brown #undef __FUNCT__ 3620c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 36357527edcSJed Brown /*@ 36428509bceSStefano Zampini PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 36557527edcSJed Brown 3669c0446d6SStefano Zampini Not collective 36757527edcSJed Brown 36857527edcSJed Brown Input Parameters: 36957527edcSJed Brown + pc - the preconditioning context 37028509bceSStefano Zampini - NeumannBoundaries - sequential IS defining the subdomain part of Neumann boundaries (in local ordering) 37157527edcSJed Brown 37257527edcSJed Brown Level: intermediate 37357527edcSJed Brown 37457527edcSJed Brown Notes: 37557527edcSJed Brown 37657527edcSJed Brown .seealso: PCBDDC 37757527edcSJed Brown @*/ 37853cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 3790c7d97c5SJed Brown { 3800c7d97c5SJed Brown PetscErrorCode ierr; 3810c7d97c5SJed Brown 3820c7d97c5SJed Brown PetscFunctionBegin; 3830c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 384674ae819SStefano Zampini PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 38553cdbc3dSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 38653cdbc3dSStefano Zampini PetscFunctionReturn(0); 38753cdbc3dSStefano Zampini } 38853cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 3891e6b0712SBarry Smith 39053cdbc3dSStefano Zampini #undef __FUNCT__ 391da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 392da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 393da1bb401SStefano Zampini { 394da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 395da1bb401SStefano Zampini 396da1bb401SStefano Zampini PetscFunctionBegin; 397da1bb401SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundaries; 398da1bb401SStefano Zampini PetscFunctionReturn(0); 399da1bb401SStefano Zampini } 4001e6b0712SBarry Smith 401da1bb401SStefano Zampini #undef __FUNCT__ 402da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 403da1bb401SStefano Zampini /*@ 40428509bceSStefano Zampini PCBDDCGetDirichletBoundaries - Get IS for local Dirichlet boundaries 405da1bb401SStefano Zampini 406da1bb401SStefano Zampini Not collective 407da1bb401SStefano Zampini 408da1bb401SStefano Zampini Input Parameters: 40928509bceSStefano Zampini . pc - the preconditioning context 410da1bb401SStefano Zampini 411da1bb401SStefano Zampini Output Parameters: 41228509bceSStefano Zampini . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 413da1bb401SStefano Zampini 414da1bb401SStefano Zampini Level: intermediate 415da1bb401SStefano Zampini 416da1bb401SStefano Zampini Notes: 417da1bb401SStefano Zampini 418da1bb401SStefano Zampini .seealso: PCBDDC 419da1bb401SStefano Zampini @*/ 420da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 421da1bb401SStefano Zampini { 422da1bb401SStefano Zampini PetscErrorCode ierr; 423da1bb401SStefano Zampini 424da1bb401SStefano Zampini PetscFunctionBegin; 425da1bb401SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 426da1bb401SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 427da1bb401SStefano Zampini PetscFunctionReturn(0); 428da1bb401SStefano Zampini } 429da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 4301e6b0712SBarry Smith 431da1bb401SStefano Zampini #undef __FUNCT__ 43253cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 43353cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 43453cdbc3dSStefano Zampini { 43553cdbc3dSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 43653cdbc3dSStefano Zampini 43753cdbc3dSStefano Zampini PetscFunctionBegin; 43853cdbc3dSStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundaries; 43953cdbc3dSStefano Zampini PetscFunctionReturn(0); 44053cdbc3dSStefano Zampini } 4411e6b0712SBarry Smith 44253cdbc3dSStefano Zampini #undef __FUNCT__ 44353cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 44453cdbc3dSStefano Zampini /*@ 44528509bceSStefano Zampini PCBDDCGetNeumannBoundaries - Get IS for local Neumann boundaries 44653cdbc3dSStefano Zampini 4479c0446d6SStefano Zampini Not collective 44853cdbc3dSStefano Zampini 44953cdbc3dSStefano Zampini Input Parameters: 45028509bceSStefano Zampini . pc - the preconditioning context 45153cdbc3dSStefano Zampini 45253cdbc3dSStefano Zampini Output Parameters: 45328509bceSStefano Zampini . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 45453cdbc3dSStefano Zampini 45553cdbc3dSStefano Zampini Level: intermediate 45653cdbc3dSStefano Zampini 45753cdbc3dSStefano Zampini Notes: 45853cdbc3dSStefano Zampini 45953cdbc3dSStefano Zampini .seealso: PCBDDC 46053cdbc3dSStefano Zampini @*/ 46153cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 46253cdbc3dSStefano Zampini { 46353cdbc3dSStefano Zampini PetscErrorCode ierr; 46453cdbc3dSStefano Zampini 46553cdbc3dSStefano Zampini PetscFunctionBegin; 46653cdbc3dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 46753cdbc3dSStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 4680c7d97c5SJed Brown PetscFunctionReturn(0); 4690c7d97c5SJed Brown } 47036e030ebSStefano Zampini /* -------------------------------------------------------------------------- */ 4711e6b0712SBarry Smith 47236e030ebSStefano Zampini #undef __FUNCT__ 473da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 4741a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 47536e030ebSStefano Zampini { 47636e030ebSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 477da1bb401SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 478da1bb401SStefano Zampini PetscErrorCode ierr; 47936e030ebSStefano Zampini 48036e030ebSStefano Zampini PetscFunctionBegin; 481674ae819SStefano Zampini /* free old CSR */ 482674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 483fb223d50SStefano Zampini /* TODO: PCBDDCGraphSetAdjacency */ 484674ae819SStefano Zampini /* get CSR into graph structure */ 485da1bb401SStefano Zampini if (copymode == PETSC_COPY_VALUES) { 486674ae819SStefano Zampini ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr); 487674ae819SStefano Zampini ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr); 488674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 489674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 490da1bb401SStefano Zampini } else if (copymode == PETSC_OWN_POINTER) { 4911a83f524SJed Brown mat_graph->xadj = (PetscInt*)xadj; 4921a83f524SJed Brown mat_graph->adjncy = (PetscInt*)adjncy; 493674ae819SStefano Zampini } 494575ad6abSStefano Zampini mat_graph->nvtxs_csr = nvtxs; 49536e030ebSStefano Zampini PetscFunctionReturn(0); 49636e030ebSStefano Zampini } 4971e6b0712SBarry Smith 49836e030ebSStefano Zampini #undef __FUNCT__ 499da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 50036e030ebSStefano Zampini /*@ 50128509bceSStefano Zampini PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix 50236e030ebSStefano Zampini 50336e030ebSStefano Zampini Not collective 50436e030ebSStefano Zampini 50536e030ebSStefano Zampini Input Parameters: 50636e030ebSStefano Zampini + pc - the preconditioning context 50728509bceSStefano Zampini . nvtxs - number of local vertices of the graph (i.e., the local size of your problem) 50828509bceSStefano Zampini . xadj, adjncy - the CSR graph 50928509bceSStefano Zampini - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 51036e030ebSStefano Zampini 51136e030ebSStefano Zampini Level: intermediate 51236e030ebSStefano Zampini 51336e030ebSStefano Zampini Notes: 51436e030ebSStefano Zampini 51528509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode 51636e030ebSStefano Zampini @*/ 5171a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 51836e030ebSStefano Zampini { 519575ad6abSStefano Zampini void (*f)(void) = 0; 52036e030ebSStefano Zampini PetscErrorCode ierr; 52136e030ebSStefano Zampini 52236e030ebSStefano Zampini PetscFunctionBegin; 52336e030ebSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 524674ae819SStefano Zampini PetscValidIntPointer(xadj,3); 525674ae819SStefano Zampini PetscValidIntPointer(xadj,4); 526674ae819SStefano Zampini if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 527674ae819SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 528da1bb401SStefano Zampini } 52936e030ebSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 530575ad6abSStefano Zampini /* free arrays if PCBDDC is not the PC type */ 531575ad6abSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 532575ad6abSStefano Zampini if (!f && copymode == PETSC_OWN_POINTER) { 533575ad6abSStefano Zampini ierr = PetscFree(xadj);CHKERRQ(ierr); 534575ad6abSStefano Zampini ierr = PetscFree(adjncy);CHKERRQ(ierr); 53536e030ebSStefano Zampini } 53636e030ebSStefano Zampini PetscFunctionReturn(0); 53736e030ebSStefano Zampini } 5389c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */ 5391e6b0712SBarry Smith 5409c0446d6SStefano Zampini #undef __FUNCT__ 5419c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 5429c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 5439c0446d6SStefano Zampini { 5449c0446d6SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 5459c0446d6SStefano Zampini PetscInt i; 5469c0446d6SStefano Zampini PetscErrorCode ierr; 5479c0446d6SStefano Zampini 5489c0446d6SStefano Zampini PetscFunctionBegin; 549da1bb401SStefano Zampini /* Destroy ISes if they were already set */ 5509c0446d6SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 5519c0446d6SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 5529c0446d6SStefano Zampini } 553d11ae9bbSstefano_zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 554da1bb401SStefano Zampini /* allocate space then set */ 5559c0446d6SStefano Zampini ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr); 5569c0446d6SStefano Zampini for (i=0;i<n_is;i++) { 557da1bb401SStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 558da1bb401SStefano Zampini pcbddc->ISForDofs[i]=ISForDofs[i]; 5599c0446d6SStefano Zampini } 5609c0446d6SStefano Zampini pcbddc->n_ISForDofs=n_is; 56160d44989SStefano Zampini pcbddc->user_provided_isfordofs = PETSC_TRUE; 5629c0446d6SStefano Zampini PetscFunctionReturn(0); 5639c0446d6SStefano Zampini } 5641e6b0712SBarry Smith 5659c0446d6SStefano Zampini #undef __FUNCT__ 5669c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting" 5679c0446d6SStefano Zampini /*@ 56828509bceSStefano Zampini PCBDDCSetDofsSplitting - Set index sets defining fields of the local Neumann matrix 5699c0446d6SStefano Zampini 5709c0446d6SStefano Zampini Not collective 5719c0446d6SStefano Zampini 5729c0446d6SStefano Zampini Input Parameters: 5739c0446d6SStefano Zampini + pc - the preconditioning context 57428509bceSStefano Zampini - n_is - number of index sets defining the fields 57528509bceSStefano Zampini . ISForDofs - array of IS describing the fields 5769c0446d6SStefano Zampini 5779c0446d6SStefano Zampini Level: intermediate 5789c0446d6SStefano Zampini 5799c0446d6SStefano Zampini Notes: 5809c0446d6SStefano Zampini 5819c0446d6SStefano Zampini .seealso: PCBDDC 5829c0446d6SStefano Zampini @*/ 5839c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 5849c0446d6SStefano Zampini { 5852b510759SStefano Zampini PetscInt i; 5869c0446d6SStefano Zampini PetscErrorCode ierr; 5879c0446d6SStefano Zampini 5889c0446d6SStefano Zampini PetscFunctionBegin; 5899c0446d6SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5902b510759SStefano Zampini for (i=0;i<n_is;i++) { 5912b510759SStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2); 5922b510759SStefano Zampini } 5939c0446d6SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 5949c0446d6SStefano Zampini PetscFunctionReturn(0); 5959c0446d6SStefano Zampini } 596da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 597534831adSStefano Zampini #undef __FUNCT__ 598534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC" 599534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 600534831adSStefano Zampini /* 601534831adSStefano Zampini PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 602534831adSStefano Zampini guess if a transformation of basis approach has been selected. 6039c0446d6SStefano Zampini 604534831adSStefano Zampini Input Parameter: 605534831adSStefano Zampini + pc - the preconditioner contex 606534831adSStefano Zampini 607534831adSStefano Zampini Application Interface Routine: PCPreSolve() 608534831adSStefano Zampini 609534831adSStefano Zampini Notes: 610534831adSStefano Zampini The interface routine PCPreSolve() is not usually called directly by 611534831adSStefano Zampini the user, but instead is called by KSPSolve(). 612534831adSStefano Zampini */ 613534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 614534831adSStefano Zampini { 615534831adSStefano Zampini PetscErrorCode ierr; 616534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 617534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 618534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 619534831adSStefano Zampini Mat temp_mat; 6203972b0daSStefano Zampini IS dirIS; 6213972b0daSStefano Zampini PetscInt dirsize,i,*is_indices; 6223972b0daSStefano Zampini PetscScalar *array_x,*array_diagonal; 6233972b0daSStefano Zampini Vec used_vec; 62492e3dcfbSStefano Zampini PetscBool guess_nonzero,flg,bddc_has_dirichlet_boundaries; 625534831adSStefano Zampini 626534831adSStefano Zampini PetscFunctionBegin; 62785c4d303SStefano Zampini /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 62885c4d303SStefano Zampini if (ksp) { 62985c4d303SStefano Zampini PetscBool iscg; 63085c4d303SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 63185c4d303SStefano Zampini if (!iscg) { 63285c4d303SStefano Zampini ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 63385c4d303SStefano Zampini } 63485c4d303SStefano Zampini } 63585c4d303SStefano Zampini /* Creates parallel work vectors used in presolve */ 63662a6ff1dSStefano Zampini if (!pcbddc->original_rhs) { 63762a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 63862a6ff1dSStefano Zampini } 63962a6ff1dSStefano Zampini if (!pcbddc->temp_solution) { 64062a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 64162a6ff1dSStefano Zampini } 6423972b0daSStefano Zampini if (x) { 6433972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 6443972b0daSStefano Zampini used_vec = x; 6453972b0daSStefano Zampini } else { 6463972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 6473972b0daSStefano Zampini used_vec = pcbddc->temp_solution; 6483972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 6493972b0daSStefano Zampini } 6503972b0daSStefano Zampini /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 6513972b0daSStefano Zampini if (ksp) { 6523972b0daSStefano Zampini ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 6533972b0daSStefano Zampini if (!guess_nonzero) { 6543972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 6553972b0daSStefano Zampini } 6563972b0daSStefano Zampini } 6573308cffdSStefano Zampini 65892e3dcfbSStefano Zampini /* TODO: remove when Dirichlet boundaries will be shared */ 65992e3dcfbSStefano Zampini ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr); 66092e3dcfbSStefano Zampini flg = PETSC_FALSE; 66192e3dcfbSStefano Zampini if (dirIS) flg = PETSC_TRUE; 66292e3dcfbSStefano Zampini ierr = MPI_Allreduce(&flg,&bddc_has_dirichlet_boundaries,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 66392e3dcfbSStefano Zampini 6643972b0daSStefano Zampini /* store the original rhs */ 6653972b0daSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 6663972b0daSStefano Zampini 6673972b0daSStefano Zampini /* Take into account zeroed rows -> change rhs and store solution removed */ 66892e3dcfbSStefano Zampini if (rhs && bddc_has_dirichlet_boundaries) { 6693972b0daSStefano Zampini ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 6703972b0daSStefano Zampini ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 6713972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6723972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6733972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6743972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6753972b0daSStefano Zampini if (dirIS) { 6763972b0daSStefano Zampini ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr); 6773972b0daSStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 6783972b0daSStefano Zampini ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 6793972b0daSStefano Zampini ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 6802fa5cd67SKarl Rupp for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 6813972b0daSStefano Zampini ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 6823972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 6833972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 6843972b0daSStefano Zampini } 6853972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6863972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 687b76ba322SStefano Zampini 6883972b0daSStefano Zampini /* remove the computed solution from the rhs */ 6893972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 6903972b0daSStefano Zampini ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 6913972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 6923308cffdSStefano Zampini } 693b76ba322SStefano Zampini 694b76ba322SStefano Zampini /* store partially computed solution and set initial guess */ 6953972b0daSStefano Zampini if (x) { 6963972b0daSStefano Zampini ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 6973972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 69885c4d303SStefano Zampini if (pcbddc->use_exact_dirichlet_trick) { 699b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 700b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 701b76ba322SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 702b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 703b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 704b76ba322SStefano Zampini if (ksp) { 705b76ba322SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 706b76ba322SStefano Zampini } 707b76ba322SStefano Zampini } 7083972b0daSStefano Zampini } 709b76ba322SStefano Zampini 71092e3dcfbSStefano Zampini /* prepare MatMult and rhs for solver */ 711674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 712b76ba322SStefano Zampini /* swap pointers for local matrices */ 713b76ba322SStefano Zampini temp_mat = matis->A; 714b76ba322SStefano Zampini matis->A = pcbddc->local_mat; 715b76ba322SStefano Zampini pcbddc->local_mat = temp_mat; 71692e3dcfbSStefano Zampini if (rhs) { 717b76ba322SStefano Zampini /* Get local rhs and apply transformation of basis */ 718b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 719b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 720b76ba322SStefano Zampini /* from original basis to modified basis */ 721b76ba322SStefano Zampini ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 722b76ba322SStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 723b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 724b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 725674ae819SStefano Zampini } 72692e3dcfbSStefano Zampini } 72792e3dcfbSStefano Zampini 72892e3dcfbSStefano Zampini /* remove nullspace if present */ 7290bdf917eSStefano Zampini if (ksp && pcbddc->NullSpace) { 730d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 731d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 732b76ba322SStefano Zampini } 7330bdf917eSStefano Zampini ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 734534831adSStefano Zampini PetscFunctionReturn(0); 735534831adSStefano Zampini } 736534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 737534831adSStefano Zampini #undef __FUNCT__ 738534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC" 739534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 740534831adSStefano Zampini /* 741534831adSStefano Zampini PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 742534831adSStefano Zampini approach has been selected. Also, restores rhs to its original state. 743534831adSStefano Zampini 744534831adSStefano Zampini Input Parameter: 745534831adSStefano Zampini + pc - the preconditioner contex 746534831adSStefano Zampini 747534831adSStefano Zampini Application Interface Routine: PCPostSolve() 748534831adSStefano Zampini 749534831adSStefano Zampini Notes: 750534831adSStefano Zampini The interface routine PCPostSolve() is not usually called directly by 751534831adSStefano Zampini the user, but instead is called by KSPSolve(). 752534831adSStefano Zampini */ 753534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 754534831adSStefano Zampini { 755534831adSStefano Zampini PetscErrorCode ierr; 756534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 757534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 758534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 759534831adSStefano Zampini Mat temp_mat; 760534831adSStefano Zampini 761534831adSStefano Zampini PetscFunctionBegin; 762674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 763534831adSStefano Zampini /* swap pointers for local matrices */ 764534831adSStefano Zampini temp_mat = matis->A; 765534831adSStefano Zampini matis->A = pcbddc->local_mat; 766534831adSStefano Zampini pcbddc->local_mat = temp_mat; 7673425bc38SStefano Zampini } 7683308cffdSStefano Zampini if (pcbddc->use_change_of_basis && x) { 769534831adSStefano Zampini /* Get Local boundary and apply transformation of basis to solution vector */ 770534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 771534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 772534831adSStefano Zampini /* from modified basis to original basis */ 773534831adSStefano Zampini ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 774534831adSStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 775534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 776534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 777534831adSStefano Zampini } 7783972b0daSStefano Zampini /* add solution removed in presolve */ 7793425bc38SStefano Zampini if (x) { 7803425bc38SStefano Zampini ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 7813425bc38SStefano Zampini } 782fb223d50SStefano Zampini /* restore rhs to its original state */ 783fb223d50SStefano Zampini if (rhs) { 784fb223d50SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 785fb223d50SStefano Zampini } 786534831adSStefano Zampini PetscFunctionReturn(0); 787534831adSStefano Zampini } 788534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 78953cdbc3dSStefano Zampini #undef __FUNCT__ 79053cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC" 7910c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 7920c7d97c5SJed Brown /* 7930c7d97c5SJed Brown PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 7940c7d97c5SJed Brown by setting data structures and options. 7950c7d97c5SJed Brown 7960c7d97c5SJed Brown Input Parameter: 79753cdbc3dSStefano Zampini + pc - the preconditioner context 7980c7d97c5SJed Brown 7990c7d97c5SJed Brown Application Interface Routine: PCSetUp() 8000c7d97c5SJed Brown 8010c7d97c5SJed Brown Notes: 8020c7d97c5SJed Brown The interface routine PCSetUp() is not usually called directly by 8030c7d97c5SJed Brown the user, but instead is called by PCApply() if necessary. 8040c7d97c5SJed Brown */ 80553cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc) 8060c7d97c5SJed Brown { 8070c7d97c5SJed Brown PetscErrorCode ierr; 8080c7d97c5SJed Brown PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 809674ae819SStefano Zampini MatStructure flag; 810674ae819SStefano Zampini PetscBool computeis,computetopography,computesolvers; 8110c7d97c5SJed Brown 8120c7d97c5SJed Brown PetscFunctionBegin; 813674ae819SStefano Zampini /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */ 814fcd409f5SStefano Zampini /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */ 8153b03a366Sstefano_zampini /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 816fcd409f5SStefano Zampini Also, BDDC directly build the Dirichlet problem */ 8173b03a366Sstefano_zampini /* Get stdout for dbg */ 818674ae819SStefano Zampini if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) { 819ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr); 820e269702eSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 8212b510759SStefano Zampini if (pcbddc->current_level) { 8222b510759SStefano Zampini ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 8232b510759SStefano Zampini } 824e269702eSStefano Zampini } 825674ae819SStefano Zampini /* first attempt to split work */ 826674ae819SStefano Zampini if (pc->setupcalled) { 827674ae819SStefano Zampini computeis = PETSC_FALSE; 828674ae819SStefano Zampini ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr); 829674ae819SStefano Zampini if (flag == SAME_PRECONDITIONER) { 830674ae819SStefano Zampini computetopography = PETSC_FALSE; 831674ae819SStefano Zampini computesolvers = PETSC_FALSE; 832674ae819SStefano Zampini } else if (flag == SAME_NONZERO_PATTERN) { 833674ae819SStefano Zampini computetopography = PETSC_FALSE; 834674ae819SStefano Zampini computesolvers = PETSC_TRUE; 835674ae819SStefano Zampini } else { /* DIFFERENT_NONZERO_PATTERN */ 836674ae819SStefano Zampini computetopography = PETSC_TRUE; 837674ae819SStefano Zampini computesolvers = PETSC_TRUE; 838674ae819SStefano Zampini } 839674ae819SStefano Zampini } else { 840674ae819SStefano Zampini computeis = PETSC_TRUE; 841674ae819SStefano Zampini computetopography = PETSC_TRUE; 842674ae819SStefano Zampini computesolvers = PETSC_TRUE; 843674ae819SStefano Zampini } 844fcd409f5SStefano Zampini /* Set up all the "iterative substructuring" common block without computing solvers */ 845674ae819SStefano Zampini if (computeis) { 846fcd409f5SStefano Zampini /* HACK INTO PCIS */ 847fcd409f5SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 848fcd409f5SStefano Zampini pcis->computesolvers = PETSC_FALSE; 849674ae819SStefano Zampini ierr = PCISSetUp(pc);CHKERRQ(ierr); 850674ae819SStefano Zampini } 851674ae819SStefano Zampini /* Analyze interface and set up local constraint and change of basis matrices */ 852674ae819SStefano Zampini if (computetopography) { 853674ae819SStefano Zampini /* reset data */ 854674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 855674ae819SStefano Zampini ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 856674ae819SStefano Zampini ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 857674ae819SStefano Zampini } 858674ae819SStefano Zampini if (computesolvers) { 859674ae819SStefano Zampini /* reset data */ 860674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 861674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 86299cc7994SStefano Zampini /* Create coarse and local stuffs */ 86399cc7994SStefano Zampini ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 864674ae819SStefano Zampini ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 8650c7d97c5SJed Brown } 8662b510759SStefano Zampini if (pcbddc->dbg_flag && pcbddc->current_level) { 8672b510759SStefano Zampini ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 8682b510759SStefano Zampini } 8690c7d97c5SJed Brown PetscFunctionReturn(0); 8700c7d97c5SJed Brown } 8710c7d97c5SJed Brown 8720c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 8730c7d97c5SJed Brown /* 8740c7d97c5SJed Brown PCApply_BDDC - Applies the BDDC preconditioner to a vector. 8750c7d97c5SJed Brown 8760c7d97c5SJed Brown Input Parameters: 8770c7d97c5SJed Brown . pc - the preconditioner context 8780c7d97c5SJed Brown . r - input vector (global) 8790c7d97c5SJed Brown 8800c7d97c5SJed Brown Output Parameter: 8810c7d97c5SJed Brown . z - output vector (global) 8820c7d97c5SJed Brown 8830c7d97c5SJed Brown Application Interface Routine: PCApply() 8840c7d97c5SJed Brown */ 8850c7d97c5SJed Brown #undef __FUNCT__ 8860c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC" 88753cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 8880c7d97c5SJed Brown { 8890c7d97c5SJed Brown PC_IS *pcis = (PC_IS*)(pc->data); 8900c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 8910c7d97c5SJed Brown PetscErrorCode ierr; 8923b03a366Sstefano_zampini const PetscScalar one = 1.0; 8933b03a366Sstefano_zampini const PetscScalar m_one = -1.0; 8942617d88aSStefano Zampini const PetscScalar zero = 0.0; 8950c7d97c5SJed Brown 8960c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN 8970c7d97c5SJed Brown NN interface preconditioner changed to BDDC 8988eeda7d8SStefano Zampini Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 8990c7d97c5SJed Brown 9000c7d97c5SJed Brown PetscFunctionBegin; 90185c4d303SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 9020c7d97c5SJed Brown /* First Dirichlet solve */ 9030c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9040c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 90553cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 9060c7d97c5SJed Brown /* 9070c7d97c5SJed Brown Assembling right hand side for BDDC operator 908674ae819SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 909674ae819SStefano Zampini - pcis->vec1_B the interface part of the global vector z 9100c7d97c5SJed Brown */ 9110c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 9120c7d97c5SJed Brown ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 9138eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 9140c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 9150c7d97c5SJed Brown ierr = VecCopy(r,z);CHKERRQ(ierr); 9160c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9170c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 918674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 919b76ba322SStefano Zampini } else { 9200bdf917eSStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 921b76ba322SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 922674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 923b76ba322SStefano Zampini } 924b76ba322SStefano Zampini 9252617d88aSStefano Zampini /* Apply interface preconditioner 9262617d88aSStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 9272617d88aSStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr); 9282617d88aSStefano Zampini 929674ae819SStefano Zampini /* Apply transpose of partition of unity operator */ 930674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 9310c7d97c5SJed Brown 9323b03a366Sstefano_zampini /* Second Dirichlet solve and assembling of output */ 9330c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9340c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9350c7d97c5SJed Brown ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 9368eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 93753cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr); 9380c7d97c5SJed Brown ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr); 9398eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 9400c7d97c5SJed Brown ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr); 9410c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9420c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9430c7d97c5SJed Brown PetscFunctionReturn(0); 9440c7d97c5SJed Brown } 945da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 946674ae819SStefano Zampini 947da1bb401SStefano Zampini #undef __FUNCT__ 948da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC" 949da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc) 950da1bb401SStefano Zampini { 951da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 952da1bb401SStefano Zampini PetscErrorCode ierr; 953da1bb401SStefano Zampini 954da1bb401SStefano Zampini PetscFunctionBegin; 955da1bb401SStefano Zampini /* free data created by PCIS */ 956da1bb401SStefano Zampini ierr = PCISDestroy(pc);CHKERRQ(ierr); 957674ae819SStefano Zampini /* free BDDC custom data */ 958674ae819SStefano Zampini ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 959674ae819SStefano Zampini /* destroy objects related to topography */ 960674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 961674ae819SStefano Zampini /* free allocated graph structure */ 962da1bb401SStefano Zampini ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 963674ae819SStefano Zampini /* free data for scaling operator */ 964674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 965674ae819SStefano Zampini /* free solvers stuff */ 966674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 96733bc96a4SStefano Zampini ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr); 96833bc96a4SStefano Zampini ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr); 9699881197aSStefano Zampini ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr); 970ac78edfcSStefano Zampini ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 97162a6ff1dSStefano Zampini /* free global vectors needed in presolve */ 97262a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 97362a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 9743425bc38SStefano Zampini /* remove functions */ 975674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 976bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 9772b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 978b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 9792b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 980bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 981bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 982bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 983bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 984bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 985bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 986bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 987bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 988bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 989bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 990674ae819SStefano Zampini /* Free the private data structure */ 991674ae819SStefano Zampini ierr = PetscFree(pc->data);CHKERRQ(ierr); 992da1bb401SStefano Zampini PetscFunctionReturn(0); 993da1bb401SStefano Zampini } 9943425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 9951e6b0712SBarry Smith 9963425bc38SStefano Zampini #undef __FUNCT__ 9973425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 9983425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 9993425bc38SStefano Zampini { 1000674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10013425bc38SStefano Zampini PC_IS* pcis; 10023425bc38SStefano Zampini PC_BDDC* pcbddc; 10033425bc38SStefano Zampini PetscErrorCode ierr; 10040c7d97c5SJed Brown 10053425bc38SStefano Zampini PetscFunctionBegin; 10063425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10073425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 10083425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 10093425bc38SStefano Zampini 10103425bc38SStefano Zampini /* change of basis for physical rhs if needed 10113425bc38SStefano Zampini It also changes the rhs in case of dirichlet boundaries */ 10123308cffdSStefano Zampini ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 10133425bc38SStefano Zampini /* store vectors for computation of fetidp final solution */ 10143425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10153425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1016fb223d50SStefano Zampini /* scale rhs since it should be unassembled */ 1017fb223d50SStefano Zampini /* TODO use counter scaling? (also below) */ 10183425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10193425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1020674ae819SStefano Zampini /* Apply partition of unity */ 10213425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1022674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 10238eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 10243425bc38SStefano Zampini /* compute partially subassembled Schur complement right-hand side */ 10253425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10263425bc38SStefano Zampini ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 10273425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 10283425bc38SStefano Zampini ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 10293425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10303425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1031674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 10323425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10333425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10343425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 10353425bc38SStefano Zampini } 10363425bc38SStefano Zampini /* BDDC rhs */ 10373425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 10388eeda7d8SStefano Zampini if (pcbddc->switch_static) { 10393425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10403425bc38SStefano Zampini } 10413425bc38SStefano Zampini /* apply BDDC */ 10423425bc38SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 10433425bc38SStefano Zampini /* Application of B_delta and assembling of rhs for fetidp fluxes */ 10443425bc38SStefano Zampini ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 10453425bc38SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 10463425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10473425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10483425bc38SStefano Zampini /* restore original rhs */ 10493425bc38SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 10503425bc38SStefano Zampini PetscFunctionReturn(0); 10513425bc38SStefano Zampini } 10521e6b0712SBarry Smith 10533425bc38SStefano Zampini #undef __FUNCT__ 10543425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 10553425bc38SStefano Zampini /*@ 105628509bceSStefano Zampini PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system 10573425bc38SStefano Zampini 10583425bc38SStefano Zampini Collective 10593425bc38SStefano Zampini 10603425bc38SStefano Zampini Input Parameters: 106128509bceSStefano Zampini + fetidp_mat - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators 106228509bceSStefano Zampini . standard_rhs - the right-hand side for your linear system 10633425bc38SStefano Zampini 10643425bc38SStefano Zampini Output Parameters: 106528509bceSStefano Zampini - fetidp_flux_rhs - the right-hand side for the FETIDP linear system 10663425bc38SStefano Zampini 10673425bc38SStefano Zampini Level: developer 10683425bc38SStefano Zampini 10693425bc38SStefano Zampini Notes: 10703425bc38SStefano Zampini 107128509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 10723425bc38SStefano Zampini @*/ 10733425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 10743425bc38SStefano Zampini { 1075674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10763425bc38SStefano Zampini PetscErrorCode ierr; 10773425bc38SStefano Zampini 10783425bc38SStefano Zampini PetscFunctionBegin; 10793425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10803425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 10813425bc38SStefano Zampini PetscFunctionReturn(0); 10823425bc38SStefano Zampini } 10833425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 10841e6b0712SBarry Smith 10853425bc38SStefano Zampini #undef __FUNCT__ 10863425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 10873425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 10883425bc38SStefano Zampini { 1089674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10903425bc38SStefano Zampini PC_IS* pcis; 10913425bc38SStefano Zampini PC_BDDC* pcbddc; 10923425bc38SStefano Zampini PetscErrorCode ierr; 10933425bc38SStefano Zampini 10943425bc38SStefano Zampini PetscFunctionBegin; 10953425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10963425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 10973425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 10983425bc38SStefano Zampini 10993425bc38SStefano Zampini /* apply B_delta^T */ 11003425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11013425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11023425bc38SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 11033425bc38SStefano Zampini /* compute rhs for BDDC application */ 11043425bc38SStefano Zampini ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 11058eeda7d8SStefano Zampini if (pcbddc->switch_static) { 11063425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 11073425bc38SStefano Zampini } 11083425bc38SStefano Zampini /* apply BDDC */ 11093425bc38SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 11103425bc38SStefano Zampini /* put values into standard global vector */ 11113425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11123425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11138eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 11143425bc38SStefano Zampini /* compute values into the interior if solved for the partially subassembled Schur complement */ 11153425bc38SStefano Zampini ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 11163425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 11173425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 11183425bc38SStefano Zampini } 11193425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11203425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11213425bc38SStefano Zampini /* final change of basis if needed 11223425bc38SStefano Zampini Is also sums the dirichlet part removed during RHS assembling */ 11233308cffdSStefano Zampini ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 11243425bc38SStefano Zampini PetscFunctionReturn(0); 11253425bc38SStefano Zampini } 11261e6b0712SBarry Smith 11273425bc38SStefano Zampini #undef __FUNCT__ 11283425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 11293425bc38SStefano Zampini /*@ 113028509bceSStefano Zampini PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system 11313425bc38SStefano Zampini 11323425bc38SStefano Zampini Collective 11333425bc38SStefano Zampini 11343425bc38SStefano Zampini Input Parameters: 113528509bceSStefano Zampini + fetidp_mat - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators 113628509bceSStefano Zampini . fetidp_flux_sol - the solution of the FETIDP linear system 11373425bc38SStefano Zampini 11383425bc38SStefano Zampini Output Parameters: 113928509bceSStefano Zampini - standard_sol - the solution defined on the physical domain 11403425bc38SStefano Zampini 11413425bc38SStefano Zampini Level: developer 11423425bc38SStefano Zampini 11433425bc38SStefano Zampini Notes: 11443425bc38SStefano Zampini 114528509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 11463425bc38SStefano Zampini @*/ 11473425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 11483425bc38SStefano Zampini { 1149674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 11503425bc38SStefano Zampini PetscErrorCode ierr; 11513425bc38SStefano Zampini 11523425bc38SStefano Zampini PetscFunctionBegin; 11533425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 11543425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 11553425bc38SStefano Zampini PetscFunctionReturn(0); 11563425bc38SStefano Zampini } 11573425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 11581e6b0712SBarry Smith 1159f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1160f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1161f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1162f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1163674ae819SStefano Zampini 11643425bc38SStefano Zampini #undef __FUNCT__ 11653425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 11663425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 11673425bc38SStefano Zampini { 1168674ae819SStefano Zampini 1169674ae819SStefano Zampini FETIDPMat_ctx fetidpmat_ctx; 11703425bc38SStefano Zampini Mat newmat; 1171674ae819SStefano Zampini FETIDPPC_ctx fetidppc_ctx; 11723425bc38SStefano Zampini PC newpc; 1173ce94432eSBarry Smith MPI_Comm comm; 11743425bc38SStefano Zampini PetscErrorCode ierr; 11753425bc38SStefano Zampini 11763425bc38SStefano Zampini PetscFunctionBegin; 1177ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 11783425bc38SStefano Zampini /* FETIDP linear matrix */ 11793425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 11803425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 11813425bc38SStefano Zampini ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 11823425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 11833425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 11843425bc38SStefano Zampini ierr = MatSetUp(newmat);CHKERRQ(ierr); 11853425bc38SStefano Zampini /* FETIDP preconditioner */ 11863425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 11873425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 11883425bc38SStefano Zampini ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 11893425bc38SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 11903425bc38SStefano Zampini ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 11913425bc38SStefano Zampini ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 11923425bc38SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 11933425bc38SStefano Zampini ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr); 11943425bc38SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 11953425bc38SStefano Zampini /* return pointers for objects created */ 11963425bc38SStefano Zampini *fetidp_mat=newmat; 11973425bc38SStefano Zampini *fetidp_pc=newpc; 11983425bc38SStefano Zampini PetscFunctionReturn(0); 11993425bc38SStefano Zampini } 12001e6b0712SBarry Smith 12013425bc38SStefano Zampini #undef __FUNCT__ 12023425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 12033425bc38SStefano Zampini /*@ 120428509bceSStefano Zampini PCBDDCCreateFETIDPOperators - Create operators for FETIDP 12053425bc38SStefano Zampini 12063425bc38SStefano Zampini Collective 12073425bc38SStefano Zampini 12083425bc38SStefano Zampini Input Parameters: 120928509bceSStefano Zampini + pc - the BDDC preconditioning context already setup 121028509bceSStefano Zampini 121128509bceSStefano Zampini Output Parameters: 121228509bceSStefano Zampini . fetidp_mat - shell FETIDP matrix object 121328509bceSStefano Zampini . fetidp_pc - shell Dirichlet preconditioner for FETIDP matrix 121428509bceSStefano Zampini 121528509bceSStefano Zampini Options Database Keys: 121628509bceSStefano Zampini - -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers 12173425bc38SStefano Zampini 12183425bc38SStefano Zampini Level: developer 12193425bc38SStefano Zampini 12203425bc38SStefano Zampini Notes: 122128509bceSStefano Zampini Currently the only operation provided for FETIDP matrix is MatMult 12223425bc38SStefano Zampini 12233425bc38SStefano Zampini .seealso: PCBDDC 12243425bc38SStefano Zampini @*/ 12253425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 12263425bc38SStefano Zampini { 12273425bc38SStefano Zampini PetscErrorCode ierr; 12283425bc38SStefano Zampini 12293425bc38SStefano Zampini PetscFunctionBegin; 12303425bc38SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12313425bc38SStefano Zampini if (pc->setupcalled) { 1232516d51a7SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1233f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 12343425bc38SStefano Zampini PetscFunctionReturn(0); 12353425bc38SStefano Zampini } 12360c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 1237da1bb401SStefano Zampini /*MC 1238da1bb401SStefano Zampini PCBDDC - Balancing Domain Decomposition by Constraints. 12390c7d97c5SJed Brown 124028509bceSStefano Zampini An implementation of the BDDC preconditioner based on 124128509bceSStefano Zampini 124228509bceSStefano Zampini .vb 124328509bceSStefano Zampini [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 124428509bceSStefano 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 124528509bceSStefano Zampini [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 124628509bceSStefano Zampini .ve 124728509bceSStefano Zampini 124828509bceSStefano Zampini The matrix to be preconditioned (Pmat) must be of type MATIS. 124928509bceSStefano Zampini 1250b6fdb6dfSStefano Zampini Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 125128509bceSStefano Zampini 125228509bceSStefano Zampini It also works with unsymmetric and indefinite problems. 125328509bceSStefano Zampini 1254b6fdb6dfSStefano 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. 1255b6fdb6dfSStefano Zampini 125628509bceSStefano Zampini Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace 125728509bceSStefano Zampini 125828509bceSStefano 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 125928509bceSStefano Zampini 126028509bceSStefano Zampini Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace. 126128509bceSStefano Zampini 1262b6fdb6dfSStefano Zampini Change of basis is performed similarly to [2] when requested. When more the one constraint is present on a single connected component (i.e. an edge or a face), a robust method based on local QR factorizations is used. 126328509bceSStefano Zampini 126428509bceSStefano Zampini The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object. 126528509bceSStefano Zampini 1266da1bb401SStefano Zampini Options Database Keys: 126728509bceSStefano Zampini 126828509bceSStefano Zampini . -pc_bddc_use_vertices <1> - use or not vertices in primal space 126928509bceSStefano Zampini . -pc_bddc_use_edges <1> - use or not edges in primal space 1270b6fdb6dfSStefano Zampini . -pc_bddc_use_faces <0> - use or not faces in primal space 127128509bceSStefano Zampini . -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only) 127228509bceSStefano Zampini . -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested 127328509bceSStefano Zampini . -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1]) 127428509bceSStefano Zampini . -pc_bddc_levels <0> - maximum number of levels for multilevel 127528509bceSStefano Zampini . -pc_bddc_coarsening_ratio - H/h ratio at the coarser level 127628509bceSStefano Zampini - -pc_bddc_check_level <0> - set verbosity level of debugging output 127728509bceSStefano Zampini 127828509bceSStefano Zampini Options for Dirichlet, Neumann or coarse solver can be set with 127928509bceSStefano Zampini .vb 128028509bceSStefano Zampini -pc_bddc_dirichlet_ 128128509bceSStefano Zampini -pc_bddc_neumann_ 128228509bceSStefano Zampini -pc_bddc_coarse_ 128328509bceSStefano Zampini .ve 128428509bceSStefano Zampini e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg 128528509bceSStefano Zampini 128628509bceSStefano Zampini When using a multilevel approach, solvers' options at the N-th level can be specified as 128728509bceSStefano Zampini .vb 128828509bceSStefano Zampini -pc_bddc_dirichlet_N_ 128928509bceSStefano Zampini -pc_bddc_neumann_N_ 129028509bceSStefano Zampini -pc_bddc_coarse_N_ 129128509bceSStefano Zampini .ve 129228509bceSStefano Zampini Note that level number ranges from the finest 0 to the coarsest N 1293da1bb401SStefano Zampini 1294da1bb401SStefano Zampini Level: intermediate 1295da1bb401SStefano Zampini 1296b6fdb6dfSStefano Zampini Developer notes: 129728509bceSStefano Zampini Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose 1298da1bb401SStefano Zampini 129928509bceSStefano Zampini New deluxe scaling operator will be available soon. 1300da1bb401SStefano Zampini 1301da1bb401SStefano Zampini Contributed by Stefano Zampini 1302da1bb401SStefano Zampini 1303da1bb401SStefano Zampini .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1304da1bb401SStefano Zampini M*/ 1305b2573a8aSBarry Smith 1306da1bb401SStefano Zampini #undef __FUNCT__ 1307da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC" 13088cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1309da1bb401SStefano Zampini { 1310da1bb401SStefano Zampini PetscErrorCode ierr; 1311da1bb401SStefano Zampini PC_BDDC *pcbddc; 1312da1bb401SStefano Zampini 1313da1bb401SStefano Zampini PetscFunctionBegin; 1314da1bb401SStefano Zampini /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1315da1bb401SStefano Zampini ierr = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr); 1316da1bb401SStefano Zampini pc->data = (void*)pcbddc; 1317da1bb401SStefano Zampini 1318da1bb401SStefano Zampini /* create PCIS data structure */ 1319da1bb401SStefano Zampini ierr = PCISCreate(pc);CHKERRQ(ierr); 1320da1bb401SStefano Zampini 132147d04d0dSStefano Zampini /* BDDC customization */ 132247d04d0dSStefano Zampini pcbddc->use_vertices = PETSC_TRUE; 132347d04d0dSStefano Zampini pcbddc->use_edges = PETSC_TRUE; 132447d04d0dSStefano Zampini pcbddc->use_faces = PETSC_FALSE; 132547d04d0dSStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 132647d04d0dSStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 132747d04d0dSStefano Zampini pcbddc->switch_static = PETSC_FALSE; 132847d04d0dSStefano Zampini pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 132947d04d0dSStefano Zampini pcbddc->dbg_flag = 0; 133047d04d0dSStefano Zampini 1331674ae819SStefano Zampini pcbddc->user_primal_vertices = 0; 13320bdf917eSStefano Zampini pcbddc->NullSpace = 0; 13333972b0daSStefano Zampini pcbddc->temp_solution = 0; 1334534831adSStefano Zampini pcbddc->original_rhs = 0; 1335534831adSStefano Zampini pcbddc->local_mat = 0; 1336534831adSStefano Zampini pcbddc->ChangeOfBasisMatrix = 0; 1337da1bb401SStefano Zampini pcbddc->coarse_vec = 0; 1338da1bb401SStefano Zampini pcbddc->coarse_rhs = 0; 1339da1bb401SStefano Zampini pcbddc->coarse_ksp = 0; 1340da1bb401SStefano Zampini pcbddc->coarse_phi_B = 0; 1341da1bb401SStefano Zampini pcbddc->coarse_phi_D = 0; 134215aaf578SStefano Zampini pcbddc->coarse_psi_B = 0; 134315aaf578SStefano Zampini pcbddc->coarse_psi_D = 0; 1344da1bb401SStefano Zampini pcbddc->vec1_P = 0; 1345da1bb401SStefano Zampini pcbddc->vec1_R = 0; 1346da1bb401SStefano Zampini pcbddc->vec2_R = 0; 1347da1bb401SStefano Zampini pcbddc->local_auxmat1 = 0; 1348da1bb401SStefano Zampini pcbddc->local_auxmat2 = 0; 1349da1bb401SStefano Zampini pcbddc->R_to_B = 0; 1350da1bb401SStefano Zampini pcbddc->R_to_D = 0; 1351da1bb401SStefano Zampini pcbddc->ksp_D = 0; 1352da1bb401SStefano Zampini pcbddc->ksp_R = 0; 1353da1bb401SStefano Zampini pcbddc->NeumannBoundaries = 0; 135460d44989SStefano Zampini pcbddc->user_provided_isfordofs = PETSC_FALSE; 135560d44989SStefano Zampini pcbddc->n_ISForDofs = 0; 1356da1bb401SStefano Zampini pcbddc->ISForDofs = 0; 1357da1bb401SStefano Zampini pcbddc->ConstraintMatrix = 0; 135885c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 135947d04d0dSStefano Zampini pcbddc->coarse_loc_to_glob = 0; 136047d04d0dSStefano Zampini pcbddc->coarsening_ratio = 8; 13614fad6a16SStefano Zampini pcbddc->current_level = 0; 13622b510759SStefano Zampini pcbddc->max_levels = 0; 1363da1bb401SStefano Zampini 1364674ae819SStefano Zampini /* create local graph structure */ 1365674ae819SStefano Zampini ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1366674ae819SStefano Zampini 1367674ae819SStefano Zampini /* scaling */ 1368674ae819SStefano Zampini pcbddc->use_deluxe_scaling = PETSC_FALSE; 1369674ae819SStefano Zampini pcbddc->work_scaling = 0; 1370da1bb401SStefano Zampini 1371da1bb401SStefano Zampini /* function pointers */ 1372da1bb401SStefano Zampini pc->ops->apply = PCApply_BDDC; 1373da1bb401SStefano Zampini pc->ops->applytranspose = 0; 1374da1bb401SStefano Zampini pc->ops->setup = PCSetUp_BDDC; 1375da1bb401SStefano Zampini pc->ops->destroy = PCDestroy_BDDC; 1376da1bb401SStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1377da1bb401SStefano Zampini pc->ops->view = 0; 1378da1bb401SStefano Zampini pc->ops->applyrichardson = 0; 1379da1bb401SStefano Zampini pc->ops->applysymmetricleft = 0; 1380da1bb401SStefano Zampini pc->ops->applysymmetricright = 0; 1381534831adSStefano Zampini pc->ops->presolve = PCPreSolve_BDDC; 1382534831adSStefano Zampini pc->ops->postsolve = PCPostSolve_BDDC; 1383da1bb401SStefano Zampini 1384da1bb401SStefano Zampini /* composing function */ 1385674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1386bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 13872b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1388b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 13892b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1390bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1391bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1392bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1393bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1394bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1395bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 1396bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1397bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1398bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1399bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1400da1bb401SStefano Zampini PetscFunctionReturn(0); 1401da1bb401SStefano Zampini } 14023425bc38SStefano Zampini 1403