153cdbc3dSStefano Zampini /* TODOLIST 2eb97c9d2SStefano Zampini 3eb97c9d2SStefano Zampini ConstraintsSetup 4eb97c9d2SStefano Zampini - tolerances for constraints as an option (take care of single precision!) 5*b9b85e73SStefano Zampini - * add support for a user which specifies a change of basis 6eb97c9d2SStefano Zampini 7eb97c9d2SStefano Zampini Solvers 8eb97c9d2SStefano Zampini - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers) 9eb97c9d2SStefano Zampini - Propagate ksp prefixes for solvers to mat objects? 10eb97c9d2SStefano Zampini - Propagate nearnullspace info among levels 11eb97c9d2SStefano Zampini 12eb97c9d2SStefano Zampini User interface 13*b9b85e73SStefano Zampini - * Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access) 14*b9b85e73SStefano Zampini - ** DofSplitting and DM attached to pc? 15eb97c9d2SStefano Zampini 16eb97c9d2SStefano Zampini Debugging output 17*b9b85e73SStefano Zampini - * Better management of verbosity levels of debugging output 18eb97c9d2SStefano Zampini 19eb97c9d2SStefano Zampini Build 20eb97c9d2SStefano Zampini - make runexe59 21eb97c9d2SStefano Zampini 22eb97c9d2SStefano Zampini Extra 23*b9b85e73SStefano Zampini - ** GetRid of PCBDDCApplySchur, use MatSchur instead 24*b9b85e73SStefano Zampini - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)? 25eb97c9d2SStefano Zampini - add support for computing h,H and related using coordinates? 26c0b83709SStefano Zampini - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap) 27eb97c9d2SStefano Zampini - Better management in PCIS code 28eb97c9d2SStefano Zampini - BDDC with MG framework? 29eb97c9d2SStefano Zampini 30eb97c9d2SStefano Zampini FETIDP 31eb97c9d2SStefano Zampini - Move FETIDP code to its own classes 32eb97c9d2SStefano Zampini 33eb97c9d2SStefano Zampini MATIS related operations contained in BDDC code 34eb97c9d2SStefano Zampini - Provide general case for subassembling 35*b9b85e73SStefano Zampini - *** Preallocation routines in MatISGetMPIAXAIJ 36eb97c9d2SStefano Zampini 3753cdbc3dSStefano Zampini */ 380c7d97c5SJed Brown 3953cdbc3dSStefano Zampini /* ---------------------------------------------------------------------------------------------------------------------------------------------- 400c7d97c5SJed Brown Implementation of BDDC preconditioner based on: 410c7d97c5SJed Brown C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 4253cdbc3dSStefano Zampini ---------------------------------------------------------------------------------------------------------------------------------------------- */ 4353cdbc3dSStefano Zampini 44ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 45ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 463b03a366Sstefano_zampini #include <petscblaslapack.h> 47674ae819SStefano Zampini 480c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 490c7d97c5SJed Brown #undef __FUNCT__ 500c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC" 510c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc) 520c7d97c5SJed Brown { 530c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 540c7d97c5SJed Brown PetscErrorCode ierr; 550c7d97c5SJed Brown 560c7d97c5SJed Brown PetscFunctionBegin; 570c7d97c5SJed Brown ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 588eeda7d8SStefano Zampini /* Verbose debugging */ 598eeda7d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 608eeda7d8SStefano Zampini /* Primal space cumstomization */ 6108a5cf49SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);CHKERRQ(ierr); 628eeda7d8SStefano 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); 638eeda7d8SStefano 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); 648eeda7d8SStefano 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); 658eeda7d8SStefano Zampini /* Change of basis */ 66*b9b85e73SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not internal change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr); 67*b9b85e73SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not internal change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr); 68674ae819SStefano Zampini if (!pcbddc->use_change_of_basis) { 69674ae819SStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 70674ae819SStefano Zampini } 718eeda7d8SStefano Zampini /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 728eeda7d8SStefano 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); 730298fd71SBarry Smith ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 742b510759SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 75323d395dSStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);CHKERRQ(ierr); 76674ae819SStefano 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); 770c7d97c5SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 780c7d97c5SJed Brown PetscFunctionReturn(0); 790c7d97c5SJed Brown } 800c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 81674ae819SStefano Zampini #undef __FUNCT__ 82*b9b85e73SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisLocalMat_BDDC" 83*b9b85e73SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisLocalMat_BDDC(PC pc, Mat change) 84*b9b85e73SStefano Zampini { 85*b9b85e73SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 86*b9b85e73SStefano Zampini PetscErrorCode ierr; 87*b9b85e73SStefano Zampini 88*b9b85e73SStefano Zampini PetscFunctionBegin; 89*b9b85e73SStefano Zampini ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr); 90*b9b85e73SStefano Zampini ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr); 91*b9b85e73SStefano Zampini pcbddc->user_ChangeOfBasisMatrix = change; 92*b9b85e73SStefano Zampini PetscFunctionReturn(0); 93*b9b85e73SStefano Zampini } 94*b9b85e73SStefano Zampini #undef __FUNCT__ 95*b9b85e73SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisLocalMat" 96*b9b85e73SStefano Zampini /*@ 97*b9b85e73SStefano Zampini PCBDDCSetChangeOfBasisLocalMat - Set user defined change of basis for local boundary dofs 98*b9b85e73SStefano Zampini 99*b9b85e73SStefano Zampini Collective on PC 100*b9b85e73SStefano Zampini 101*b9b85e73SStefano Zampini Input Parameters: 102*b9b85e73SStefano Zampini + pc - the preconditioning context 103*b9b85e73SStefano Zampini - change - the local change of basis matrix, either in local (internal + boundary) or in local boundary numbering 104*b9b85e73SStefano Zampini 105*b9b85e73SStefano Zampini Level: intermediate 106*b9b85e73SStefano Zampini 107*b9b85e73SStefano Zampini Notes: 108*b9b85e73SStefano Zampini 109*b9b85e73SStefano Zampini .seealso: PCBDDC 110*b9b85e73SStefano Zampini @*/ 111*b9b85e73SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisLocalMat(PC pc, Mat change) 112*b9b85e73SStefano Zampini { 113*b9b85e73SStefano Zampini PetscErrorCode ierr; 114*b9b85e73SStefano Zampini 115*b9b85e73SStefano Zampini PetscFunctionBegin; 116*b9b85e73SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 117*b9b85e73SStefano Zampini PetscValidHeaderSpecific(change,MAT_CLASSID,2); 118*b9b85e73SStefano Zampini PetscValidLogicalCollectiveBool(pc,change,2); 119*b9b85e73SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisLocalMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr); 120*b9b85e73SStefano Zampini PetscFunctionReturn(0); 121*b9b85e73SStefano Zampini } 122*b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */ 123*b9b85e73SStefano Zampini #undef __FUNCT__ 124674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 125674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 126674ae819SStefano Zampini { 127674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 128674ae819SStefano Zampini PetscErrorCode ierr; 1291e6b0712SBarry Smith 130674ae819SStefano Zampini PetscFunctionBegin; 131674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 132674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 133674ae819SStefano Zampini pcbddc->user_primal_vertices = PrimalVertices; 134674ae819SStefano Zampini PetscFunctionReturn(0); 135674ae819SStefano Zampini } 136674ae819SStefano Zampini #undef __FUNCT__ 137674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 138674ae819SStefano Zampini /*@ 13928509bceSStefano Zampini PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 140674ae819SStefano Zampini 141674ae819SStefano Zampini Not collective 142674ae819SStefano Zampini 143674ae819SStefano Zampini Input Parameters: 144674ae819SStefano Zampini + pc - the preconditioning context 14528509bceSStefano Zampini - PrimalVertices - index set of primal vertices in local numbering 146674ae819SStefano Zampini 147674ae819SStefano Zampini Level: intermediate 148674ae819SStefano Zampini 149674ae819SStefano Zampini Notes: 150674ae819SStefano Zampini 151674ae819SStefano Zampini .seealso: PCBDDC 152674ae819SStefano Zampini @*/ 153674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 154674ae819SStefano Zampini { 155674ae819SStefano Zampini PetscErrorCode ierr; 156674ae819SStefano Zampini 157674ae819SStefano Zampini PetscFunctionBegin; 158674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 159674ae819SStefano Zampini PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 160674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 161674ae819SStefano Zampini PetscFunctionReturn(0); 162674ae819SStefano Zampini } 163674ae819SStefano Zampini /* -------------------------------------------------------------------------- */ 1640c7d97c5SJed Brown #undef __FUNCT__ 1654fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 1664fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 1674fad6a16SStefano Zampini { 1684fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1694fad6a16SStefano Zampini 1704fad6a16SStefano Zampini PetscFunctionBegin; 1714fad6a16SStefano Zampini pcbddc->coarsening_ratio = k; 1724fad6a16SStefano Zampini PetscFunctionReturn(0); 1734fad6a16SStefano Zampini } 1741e6b0712SBarry Smith 1754fad6a16SStefano Zampini #undef __FUNCT__ 1764fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio" 1774fad6a16SStefano Zampini /*@ 17828509bceSStefano Zampini PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 1794fad6a16SStefano Zampini 1804fad6a16SStefano Zampini Logically collective on PC 1814fad6a16SStefano Zampini 1824fad6a16SStefano Zampini Input Parameters: 1834fad6a16SStefano Zampini + pc - the preconditioning context 18428509bceSStefano Zampini - k - coarsening ratio (H/h at the coarser level) 1854fad6a16SStefano Zampini 18628509bceSStefano Zampini Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 1874fad6a16SStefano Zampini 1884fad6a16SStefano Zampini Level: intermediate 1894fad6a16SStefano Zampini 1904fad6a16SStefano Zampini Notes: 1914fad6a16SStefano Zampini 1924fad6a16SStefano Zampini .seealso: PCBDDC 1934fad6a16SStefano Zampini @*/ 1944fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 1954fad6a16SStefano Zampini { 1964fad6a16SStefano Zampini PetscErrorCode ierr; 1974fad6a16SStefano Zampini 1984fad6a16SStefano Zampini PetscFunctionBegin; 1994fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2002b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,k,2); 2014fad6a16SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 2024fad6a16SStefano Zampini PetscFunctionReturn(0); 2034fad6a16SStefano Zampini } 2042b510759SStefano Zampini 205b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 2062b510759SStefano Zampini #undef __FUNCT__ 207b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 208b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 209b8ffe317SStefano Zampini { 210b8ffe317SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 211b8ffe317SStefano Zampini 212b8ffe317SStefano Zampini PetscFunctionBegin; 21385c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = flg; 214b8ffe317SStefano Zampini PetscFunctionReturn(0); 215b8ffe317SStefano Zampini } 216b8ffe317SStefano Zampini 217b8ffe317SStefano Zampini #undef __FUNCT__ 218b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 219b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 2202b510759SStefano Zampini { 2212b510759SStefano Zampini PetscErrorCode ierr; 2222b510759SStefano Zampini 2232b510759SStefano Zampini PetscFunctionBegin; 2242b510759SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 225b8ffe317SStefano Zampini PetscValidLogicalCollectiveBool(pc,flg,2); 226b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 2272b510759SStefano Zampini PetscFunctionReturn(0); 2282b510759SStefano Zampini } 2291e6b0712SBarry Smith 2304fad6a16SStefano Zampini #undef __FUNCT__ 2312b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC" 2322b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 2334fad6a16SStefano Zampini { 2344fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2354fad6a16SStefano Zampini 2364fad6a16SStefano Zampini PetscFunctionBegin; 2372b510759SStefano Zampini pcbddc->current_level = level; 2384fad6a16SStefano Zampini PetscFunctionReturn(0); 2394fad6a16SStefano Zampini } 2401e6b0712SBarry Smith 2414fad6a16SStefano Zampini #undef __FUNCT__ 242b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel" 243b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 244b8ffe317SStefano Zampini { 245b8ffe317SStefano Zampini PetscErrorCode ierr; 246b8ffe317SStefano Zampini 247b8ffe317SStefano Zampini PetscFunctionBegin; 248b8ffe317SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 249b8ffe317SStefano Zampini PetscValidLogicalCollectiveInt(pc,level,2); 250b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 251b8ffe317SStefano Zampini PetscFunctionReturn(0); 252b8ffe317SStefano Zampini } 253b8ffe317SStefano Zampini 254b8ffe317SStefano Zampini #undef __FUNCT__ 2552b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC" 2562b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 2572b510759SStefano Zampini { 2582b510759SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2592b510759SStefano Zampini 2602b510759SStefano Zampini PetscFunctionBegin; 2612b510759SStefano Zampini pcbddc->max_levels = levels; 2622b510759SStefano Zampini PetscFunctionReturn(0); 2632b510759SStefano Zampini } 2642b510759SStefano Zampini 2652b510759SStefano Zampini #undef __FUNCT__ 2662b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels" 2674fad6a16SStefano Zampini /*@ 26828509bceSStefano Zampini PCBDDCSetLevels - Sets the maximum number of levels for multilevel 2694fad6a16SStefano Zampini 2704fad6a16SStefano Zampini Logically collective on PC 2714fad6a16SStefano Zampini 2724fad6a16SStefano Zampini Input Parameters: 2734fad6a16SStefano Zampini + pc - the preconditioning context 27428509bceSStefano Zampini - levels - the maximum number of levels (max 9) 2754fad6a16SStefano Zampini 27628509bceSStefano Zampini Default value is 0, i.e. traditional one-level BDDC 2774fad6a16SStefano Zampini 2784fad6a16SStefano Zampini Level: intermediate 2794fad6a16SStefano Zampini 2804fad6a16SStefano Zampini Notes: 2814fad6a16SStefano Zampini 2824fad6a16SStefano Zampini .seealso: PCBDDC 2834fad6a16SStefano Zampini @*/ 2842b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 2854fad6a16SStefano Zampini { 2864fad6a16SStefano Zampini PetscErrorCode ierr; 2874fad6a16SStefano Zampini 2884fad6a16SStefano Zampini PetscFunctionBegin; 2894fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2902b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,levels,2); 291312be037SStefano Zampini if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n"); 2922b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 2934fad6a16SStefano Zampini PetscFunctionReturn(0); 2944fad6a16SStefano Zampini } 2954fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */ 2961e6b0712SBarry Smith 2974fad6a16SStefano Zampini #undef __FUNCT__ 2980bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 2990bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 3000bdf917eSStefano Zampini { 3010bdf917eSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3020bdf917eSStefano Zampini PetscErrorCode ierr; 3030bdf917eSStefano Zampini 3040bdf917eSStefano Zampini PetscFunctionBegin; 3050bdf917eSStefano Zampini ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 3060bdf917eSStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 3070bdf917eSStefano Zampini pcbddc->NullSpace = NullSpace; 3080bdf917eSStefano Zampini PetscFunctionReturn(0); 3090bdf917eSStefano Zampini } 3101e6b0712SBarry Smith 3110bdf917eSStefano Zampini #undef __FUNCT__ 3120bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace" 3130bdf917eSStefano Zampini /*@ 31428509bceSStefano Zampini PCBDDCSetNullSpace - Set nullspace for BDDC operator 3150bdf917eSStefano Zampini 3160bdf917eSStefano Zampini Logically collective on PC and MatNullSpace 3170bdf917eSStefano Zampini 3180bdf917eSStefano Zampini Input Parameters: 3190bdf917eSStefano Zampini + pc - the preconditioning context 32028509bceSStefano Zampini - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 3210bdf917eSStefano Zampini 3220bdf917eSStefano Zampini Level: intermediate 3230bdf917eSStefano Zampini 3240bdf917eSStefano Zampini Notes: 3250bdf917eSStefano Zampini 3260bdf917eSStefano Zampini .seealso: PCBDDC 3270bdf917eSStefano Zampini @*/ 3280bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 3290bdf917eSStefano Zampini { 3300bdf917eSStefano Zampini PetscErrorCode ierr; 3310bdf917eSStefano Zampini 3320bdf917eSStefano Zampini PetscFunctionBegin; 3330bdf917eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 334674ae819SStefano Zampini PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 3352b510759SStefano Zampini PetscCheckSameComm(pc,1,NullSpace,2); 3360bdf917eSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 3370bdf917eSStefano Zampini PetscFunctionReturn(0); 3380bdf917eSStefano Zampini } 3390bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */ 3401e6b0712SBarry Smith 3410bdf917eSStefano Zampini #undef __FUNCT__ 3423b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 3433b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 3443b03a366Sstefano_zampini { 3453b03a366Sstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3463b03a366Sstefano_zampini PetscErrorCode ierr; 3473b03a366Sstefano_zampini 3483b03a366Sstefano_zampini PetscFunctionBegin; 349785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 350785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 3513b03a366Sstefano_zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 35236e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 35336e030ebSStefano Zampini pcbddc->DirichletBoundaries = DirichletBoundaries; 354fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 3553b03a366Sstefano_zampini PetscFunctionReturn(0); 3563b03a366Sstefano_zampini } 3571e6b0712SBarry Smith 3583b03a366Sstefano_zampini #undef __FUNCT__ 3593b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 3603b03a366Sstefano_zampini /*@ 36128509bceSStefano Zampini PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 3623b03a366Sstefano_zampini 363785d1243SStefano Zampini Collective 3643b03a366Sstefano_zampini 3653b03a366Sstefano_zampini Input Parameters: 3663b03a366Sstefano_zampini + pc - the preconditioning context 367785d1243SStefano Zampini - DirichletBoundaries - parallel IS defining the Dirichlet boundaries 3683b03a366Sstefano_zampini 3693b03a366Sstefano_zampini Level: intermediate 3703b03a366Sstefano_zampini 371785d1243SStefano Zampini Notes: Any process can list any global node 3723b03a366Sstefano_zampini 3733b03a366Sstefano_zampini .seealso: PCBDDC 3743b03a366Sstefano_zampini @*/ 3753b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 3763b03a366Sstefano_zampini { 3773b03a366Sstefano_zampini PetscErrorCode ierr; 3783b03a366Sstefano_zampini 3793b03a366Sstefano_zampini PetscFunctionBegin; 3803b03a366Sstefano_zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 381674ae819SStefano Zampini PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 382785d1243SStefano Zampini PetscCheckSameComm(pc,1,DirichletBoundaries,2); 3833b03a366Sstefano_zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 3843b03a366Sstefano_zampini PetscFunctionReturn(0); 3853b03a366Sstefano_zampini } 3863b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */ 3871e6b0712SBarry Smith 3883b03a366Sstefano_zampini #undef __FUNCT__ 38982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC" 39082ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries) 3910c7d97c5SJed Brown { 3920c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3930c7d97c5SJed Brown PetscErrorCode ierr; 3940c7d97c5SJed Brown 3950c7d97c5SJed Brown PetscFunctionBegin; 396785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 397785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 3980c7d97c5SJed Brown ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 3990c7d97c5SJed Brown ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 400785d1243SStefano Zampini pcbddc->DirichletBoundariesLocal = DirichletBoundaries; 4017d24d155SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4020c7d97c5SJed Brown PetscFunctionReturn(0); 4030c7d97c5SJed Brown } 4040c7d97c5SJed Brown 4050c7d97c5SJed Brown #undef __FUNCT__ 40682ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal" 4079c0446d6SStefano Zampini /*@ 40882ba6b80SStefano Zampini PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering. 4099c0446d6SStefano Zampini 410785d1243SStefano Zampini Collective 41153cdbc3dSStefano Zampini 41253cdbc3dSStefano Zampini Input Parameters: 41353cdbc3dSStefano Zampini + pc - the preconditioning context 41482ba6b80SStefano Zampini - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering) 41553cdbc3dSStefano Zampini 41653cdbc3dSStefano Zampini Level: intermediate 41753cdbc3dSStefano Zampini 4189c0446d6SStefano Zampini Notes: 41953cdbc3dSStefano Zampini 42053cdbc3dSStefano Zampini .seealso: PCBDDC 42153cdbc3dSStefano Zampini @*/ 42282ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries) 4230c7d97c5SJed Brown { 4240c7d97c5SJed Brown PetscErrorCode ierr; 4250c7d97c5SJed Brown 4260c7d97c5SJed Brown PetscFunctionBegin; 4270c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4280c7d97c5SJed Brown PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 42982ba6b80SStefano Zampini PetscCheckSameComm(pc,1,DirichletBoundaries,2); 43082ba6b80SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 4310c7d97c5SJed Brown PetscFunctionReturn(0); 4320c7d97c5SJed Brown } 4330c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 4340c7d97c5SJed Brown 4350c7d97c5SJed Brown #undef __FUNCT__ 4360c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 43753cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 4380c7d97c5SJed Brown { 4390c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 44053cdbc3dSStefano Zampini PetscErrorCode ierr; 4410c7d97c5SJed Brown 4420c7d97c5SJed Brown PetscFunctionBegin; 443785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 444785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 44553cdbc3dSStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 44636e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 44736e030ebSStefano Zampini pcbddc->NeumannBoundaries = NeumannBoundaries; 448fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4490c7d97c5SJed Brown PetscFunctionReturn(0); 4500c7d97c5SJed Brown } 4511e6b0712SBarry Smith 4520c7d97c5SJed Brown #undef __FUNCT__ 4530c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 45457527edcSJed Brown /*@ 45528509bceSStefano Zampini PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 45657527edcSJed Brown 457785d1243SStefano Zampini Collective 45857527edcSJed Brown 45957527edcSJed Brown Input Parameters: 46057527edcSJed Brown + pc - the preconditioning context 461785d1243SStefano Zampini - NeumannBoundaries - parallel IS defining the Neumann boundaries 46257527edcSJed Brown 46357527edcSJed Brown Level: intermediate 46457527edcSJed Brown 465785d1243SStefano Zampini Notes: Any process can list any global node 46657527edcSJed Brown 46757527edcSJed Brown .seealso: PCBDDC 46857527edcSJed Brown @*/ 46953cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 4700c7d97c5SJed Brown { 4710c7d97c5SJed Brown PetscErrorCode ierr; 4720c7d97c5SJed Brown 4730c7d97c5SJed Brown PetscFunctionBegin; 4740c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 475674ae819SStefano Zampini PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 476785d1243SStefano Zampini PetscCheckSameComm(pc,1,NeumannBoundaries,2); 47753cdbc3dSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 47853cdbc3dSStefano Zampini PetscFunctionReturn(0); 47953cdbc3dSStefano Zampini } 48053cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 4811e6b0712SBarry Smith 48253cdbc3dSStefano Zampini #undef __FUNCT__ 48382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC" 48482ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries) 4850c7d97c5SJed Brown { 4860c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 4870c7d97c5SJed Brown PetscErrorCode ierr; 4880c7d97c5SJed Brown 4890c7d97c5SJed Brown PetscFunctionBegin; 490785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 491785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 4920c7d97c5SJed Brown ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 4930c7d97c5SJed Brown ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 494785d1243SStefano Zampini pcbddc->NeumannBoundariesLocal = NeumannBoundaries; 4957d24d155SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4960c7d97c5SJed Brown PetscFunctionReturn(0); 4970c7d97c5SJed Brown } 4980c7d97c5SJed Brown 4990c7d97c5SJed Brown #undef __FUNCT__ 50082ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal" 5010c7d97c5SJed Brown /*@ 50282ba6b80SStefano Zampini PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering. 5030c7d97c5SJed Brown 504785d1243SStefano Zampini Collective 5050c7d97c5SJed Brown 5060c7d97c5SJed Brown Input Parameters: 5070c7d97c5SJed Brown + pc - the preconditioning context 50882ba6b80SStefano Zampini - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering) 5090c7d97c5SJed Brown 5100c7d97c5SJed Brown Level: intermediate 5110c7d97c5SJed Brown 5120c7d97c5SJed Brown Notes: 5130c7d97c5SJed Brown 5140c7d97c5SJed Brown .seealso: PCBDDC 5150c7d97c5SJed Brown @*/ 51682ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries) 5170c7d97c5SJed Brown { 5180c7d97c5SJed Brown PetscErrorCode ierr; 5190c7d97c5SJed Brown 5200c7d97c5SJed Brown PetscFunctionBegin; 5210c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5220c7d97c5SJed Brown PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 52382ba6b80SStefano Zampini PetscCheckSameComm(pc,1,NeumannBoundaries,2); 52482ba6b80SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 52553cdbc3dSStefano Zampini PetscFunctionReturn(0); 52653cdbc3dSStefano Zampini } 52753cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 52853cdbc3dSStefano Zampini 52953cdbc3dSStefano Zampini #undef __FUNCT__ 530da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 531da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 532da1bb401SStefano Zampini { 533da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 534da1bb401SStefano Zampini 535da1bb401SStefano Zampini PetscFunctionBegin; 536da1bb401SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundaries; 537da1bb401SStefano Zampini PetscFunctionReturn(0); 538da1bb401SStefano Zampini } 5391e6b0712SBarry Smith 540da1bb401SStefano Zampini #undef __FUNCT__ 541da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 542da1bb401SStefano Zampini /*@ 543785d1243SStefano Zampini PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries 544da1bb401SStefano Zampini 545785d1243SStefano Zampini Collective 546785d1243SStefano Zampini 547785d1243SStefano Zampini Input Parameters: 548785d1243SStefano Zampini . pc - the preconditioning context 549785d1243SStefano Zampini 550785d1243SStefano Zampini Output Parameters: 551785d1243SStefano Zampini . DirichletBoundaries - index set defining the Dirichlet boundaries 552785d1243SStefano Zampini 553785d1243SStefano Zampini Level: intermediate 554785d1243SStefano Zampini 555785d1243SStefano Zampini Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries 556785d1243SStefano Zampini 557785d1243SStefano Zampini .seealso: PCBDDC 558785d1243SStefano Zampini @*/ 559785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 560785d1243SStefano Zampini { 561785d1243SStefano Zampini PetscErrorCode ierr; 562785d1243SStefano Zampini 563785d1243SStefano Zampini PetscFunctionBegin; 564785d1243SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 565785d1243SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 566785d1243SStefano Zampini PetscFunctionReturn(0); 567785d1243SStefano Zampini } 568785d1243SStefano Zampini /* -------------------------------------------------------------------------- */ 569785d1243SStefano Zampini 570785d1243SStefano Zampini #undef __FUNCT__ 571785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC" 572785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries) 573785d1243SStefano Zampini { 574785d1243SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 575785d1243SStefano Zampini 576785d1243SStefano Zampini PetscFunctionBegin; 577785d1243SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundariesLocal; 578785d1243SStefano Zampini PetscFunctionReturn(0); 579785d1243SStefano Zampini } 580785d1243SStefano Zampini 581785d1243SStefano Zampini #undef __FUNCT__ 58282ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal" 583da1bb401SStefano Zampini /*@ 58482ba6b80SStefano Zampini PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering) 585da1bb401SStefano Zampini 586785d1243SStefano Zampini Collective 587da1bb401SStefano Zampini 588da1bb401SStefano Zampini Input Parameters: 58928509bceSStefano Zampini . pc - the preconditioning context 590da1bb401SStefano Zampini 591da1bb401SStefano Zampini Output Parameters: 59228509bceSStefano Zampini . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 593da1bb401SStefano Zampini 594da1bb401SStefano Zampini Level: intermediate 595da1bb401SStefano Zampini 596da1bb401SStefano Zampini Notes: 597da1bb401SStefano Zampini 598da1bb401SStefano Zampini .seealso: PCBDDC 599da1bb401SStefano Zampini @*/ 60082ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries) 601da1bb401SStefano Zampini { 602da1bb401SStefano Zampini PetscErrorCode ierr; 603da1bb401SStefano Zampini 604da1bb401SStefano Zampini PetscFunctionBegin; 605da1bb401SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 60682ba6b80SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 607da1bb401SStefano Zampini PetscFunctionReturn(0); 608da1bb401SStefano Zampini } 609da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 6101e6b0712SBarry Smith 611da1bb401SStefano Zampini #undef __FUNCT__ 61253cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 61353cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 61453cdbc3dSStefano Zampini { 61553cdbc3dSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 61653cdbc3dSStefano Zampini 61753cdbc3dSStefano Zampini PetscFunctionBegin; 61853cdbc3dSStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundaries; 61953cdbc3dSStefano Zampini PetscFunctionReturn(0); 62053cdbc3dSStefano Zampini } 6211e6b0712SBarry Smith 62253cdbc3dSStefano Zampini #undef __FUNCT__ 62353cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 62453cdbc3dSStefano Zampini /*@ 625785d1243SStefano Zampini PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries 62653cdbc3dSStefano Zampini 627785d1243SStefano Zampini Collective 628785d1243SStefano Zampini 629785d1243SStefano Zampini Input Parameters: 630785d1243SStefano Zampini . pc - the preconditioning context 631785d1243SStefano Zampini 632785d1243SStefano Zampini Output Parameters: 633785d1243SStefano Zampini . NeumannBoundaries - index set defining the Neumann boundaries 634785d1243SStefano Zampini 635785d1243SStefano Zampini Level: intermediate 636785d1243SStefano Zampini 637785d1243SStefano Zampini Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries 638785d1243SStefano Zampini 639785d1243SStefano Zampini .seealso: PCBDDC 640785d1243SStefano Zampini @*/ 641785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 642785d1243SStefano Zampini { 643785d1243SStefano Zampini PetscErrorCode ierr; 644785d1243SStefano Zampini 645785d1243SStefano Zampini PetscFunctionBegin; 646785d1243SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 647785d1243SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 648785d1243SStefano Zampini PetscFunctionReturn(0); 649785d1243SStefano Zampini } 650785d1243SStefano Zampini /* -------------------------------------------------------------------------- */ 651785d1243SStefano Zampini 652785d1243SStefano Zampini #undef __FUNCT__ 653785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC" 654785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries) 655785d1243SStefano Zampini { 656785d1243SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 657785d1243SStefano Zampini 658785d1243SStefano Zampini PetscFunctionBegin; 659785d1243SStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundariesLocal; 660785d1243SStefano Zampini PetscFunctionReturn(0); 661785d1243SStefano Zampini } 662785d1243SStefano Zampini 663785d1243SStefano Zampini #undef __FUNCT__ 66482ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal" 66553cdbc3dSStefano Zampini /*@ 66682ba6b80SStefano Zampini PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering) 66753cdbc3dSStefano Zampini 668785d1243SStefano Zampini Collective 66953cdbc3dSStefano Zampini 67053cdbc3dSStefano Zampini Input Parameters: 67128509bceSStefano Zampini . pc - the preconditioning context 67253cdbc3dSStefano Zampini 67353cdbc3dSStefano Zampini Output Parameters: 67428509bceSStefano Zampini . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 67553cdbc3dSStefano Zampini 67653cdbc3dSStefano Zampini Level: intermediate 67753cdbc3dSStefano Zampini 67853cdbc3dSStefano Zampini Notes: 67953cdbc3dSStefano Zampini 68053cdbc3dSStefano Zampini .seealso: PCBDDC 68153cdbc3dSStefano Zampini @*/ 68282ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries) 68353cdbc3dSStefano Zampini { 68453cdbc3dSStefano Zampini PetscErrorCode ierr; 68553cdbc3dSStefano Zampini 68653cdbc3dSStefano Zampini PetscFunctionBegin; 68753cdbc3dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 68882ba6b80SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 6890c7d97c5SJed Brown PetscFunctionReturn(0); 6900c7d97c5SJed Brown } 69136e030ebSStefano Zampini /* -------------------------------------------------------------------------- */ 6921e6b0712SBarry Smith 69336e030ebSStefano Zampini #undef __FUNCT__ 694da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 6951a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 69636e030ebSStefano Zampini { 69736e030ebSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 698da1bb401SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 699da1bb401SStefano Zampini PetscErrorCode ierr; 70036e030ebSStefano Zampini 70136e030ebSStefano Zampini PetscFunctionBegin; 702674ae819SStefano Zampini /* free old CSR */ 703674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 704fb223d50SStefano Zampini /* TODO: PCBDDCGraphSetAdjacency */ 705674ae819SStefano Zampini /* get CSR into graph structure */ 706da1bb401SStefano Zampini if (copymode == PETSC_COPY_VALUES) { 707785e854fSJed Brown ierr = PetscMalloc1((nvtxs+1),&mat_graph->xadj);CHKERRQ(ierr); 708785e854fSJed Brown ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr); 709674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 710674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 711da1bb401SStefano Zampini } else if (copymode == PETSC_OWN_POINTER) { 7121a83f524SJed Brown mat_graph->xadj = (PetscInt*)xadj; 7131a83f524SJed Brown mat_graph->adjncy = (PetscInt*)adjncy; 714674ae819SStefano Zampini } 715575ad6abSStefano Zampini mat_graph->nvtxs_csr = nvtxs; 71636e030ebSStefano Zampini PetscFunctionReturn(0); 71736e030ebSStefano Zampini } 7181e6b0712SBarry Smith 71936e030ebSStefano Zampini #undef __FUNCT__ 720da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 72136e030ebSStefano Zampini /*@ 72228509bceSStefano Zampini PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix 72336e030ebSStefano Zampini 72436e030ebSStefano Zampini Not collective 72536e030ebSStefano Zampini 72636e030ebSStefano Zampini Input Parameters: 72736e030ebSStefano Zampini + pc - the preconditioning context 72828509bceSStefano Zampini . nvtxs - number of local vertices of the graph (i.e., the local size of your problem) 72928509bceSStefano Zampini . xadj, adjncy - the CSR graph 73028509bceSStefano Zampini - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 73136e030ebSStefano Zampini 73236e030ebSStefano Zampini Level: intermediate 73336e030ebSStefano Zampini 73436e030ebSStefano Zampini Notes: 73536e030ebSStefano Zampini 73628509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode 73736e030ebSStefano Zampini @*/ 7381a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 73936e030ebSStefano Zampini { 740575ad6abSStefano Zampini void (*f)(void) = 0; 74136e030ebSStefano Zampini PetscErrorCode ierr; 74236e030ebSStefano Zampini 74336e030ebSStefano Zampini PetscFunctionBegin; 74436e030ebSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 745674ae819SStefano Zampini PetscValidIntPointer(xadj,3); 746674ae819SStefano Zampini PetscValidIntPointer(xadj,4); 747674ae819SStefano Zampini if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 748674ae819SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 749da1bb401SStefano Zampini } 75036e030ebSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 751575ad6abSStefano Zampini /* free arrays if PCBDDC is not the PC type */ 752575ad6abSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 753575ad6abSStefano Zampini if (!f && copymode == PETSC_OWN_POINTER) { 754575ad6abSStefano Zampini ierr = PetscFree(xadj);CHKERRQ(ierr); 755575ad6abSStefano Zampini ierr = PetscFree(adjncy);CHKERRQ(ierr); 75636e030ebSStefano Zampini } 75736e030ebSStefano Zampini PetscFunctionReturn(0); 75836e030ebSStefano Zampini } 7599c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */ 7601e6b0712SBarry Smith 7619c0446d6SStefano Zampini #undef __FUNCT__ 76263602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC" 76363602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 76463602bcaSStefano Zampini { 76563602bcaSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 76663602bcaSStefano Zampini PetscInt i; 76763602bcaSStefano Zampini PetscErrorCode ierr; 76863602bcaSStefano Zampini 76963602bcaSStefano Zampini PetscFunctionBegin; 77063602bcaSStefano Zampini /* Destroy ISes if they were already set */ 77163602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 77263602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 77363602bcaSStefano Zampini } 77463602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 77563602bcaSStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 77663602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 77763602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 77863602bcaSStefano Zampini } 77963602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 78063602bcaSStefano Zampini pcbddc->n_ISForDofs = 0; 78163602bcaSStefano Zampini /* allocate space then set */ 78263602bcaSStefano Zampini ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr); 78363602bcaSStefano Zampini for (i=0;i<n_is;i++) { 78463602bcaSStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 78563602bcaSStefano Zampini pcbddc->ISForDofsLocal[i]=ISForDofs[i]; 78663602bcaSStefano Zampini } 78763602bcaSStefano Zampini pcbddc->n_ISForDofsLocal=n_is; 78863602bcaSStefano Zampini if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 7891a8a16d1SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 79063602bcaSStefano Zampini PetscFunctionReturn(0); 79163602bcaSStefano Zampini } 79263602bcaSStefano Zampini 79363602bcaSStefano Zampini #undef __FUNCT__ 79463602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal" 79563602bcaSStefano Zampini /*@ 79663602bcaSStefano Zampini PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix 79763602bcaSStefano Zampini 79863602bcaSStefano Zampini Collective 79963602bcaSStefano Zampini 80063602bcaSStefano Zampini Input Parameters: 80163602bcaSStefano Zampini + pc - the preconditioning context 80263602bcaSStefano Zampini - n_is - number of index sets defining the fields 80363602bcaSStefano Zampini . ISForDofs - array of IS describing the fields in local ordering 80463602bcaSStefano Zampini 80563602bcaSStefano Zampini Level: intermediate 80663602bcaSStefano Zampini 80763602bcaSStefano Zampini Notes: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to a different field. 80863602bcaSStefano Zampini 80963602bcaSStefano Zampini .seealso: PCBDDC 81063602bcaSStefano Zampini @*/ 81163602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[]) 81263602bcaSStefano Zampini { 81363602bcaSStefano Zampini PetscInt i; 81463602bcaSStefano Zampini PetscErrorCode ierr; 81563602bcaSStefano Zampini 81663602bcaSStefano Zampini PetscFunctionBegin; 81763602bcaSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 81863602bcaSStefano Zampini PetscValidLogicalCollectiveInt(pc,n_is,2); 81963602bcaSStefano Zampini for (i=0;i<n_is;i++) { 82063602bcaSStefano Zampini PetscCheckSameComm(pc,1,ISForDofs[i],3); 82163602bcaSStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 82263602bcaSStefano Zampini } 82363602bcaSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 82463602bcaSStefano Zampini PetscFunctionReturn(0); 82563602bcaSStefano Zampini } 82663602bcaSStefano Zampini /* -------------------------------------------------------------------------- */ 82763602bcaSStefano Zampini 82863602bcaSStefano Zampini #undef __FUNCT__ 8299c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 8309c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 8319c0446d6SStefano Zampini { 8329c0446d6SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 8339c0446d6SStefano Zampini PetscInt i; 8349c0446d6SStefano Zampini PetscErrorCode ierr; 8359c0446d6SStefano Zampini 8369c0446d6SStefano Zampini PetscFunctionBegin; 837da1bb401SStefano Zampini /* Destroy ISes if they were already set */ 8389c0446d6SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 8399c0446d6SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 8409c0446d6SStefano Zampini } 841d11ae9bbSstefano_zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 84263602bcaSStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 84363602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 84463602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 84563602bcaSStefano Zampini } 84663602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 84763602bcaSStefano Zampini pcbddc->n_ISForDofsLocal = 0; 848da1bb401SStefano Zampini /* allocate space then set */ 849785e854fSJed Brown ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr); 8509c0446d6SStefano Zampini for (i=0;i<n_is;i++) { 851da1bb401SStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 852da1bb401SStefano Zampini pcbddc->ISForDofs[i]=ISForDofs[i]; 8539c0446d6SStefano Zampini } 8549c0446d6SStefano Zampini pcbddc->n_ISForDofs=n_is; 85563602bcaSStefano Zampini if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 8561a8a16d1SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 8579c0446d6SStefano Zampini PetscFunctionReturn(0); 8589c0446d6SStefano Zampini } 8591e6b0712SBarry Smith 8609c0446d6SStefano Zampini #undef __FUNCT__ 8619c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting" 8629c0446d6SStefano Zampini /*@ 86363602bcaSStefano Zampini PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix 8649c0446d6SStefano Zampini 86563602bcaSStefano Zampini Collective 8669c0446d6SStefano Zampini 8679c0446d6SStefano Zampini Input Parameters: 8689c0446d6SStefano Zampini + pc - the preconditioning context 86928509bceSStefano Zampini - n_is - number of index sets defining the fields 87063602bcaSStefano Zampini . ISForDofs - array of IS describing the fields in global ordering 8719c0446d6SStefano Zampini 8729c0446d6SStefano Zampini Level: intermediate 8739c0446d6SStefano Zampini 87463602bcaSStefano Zampini Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field. 8759c0446d6SStefano Zampini 8769c0446d6SStefano Zampini .seealso: PCBDDC 8779c0446d6SStefano Zampini @*/ 8789c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 8799c0446d6SStefano Zampini { 8802b510759SStefano Zampini PetscInt i; 8819c0446d6SStefano Zampini PetscErrorCode ierr; 8829c0446d6SStefano Zampini 8839c0446d6SStefano Zampini PetscFunctionBegin; 8849c0446d6SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 88563602bcaSStefano Zampini PetscValidLogicalCollectiveInt(pc,n_is,2); 8862b510759SStefano Zampini for (i=0;i<n_is;i++) { 88763602bcaSStefano Zampini PetscCheckSameComm(pc,1,ISForDofs[i],3); 88863602bcaSStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 8892b510759SStefano Zampini } 8909c0446d6SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 8919c0446d6SStefano Zampini PetscFunctionReturn(0); 8929c0446d6SStefano Zampini } 893da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 894534831adSStefano Zampini #undef __FUNCT__ 895534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC" 896534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 897534831adSStefano Zampini /* 898534831adSStefano Zampini PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 899534831adSStefano Zampini guess if a transformation of basis approach has been selected. 9009c0446d6SStefano Zampini 901534831adSStefano Zampini Input Parameter: 902534831adSStefano Zampini + pc - the preconditioner contex 903534831adSStefano Zampini 904534831adSStefano Zampini Application Interface Routine: PCPreSolve() 905534831adSStefano Zampini 906534831adSStefano Zampini Notes: 907534831adSStefano Zampini The interface routine PCPreSolve() is not usually called directly by 908534831adSStefano Zampini the user, but instead is called by KSPSolve(). 909534831adSStefano Zampini */ 910534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 911534831adSStefano Zampini { 912534831adSStefano Zampini PetscErrorCode ierr; 913534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 914534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 915534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 916534831adSStefano Zampini Mat temp_mat; 9173972b0daSStefano Zampini IS dirIS; 9183972b0daSStefano Zampini Vec used_vec; 91982ba6b80SStefano Zampini PetscBool guess_nonzero; 920534831adSStefano Zampini 921534831adSStefano Zampini PetscFunctionBegin; 92285c4d303SStefano Zampini /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 92385c4d303SStefano Zampini if (ksp) { 92485c4d303SStefano Zampini PetscBool iscg; 92585c4d303SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 92685c4d303SStefano Zampini if (!iscg) { 92785c4d303SStefano Zampini ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 92885c4d303SStefano Zampini } 92985c4d303SStefano Zampini } 93085c4d303SStefano Zampini /* Creates parallel work vectors used in presolve */ 93162a6ff1dSStefano Zampini if (!pcbddc->original_rhs) { 93262a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 93362a6ff1dSStefano Zampini } 93462a6ff1dSStefano Zampini if (!pcbddc->temp_solution) { 93562a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 93662a6ff1dSStefano Zampini } 9373972b0daSStefano Zampini if (x) { 9383972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 9393972b0daSStefano Zampini used_vec = x; 9403972b0daSStefano Zampini } else { 9413972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 9423972b0daSStefano Zampini used_vec = pcbddc->temp_solution; 9433972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 9443972b0daSStefano Zampini } 9453972b0daSStefano Zampini /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 9463972b0daSStefano Zampini if (ksp) { 9473972b0daSStefano Zampini ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 9483972b0daSStefano Zampini if (!guess_nonzero) { 9493972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 9503972b0daSStefano Zampini } 9513972b0daSStefano Zampini } 9523308cffdSStefano Zampini 9533972b0daSStefano Zampini /* store the original rhs */ 9543972b0daSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 9553972b0daSStefano Zampini 9563972b0daSStefano Zampini /* Take into account zeroed rows -> change rhs and store solution removed */ 957785d1243SStefano Zampini /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */ 95882ba6b80SStefano Zampini ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr); 95982ba6b80SStefano Zampini if (rhs && dirIS) { 960785d1243SStefano Zampini PetscInt dirsize,i,*is_indices; 961785d1243SStefano Zampini PetscScalar *array_x,*array_diagonal; 962785d1243SStefano Zampini 9633972b0daSStefano Zampini ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 9643972b0daSStefano Zampini ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 9653972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9663972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9673972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9683972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 96982ba6b80SStefano Zampini ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr); 9703972b0daSStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 9713972b0daSStefano Zampini ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 9723972b0daSStefano Zampini ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 9732fa5cd67SKarl Rupp for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 9743972b0daSStefano Zampini ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 9753972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 9763972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 9773972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9783972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 979b76ba322SStefano Zampini 9803972b0daSStefano Zampini /* remove the computed solution from the rhs */ 9813972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 9823972b0daSStefano Zampini ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 9833972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 9843308cffdSStefano Zampini } 985b76ba322SStefano Zampini 986b76ba322SStefano Zampini /* store partially computed solution and set initial guess */ 9873972b0daSStefano Zampini if (x) { 9883972b0daSStefano Zampini ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 9893972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 99085c4d303SStefano Zampini if (pcbddc->use_exact_dirichlet_trick) { 991b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 992b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 993b76ba322SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 994b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 995b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 996b76ba322SStefano Zampini if (ksp) { 997b76ba322SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 998b76ba322SStefano Zampini } 999b76ba322SStefano Zampini } 10003972b0daSStefano Zampini } 1001b76ba322SStefano Zampini 100292e3dcfbSStefano Zampini /* prepare MatMult and rhs for solver */ 1003*b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix) { 1004b76ba322SStefano Zampini /* swap pointers for local matrices */ 1005b76ba322SStefano Zampini temp_mat = matis->A; 1006b76ba322SStefano Zampini matis->A = pcbddc->local_mat; 1007b76ba322SStefano Zampini pcbddc->local_mat = temp_mat; 100892e3dcfbSStefano Zampini if (rhs) { 1009b76ba322SStefano Zampini /* Get local rhs and apply transformation of basis */ 1010b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1011b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1012b76ba322SStefano Zampini /* from original basis to modified basis */ 1013b76ba322SStefano Zampini ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 1014b76ba322SStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 1015b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1016b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1017674ae819SStefano Zampini } 101892e3dcfbSStefano Zampini } 101992e3dcfbSStefano Zampini 102092e3dcfbSStefano Zampini /* remove nullspace if present */ 10210bdf917eSStefano Zampini if (ksp && pcbddc->NullSpace) { 1022d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 1023d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 1024b76ba322SStefano Zampini } 10250bdf917eSStefano Zampini ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 1026534831adSStefano Zampini PetscFunctionReturn(0); 1027534831adSStefano Zampini } 1028534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 1029534831adSStefano Zampini #undef __FUNCT__ 1030534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC" 1031534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 1032534831adSStefano Zampini /* 1033534831adSStefano Zampini PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 1034534831adSStefano Zampini approach has been selected. Also, restores rhs to its original state. 1035534831adSStefano Zampini 1036534831adSStefano Zampini Input Parameter: 1037534831adSStefano Zampini + pc - the preconditioner contex 1038534831adSStefano Zampini 1039534831adSStefano Zampini Application Interface Routine: PCPostSolve() 1040534831adSStefano Zampini 1041534831adSStefano Zampini Notes: 1042534831adSStefano Zampini The interface routine PCPostSolve() is not usually called directly by 1043534831adSStefano Zampini the user, but instead is called by KSPSolve(). 1044534831adSStefano Zampini */ 1045534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1046534831adSStefano Zampini { 1047534831adSStefano Zampini PetscErrorCode ierr; 1048534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1049534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 1050534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 1051534831adSStefano Zampini Mat temp_mat; 1052534831adSStefano Zampini 1053534831adSStefano Zampini PetscFunctionBegin; 1054*b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix) { 1055534831adSStefano Zampini /* swap pointers for local matrices */ 1056534831adSStefano Zampini temp_mat = matis->A; 1057534831adSStefano Zampini matis->A = pcbddc->local_mat; 1058534831adSStefano Zampini pcbddc->local_mat = temp_mat; 10593425bc38SStefano Zampini } 1060*b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix && x) { 1061534831adSStefano Zampini /* Get Local boundary and apply transformation of basis to solution vector */ 1062534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1063534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1064534831adSStefano Zampini /* from modified basis to original basis */ 1065534831adSStefano Zampini ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 1066534831adSStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 1067534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1068534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1069534831adSStefano Zampini } 10703972b0daSStefano Zampini /* add solution removed in presolve */ 10713425bc38SStefano Zampini if (x) { 10723425bc38SStefano Zampini ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 10733425bc38SStefano Zampini } 1074fb223d50SStefano Zampini /* restore rhs to its original state */ 1075fb223d50SStefano Zampini if (rhs) { 1076fb223d50SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 1077fb223d50SStefano Zampini } 1078534831adSStefano Zampini PetscFunctionReturn(0); 1079534831adSStefano Zampini } 1080534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 108153cdbc3dSStefano Zampini #undef __FUNCT__ 108253cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC" 10830c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 10840c7d97c5SJed Brown /* 10850c7d97c5SJed Brown PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 10860c7d97c5SJed Brown by setting data structures and options. 10870c7d97c5SJed Brown 10880c7d97c5SJed Brown Input Parameter: 108953cdbc3dSStefano Zampini + pc - the preconditioner context 10900c7d97c5SJed Brown 10910c7d97c5SJed Brown Application Interface Routine: PCSetUp() 10920c7d97c5SJed Brown 10930c7d97c5SJed Brown Notes: 10940c7d97c5SJed Brown The interface routine PCSetUp() is not usually called directly by 10950c7d97c5SJed Brown the user, but instead is called by PCApply() if necessary. 10960c7d97c5SJed Brown */ 109753cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc) 10980c7d97c5SJed Brown { 10990c7d97c5SJed Brown PetscErrorCode ierr; 11000c7d97c5SJed Brown PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 1101f4ddd8eeSStefano Zampini MatNullSpace nearnullspace; 1102674ae819SStefano Zampini PetscBool computeis,computetopography,computesolvers; 1103165b64e2SStefano Zampini PetscBool new_nearnullspace_provided; 11040c7d97c5SJed Brown 11050c7d97c5SJed Brown PetscFunctionBegin; 1106f4ddd8eeSStefano Zampini /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 11073b03a366Sstefano_zampini /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1108fcd409f5SStefano Zampini Also, BDDC directly build the Dirichlet problem */ 1109f4ddd8eeSStefano Zampini /* split work */ 1110674ae819SStefano Zampini if (pc->setupcalled) { 1111674ae819SStefano Zampini computeis = PETSC_FALSE; 1112d1e9a80fSBarry Smith if (pc->flag == SAME_NONZERO_PATTERN) { 1113674ae819SStefano Zampini computetopography = PETSC_FALSE; 1114674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1115674ae819SStefano Zampini } else { /* DIFFERENT_NONZERO_PATTERN */ 1116674ae819SStefano Zampini computetopography = PETSC_TRUE; 1117674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1118674ae819SStefano Zampini } 1119674ae819SStefano Zampini } else { 1120674ae819SStefano Zampini computeis = PETSC_TRUE; 1121674ae819SStefano Zampini computetopography = PETSC_TRUE; 1122674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1123674ae819SStefano Zampini } 1124fb180af4SStefano Zampini if (pcbddc->recompute_topography) { 1125fb180af4SStefano Zampini computetopography = PETSC_TRUE; 1126fb180af4SStefano Zampini } 1127f4ddd8eeSStefano Zampini 1128f4ddd8eeSStefano Zampini /* Get stdout for dbg */ 112970cf5478SStefano Zampini if (pcbddc->dbg_flag) { 113070cf5478SStefano Zampini if (!pcbddc->dbg_viewer) { 113158a03d70SStefano Zampini pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)); 1132f4ddd8eeSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 113370cf5478SStefano Zampini } 113458a03d70SStefano Zampini ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1135f4ddd8eeSStefano Zampini } 1136f4ddd8eeSStefano Zampini 1137fcd409f5SStefano Zampini /* Set up all the "iterative substructuring" common block without computing solvers */ 1138674ae819SStefano Zampini if (computeis) { 1139fcd409f5SStefano Zampini /* HACK INTO PCIS */ 1140fcd409f5SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 1141fcd409f5SStefano Zampini pcis->computesolvers = PETSC_FALSE; 1142674ae819SStefano Zampini ierr = PCISSetUp(pc);CHKERRQ(ierr); 114339e2fb2aSStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr); 1144674ae819SStefano Zampini } 1145f4ddd8eeSStefano Zampini 1146*b9b85e73SStefano Zampini /* check user defined change of basis (if any) */ 1147*b9b85e73SStefano Zampini if (pcbddc->user_ChangeOfBasisMatrix) { 1148*b9b85e73SStefano Zampini PC_IS* pcis= (PC_IS*)pc->data; 1149*b9b85e73SStefano Zampini PetscInt n,m; 1150*b9b85e73SStefano Zampini ierr = MatGetSize(pcbddc->user_ChangeOfBasisMatrix,&n,&m);CHKERRQ(ierr); 1151*b9b85e73SStefano Zampini if (n != m) { 1152*b9b85e73SStefano Zampini SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Change of basis matrix should be square %d != %d\n",n,m); 1153*b9b85e73SStefano Zampini } else if (n != pcis->n_B && n != pcis->n) { 1154*b9b85e73SStefano Zampini SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Size of change of basis matrix %d differs from boundary size %d and local size %d\n",n,pcis->n_B,pcis->n); 1155*b9b85e73SStefano Zampini } 1156*b9b85e73SStefano Zampini /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */ 1157*b9b85e73SStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 1158*b9b85e73SStefano Zampini } 1159f4ddd8eeSStefano Zampini /* Analyze interface */ 1160674ae819SStefano Zampini if (computetopography) { 1161674ae819SStefano Zampini ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 1162fb8d54d4SStefano Zampini } 1163fb8d54d4SStefano Zampini 1164f4ddd8eeSStefano Zampini /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1165fb8d54d4SStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 1166f4ddd8eeSStefano Zampini ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr); 1167f4ddd8eeSStefano Zampini if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1168f4ddd8eeSStefano Zampini if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1169f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1170f4ddd8eeSStefano Zampini } else { 1171f4ddd8eeSStefano Zampini /* determine if the two nullspaces are different (should be lightweight) */ 1172f4ddd8eeSStefano Zampini if (nearnullspace != pcbddc->onearnullspace) { 1173f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1174165b64e2SStefano Zampini } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1175f4ddd8eeSStefano Zampini PetscInt i; 1176165b64e2SStefano Zampini const Vec *nearnullvecs; 1177165b64e2SStefano Zampini PetscObjectState state; 1178165b64e2SStefano Zampini PetscInt nnsp_size; 1179165b64e2SStefano Zampini ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 1180f4ddd8eeSStefano Zampini for (i=0;i<nnsp_size;i++) { 1181f4ddd8eeSStefano Zampini ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr); 1182165b64e2SStefano Zampini if (pcbddc->onearnullvecs_state[i] != state) { 1183f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1184f4ddd8eeSStefano Zampini break; 1185f4ddd8eeSStefano Zampini } 1186f4ddd8eeSStefano Zampini } 1187f4ddd8eeSStefano Zampini } 1188f4ddd8eeSStefano Zampini } 1189f4ddd8eeSStefano Zampini } else { 1190f4ddd8eeSStefano Zampini if (!nearnullspace) { /* both nearnullspaces are null */ 1191f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 1192f4ddd8eeSStefano Zampini } else { /* nearnullspace attached later */ 1193f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1194f4ddd8eeSStefano Zampini } 1195f4ddd8eeSStefano Zampini } 1196f4ddd8eeSStefano Zampini 1197f4ddd8eeSStefano Zampini /* Setup constraints and related work vectors */ 1198727cdba6SStefano Zampini /* reset primal space flags */ 1199f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1200727cdba6SStefano Zampini pcbddc->new_primal_space_local = PETSC_FALSE; 1201fb8d54d4SStefano Zampini if (computetopography || new_nearnullspace_provided) { 1202727cdba6SStefano Zampini /* It also sets the primal space flags */ 1203674ae819SStefano Zampini ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 1204e7b262bdSStefano Zampini /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1205f4ddd8eeSStefano Zampini ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr); 1206674ae819SStefano Zampini } 1207fb8d54d4SStefano Zampini 1208f4ddd8eeSStefano Zampini if (computesolvers || pcbddc->new_primal_space) { 1209674ae819SStefano Zampini /* reset data */ 1210674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 121199cc7994SStefano Zampini /* Create coarse and local stuffs */ 121299cc7994SStefano Zampini ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 1213674ae819SStefano Zampini ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 12140c7d97c5SJed Brown } 121570cf5478SStefano Zampini 121658a03d70SStefano Zampini if (pcbddc->dbg_flag) { 121758a03d70SStefano Zampini ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 12182b510759SStefano Zampini } 12190c7d97c5SJed Brown PetscFunctionReturn(0); 12200c7d97c5SJed Brown } 12210c7d97c5SJed Brown 12220c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 12230c7d97c5SJed Brown /* 122450efa1b5SStefano Zampini PCApply_BDDC - Applies the BDDC operator to a vector. 12250c7d97c5SJed Brown 12260c7d97c5SJed Brown Input Parameters: 12270c7d97c5SJed Brown . pc - the preconditioner context 12280c7d97c5SJed Brown . r - input vector (global) 12290c7d97c5SJed Brown 12300c7d97c5SJed Brown Output Parameter: 12310c7d97c5SJed Brown . z - output vector (global) 12320c7d97c5SJed Brown 12330c7d97c5SJed Brown Application Interface Routine: PCApply() 12340c7d97c5SJed Brown */ 12350c7d97c5SJed Brown #undef __FUNCT__ 12360c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC" 123753cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 12380c7d97c5SJed Brown { 12390c7d97c5SJed Brown PC_IS *pcis = (PC_IS*)(pc->data); 12400c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 12410c7d97c5SJed Brown PetscErrorCode ierr; 12423b03a366Sstefano_zampini const PetscScalar one = 1.0; 12433b03a366Sstefano_zampini const PetscScalar m_one = -1.0; 12442617d88aSStefano Zampini const PetscScalar zero = 0.0; 12450c7d97c5SJed Brown 12460c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN 12470c7d97c5SJed Brown NN interface preconditioner changed to BDDC 12488eeda7d8SStefano Zampini Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 12490c7d97c5SJed Brown 12500c7d97c5SJed Brown PetscFunctionBegin; 125185c4d303SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 12520c7d97c5SJed Brown /* First Dirichlet solve */ 12530c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12540c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 125553cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 12560c7d97c5SJed Brown /* 12570c7d97c5SJed Brown Assembling right hand side for BDDC operator 1258674ae819SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 1259674ae819SStefano Zampini - pcis->vec1_B the interface part of the global vector z 12600c7d97c5SJed Brown */ 12610c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 12620c7d97c5SJed Brown ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 12638eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 12640c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 12650c7d97c5SJed Brown ierr = VecCopy(r,z);CHKERRQ(ierr); 12660c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12670c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1268674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1269b76ba322SStefano Zampini } else { 12700bdf917eSStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 1271b76ba322SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 1272674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1273b76ba322SStefano Zampini } 1274b76ba322SStefano Zampini 12752617d88aSStefano Zampini /* Apply interface preconditioner 12762617d88aSStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1277dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr); 12782617d88aSStefano Zampini 1279674ae819SStefano Zampini /* Apply transpose of partition of unity operator */ 1280674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 12810c7d97c5SJed Brown 12823b03a366Sstefano_zampini /* Second Dirichlet solve and assembling of output */ 12830c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12840c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12850c7d97c5SJed Brown ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 12868eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1287df187020SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1288df187020SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1289df187020SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1290df187020SStefano Zampini ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 12910c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12920c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12930c7d97c5SJed Brown PetscFunctionReturn(0); 12940c7d97c5SJed Brown } 129550efa1b5SStefano Zampini 129650efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */ 129750efa1b5SStefano Zampini /* 129850efa1b5SStefano Zampini PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector. 129950efa1b5SStefano Zampini 130050efa1b5SStefano Zampini Input Parameters: 130150efa1b5SStefano Zampini . pc - the preconditioner context 130250efa1b5SStefano Zampini . r - input vector (global) 130350efa1b5SStefano Zampini 130450efa1b5SStefano Zampini Output Parameter: 130550efa1b5SStefano Zampini . z - output vector (global) 130650efa1b5SStefano Zampini 130750efa1b5SStefano Zampini Application Interface Routine: PCApplyTranspose() 130850efa1b5SStefano Zampini */ 130950efa1b5SStefano Zampini #undef __FUNCT__ 131050efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC" 131150efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z) 131250efa1b5SStefano Zampini { 131350efa1b5SStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 131450efa1b5SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 131550efa1b5SStefano Zampini PetscErrorCode ierr; 131650efa1b5SStefano Zampini const PetscScalar one = 1.0; 131750efa1b5SStefano Zampini const PetscScalar m_one = -1.0; 131850efa1b5SStefano Zampini const PetscScalar zero = 0.0; 131950efa1b5SStefano Zampini 132050efa1b5SStefano Zampini PetscFunctionBegin; 132150efa1b5SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 132250efa1b5SStefano Zampini /* First Dirichlet solve */ 132350efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 132450efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 132550efa1b5SStefano Zampini ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 132650efa1b5SStefano Zampini /* 132750efa1b5SStefano Zampini Assembling right hand side for BDDC operator 132850efa1b5SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 132950efa1b5SStefano Zampini - pcis->vec1_B the interface part of the global vector z 133050efa1b5SStefano Zampini */ 133150efa1b5SStefano Zampini ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 133250efa1b5SStefano Zampini ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 133350efa1b5SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 133450efa1b5SStefano Zampini ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 133550efa1b5SStefano Zampini ierr = VecCopy(r,z);CHKERRQ(ierr); 133650efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 133750efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 133850efa1b5SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 133950efa1b5SStefano Zampini } else { 134050efa1b5SStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 134150efa1b5SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 134250efa1b5SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 134350efa1b5SStefano Zampini } 134450efa1b5SStefano Zampini 134550efa1b5SStefano Zampini /* Apply interface preconditioner 134650efa1b5SStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1347dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr); 134850efa1b5SStefano Zampini 134950efa1b5SStefano Zampini /* Apply transpose of partition of unity operator */ 135050efa1b5SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 135150efa1b5SStefano Zampini 135250efa1b5SStefano Zampini /* Second Dirichlet solve and assembling of output */ 135350efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 135450efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 135550efa1b5SStefano Zampini ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 135650efa1b5SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1357b0147a47SStefano Zampini ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1358b0147a47SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1359b0147a47SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1360b0147a47SStefano Zampini ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 136150efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 136250efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 136350efa1b5SStefano Zampini PetscFunctionReturn(0); 136450efa1b5SStefano Zampini } 1365da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 1366674ae819SStefano Zampini 1367da1bb401SStefano Zampini #undef __FUNCT__ 1368da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC" 1369da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc) 1370da1bb401SStefano Zampini { 1371da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1372da1bb401SStefano Zampini PetscErrorCode ierr; 1373da1bb401SStefano Zampini 1374da1bb401SStefano Zampini PetscFunctionBegin; 1375da1bb401SStefano Zampini /* free data created by PCIS */ 1376da1bb401SStefano Zampini ierr = PCISDestroy(pc);CHKERRQ(ierr); 1377674ae819SStefano Zampini /* free BDDC custom data */ 1378674ae819SStefano Zampini ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 1379674ae819SStefano Zampini /* destroy objects related to topography */ 1380674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 1381674ae819SStefano Zampini /* free allocated graph structure */ 1382da1bb401SStefano Zampini ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 1383674ae819SStefano Zampini /* free data for scaling operator */ 1384674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 1385674ae819SStefano Zampini /* free solvers stuff */ 1386674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 138762a6ff1dSStefano Zampini /* free global vectors needed in presolve */ 138862a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 138962a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 139039e2fb2aSStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr); 13913425bc38SStefano Zampini /* remove functions */ 1392*b9b85e73SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisLocalMat_C",NULL);CHKERRQ(ierr); 1393674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 1394bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 13952b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 1396b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 13972b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 1398bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 1399bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 140082ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1401bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 140282ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1403bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 140482ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1405bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1406785d1243SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1407bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 140863602bcaSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr); 1409bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 1410bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 1411bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 1412bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 1413674ae819SStefano Zampini /* Free the private data structure */ 1414674ae819SStefano Zampini ierr = PetscFree(pc->data);CHKERRQ(ierr); 1415da1bb401SStefano Zampini PetscFunctionReturn(0); 1416da1bb401SStefano Zampini } 14173425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 14181e6b0712SBarry Smith 14193425bc38SStefano Zampini #undef __FUNCT__ 14203425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 14213425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 14223425bc38SStefano Zampini { 1423674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 14243425bc38SStefano Zampini PC_IS* pcis; 14253425bc38SStefano Zampini PC_BDDC* pcbddc; 14263425bc38SStefano Zampini PetscErrorCode ierr; 14270c7d97c5SJed Brown 14283425bc38SStefano Zampini PetscFunctionBegin; 14293425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 14303425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 14313425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 14323425bc38SStefano Zampini 14333425bc38SStefano Zampini /* change of basis for physical rhs if needed 14343425bc38SStefano Zampini It also changes the rhs in case of dirichlet boundaries */ 14353308cffdSStefano Zampini ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 14363425bc38SStefano Zampini /* store vectors for computation of fetidp final solution */ 14373425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14383425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1439fb223d50SStefano Zampini /* scale rhs since it should be unassembled */ 1440fb223d50SStefano Zampini /* TODO use counter scaling? (also below) */ 14413425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14423425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1443674ae819SStefano Zampini /* Apply partition of unity */ 14443425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1445674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 14468eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 14473425bc38SStefano Zampini /* compute partially subassembled Schur complement right-hand side */ 14483425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 14493425bc38SStefano Zampini ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 14503425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 14513425bc38SStefano Zampini ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 14523425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14533425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1454674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 14553425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14563425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14573425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 14583425bc38SStefano Zampini } 14593425bc38SStefano Zampini /* BDDC rhs */ 14603425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 14618eeda7d8SStefano Zampini if (pcbddc->switch_static) { 14623425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 14633425bc38SStefano Zampini } 14643425bc38SStefano Zampini /* apply BDDC */ 1465dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 14663425bc38SStefano Zampini /* Application of B_delta and assembling of rhs for fetidp fluxes */ 14673425bc38SStefano Zampini ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 14683425bc38SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 14693425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14703425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14713425bc38SStefano Zampini /* restore original rhs */ 14723425bc38SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 14733425bc38SStefano Zampini PetscFunctionReturn(0); 14743425bc38SStefano Zampini } 14751e6b0712SBarry Smith 14763425bc38SStefano Zampini #undef __FUNCT__ 14773425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 14783425bc38SStefano Zampini /*@ 147928509bceSStefano Zampini PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system 14803425bc38SStefano Zampini 14813425bc38SStefano Zampini Collective 14823425bc38SStefano Zampini 14833425bc38SStefano Zampini Input Parameters: 148428509bceSStefano Zampini + fetidp_mat - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators 148528509bceSStefano Zampini . standard_rhs - the right-hand side for your linear system 14863425bc38SStefano Zampini 14873425bc38SStefano Zampini Output Parameters: 148828509bceSStefano Zampini - fetidp_flux_rhs - the right-hand side for the FETIDP linear system 14893425bc38SStefano Zampini 14903425bc38SStefano Zampini Level: developer 14913425bc38SStefano Zampini 14923425bc38SStefano Zampini Notes: 14933425bc38SStefano Zampini 149428509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 14953425bc38SStefano Zampini @*/ 14963425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 14973425bc38SStefano Zampini { 1498674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 14993425bc38SStefano Zampini PetscErrorCode ierr; 15003425bc38SStefano Zampini 15013425bc38SStefano Zampini PetscFunctionBegin; 15023425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15033425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 15043425bc38SStefano Zampini PetscFunctionReturn(0); 15053425bc38SStefano Zampini } 15063425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 15071e6b0712SBarry Smith 15083425bc38SStefano Zampini #undef __FUNCT__ 15093425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 15103425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 15113425bc38SStefano Zampini { 1512674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 15133425bc38SStefano Zampini PC_IS* pcis; 15143425bc38SStefano Zampini PC_BDDC* pcbddc; 15153425bc38SStefano Zampini PetscErrorCode ierr; 15163425bc38SStefano Zampini 15173425bc38SStefano Zampini PetscFunctionBegin; 15183425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15193425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 15203425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 15213425bc38SStefano Zampini 15223425bc38SStefano Zampini /* apply B_delta^T */ 15233425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15243425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15253425bc38SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 15263425bc38SStefano Zampini /* compute rhs for BDDC application */ 15273425bc38SStefano Zampini ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 15288eeda7d8SStefano Zampini if (pcbddc->switch_static) { 15293425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 15303425bc38SStefano Zampini } 15313425bc38SStefano Zampini /* apply BDDC */ 1532dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 15333425bc38SStefano Zampini /* put values into standard global vector */ 15343425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15353425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15368eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 15373425bc38SStefano Zampini /* compute values into the interior if solved for the partially subassembled Schur complement */ 15383425bc38SStefano Zampini ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 15393425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 15403425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 15413425bc38SStefano Zampini } 15423425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15433425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15443425bc38SStefano Zampini /* final change of basis if needed 15453425bc38SStefano Zampini Is also sums the dirichlet part removed during RHS assembling */ 15463308cffdSStefano Zampini ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 15473425bc38SStefano Zampini PetscFunctionReturn(0); 15483425bc38SStefano Zampini } 15491e6b0712SBarry Smith 15503425bc38SStefano Zampini #undef __FUNCT__ 15513425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 15523425bc38SStefano Zampini /*@ 155328509bceSStefano Zampini PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system 15543425bc38SStefano Zampini 15553425bc38SStefano Zampini Collective 15563425bc38SStefano Zampini 15573425bc38SStefano Zampini Input Parameters: 155828509bceSStefano Zampini + fetidp_mat - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators 155928509bceSStefano Zampini . fetidp_flux_sol - the solution of the FETIDP linear system 15603425bc38SStefano Zampini 15613425bc38SStefano Zampini Output Parameters: 156228509bceSStefano Zampini - standard_sol - the solution defined on the physical domain 15633425bc38SStefano Zampini 15643425bc38SStefano Zampini Level: developer 15653425bc38SStefano Zampini 15663425bc38SStefano Zampini Notes: 15673425bc38SStefano Zampini 156828509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 15693425bc38SStefano Zampini @*/ 15703425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 15713425bc38SStefano Zampini { 1572674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 15733425bc38SStefano Zampini PetscErrorCode ierr; 15743425bc38SStefano Zampini 15753425bc38SStefano Zampini PetscFunctionBegin; 15763425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15773425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 15783425bc38SStefano Zampini PetscFunctionReturn(0); 15793425bc38SStefano Zampini } 15803425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 15811e6b0712SBarry Smith 1582f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1583edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec); 1584f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1585f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1586edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec); 1587f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1588674ae819SStefano Zampini 15893425bc38SStefano Zampini #undef __FUNCT__ 15903425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 15913425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 15923425bc38SStefano Zampini { 1593674ae819SStefano Zampini 1594674ae819SStefano Zampini FETIDPMat_ctx fetidpmat_ctx; 15953425bc38SStefano Zampini Mat newmat; 1596674ae819SStefano Zampini FETIDPPC_ctx fetidppc_ctx; 15973425bc38SStefano Zampini PC newpc; 1598ce94432eSBarry Smith MPI_Comm comm; 15993425bc38SStefano Zampini PetscErrorCode ierr; 16003425bc38SStefano Zampini 16013425bc38SStefano Zampini PetscFunctionBegin; 1602ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 16033425bc38SStefano Zampini /* FETIDP linear matrix */ 16043425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 16053425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 16063425bc38SStefano Zampini ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 16073425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 1608edf7251bSStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr); 16093425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 16103425bc38SStefano Zampini ierr = MatSetUp(newmat);CHKERRQ(ierr); 16113425bc38SStefano Zampini /* FETIDP preconditioner */ 16123425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 16133425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 16143425bc38SStefano Zampini ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 16153425bc38SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 16163425bc38SStefano Zampini ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 16173425bc38SStefano Zampini ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 1618edf7251bSStefano Zampini ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr); 16193425bc38SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 162023ee1639SBarry Smith ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr); 16213425bc38SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 16223425bc38SStefano Zampini /* return pointers for objects created */ 16233425bc38SStefano Zampini *fetidp_mat=newmat; 16243425bc38SStefano Zampini *fetidp_pc=newpc; 16253425bc38SStefano Zampini PetscFunctionReturn(0); 16263425bc38SStefano Zampini } 16271e6b0712SBarry Smith 16283425bc38SStefano Zampini #undef __FUNCT__ 16293425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 16303425bc38SStefano Zampini /*@ 163128509bceSStefano Zampini PCBDDCCreateFETIDPOperators - Create operators for FETIDP 16323425bc38SStefano Zampini 16333425bc38SStefano Zampini Collective 16343425bc38SStefano Zampini 16353425bc38SStefano Zampini Input Parameters: 163628509bceSStefano Zampini + pc - the BDDC preconditioning context already setup 163728509bceSStefano Zampini 163828509bceSStefano Zampini Output Parameters: 163928509bceSStefano Zampini . fetidp_mat - shell FETIDP matrix object 164028509bceSStefano Zampini . fetidp_pc - shell Dirichlet preconditioner for FETIDP matrix 164128509bceSStefano Zampini 164228509bceSStefano Zampini Options Database Keys: 164328509bceSStefano Zampini - -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers 16443425bc38SStefano Zampini 16453425bc38SStefano Zampini Level: developer 16463425bc38SStefano Zampini 16473425bc38SStefano Zampini Notes: 164828509bceSStefano Zampini Currently the only operation provided for FETIDP matrix is MatMult 16493425bc38SStefano Zampini 16503425bc38SStefano Zampini .seealso: PCBDDC 16513425bc38SStefano Zampini @*/ 16523425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 16533425bc38SStefano Zampini { 16543425bc38SStefano Zampini PetscErrorCode ierr; 16553425bc38SStefano Zampini 16563425bc38SStefano Zampini PetscFunctionBegin; 16573425bc38SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16583425bc38SStefano Zampini if (pc->setupcalled) { 1659516d51a7SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1660f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 16613425bc38SStefano Zampini PetscFunctionReturn(0); 16623425bc38SStefano Zampini } 16630c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 1664da1bb401SStefano Zampini /*MC 1665da1bb401SStefano Zampini PCBDDC - Balancing Domain Decomposition by Constraints. 16660c7d97c5SJed Brown 166728509bceSStefano Zampini An implementation of the BDDC preconditioner based on 166828509bceSStefano Zampini 166928509bceSStefano Zampini .vb 167028509bceSStefano Zampini [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 167128509bceSStefano 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 167228509bceSStefano Zampini [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 167328509bceSStefano Zampini .ve 167428509bceSStefano Zampini 167528509bceSStefano Zampini The matrix to be preconditioned (Pmat) must be of type MATIS. 167628509bceSStefano Zampini 1677b6fdb6dfSStefano Zampini Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 167828509bceSStefano Zampini 167928509bceSStefano Zampini It also works with unsymmetric and indefinite problems. 168028509bceSStefano Zampini 1681b6fdb6dfSStefano 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. 1682b6fdb6dfSStefano Zampini 168328509bceSStefano Zampini Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace 168428509bceSStefano Zampini 168528509bceSStefano 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 168628509bceSStefano Zampini 168728509bceSStefano Zampini Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace. 168828509bceSStefano Zampini 1689b6fdb6dfSStefano 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. 169028509bceSStefano Zampini 169128509bceSStefano Zampini The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object. 169228509bceSStefano Zampini 1693da1bb401SStefano Zampini Options Database Keys: 169428509bceSStefano Zampini 169528509bceSStefano Zampini . -pc_bddc_use_vertices <1> - use or not vertices in primal space 169628509bceSStefano Zampini . -pc_bddc_use_edges <1> - use or not edges in primal space 1697b6fdb6dfSStefano Zampini . -pc_bddc_use_faces <0> - use or not faces in primal space 169828509bceSStefano Zampini . -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only) 169928509bceSStefano Zampini . -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested 170028509bceSStefano Zampini . -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1]) 170128509bceSStefano Zampini . -pc_bddc_levels <0> - maximum number of levels for multilevel 170228509bceSStefano Zampini . -pc_bddc_coarsening_ratio - H/h ratio at the coarser level 170328509bceSStefano Zampini - -pc_bddc_check_level <0> - set verbosity level of debugging output 170428509bceSStefano Zampini 170528509bceSStefano Zampini Options for Dirichlet, Neumann or coarse solver can be set with 170628509bceSStefano Zampini .vb 170728509bceSStefano Zampini -pc_bddc_dirichlet_ 170828509bceSStefano Zampini -pc_bddc_neumann_ 170928509bceSStefano Zampini -pc_bddc_coarse_ 171028509bceSStefano Zampini .ve 171128509bceSStefano Zampini e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg 171228509bceSStefano Zampini 171328509bceSStefano Zampini When using a multilevel approach, solvers' options at the N-th level can be specified as 171428509bceSStefano Zampini .vb 1715312be037SStefano Zampini -pc_bddc_dirichlet_lN_ 1716312be037SStefano Zampini -pc_bddc_neumann_lN_ 1717312be037SStefano Zampini -pc_bddc_coarse_lN_ 171828509bceSStefano Zampini .ve 1719312be037SStefano Zampini Note that level number ranges from the finest 0 to the coarsest N. 1720da1bb401SStefano Zampini 1721da1bb401SStefano Zampini Level: intermediate 1722da1bb401SStefano Zampini 1723b6fdb6dfSStefano Zampini Developer notes: 1724da1bb401SStefano Zampini 172528509bceSStefano Zampini New deluxe scaling operator will be available soon. 1726da1bb401SStefano Zampini 1727da1bb401SStefano Zampini Contributed by Stefano Zampini 1728da1bb401SStefano Zampini 1729da1bb401SStefano Zampini .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1730da1bb401SStefano Zampini M*/ 1731b2573a8aSBarry Smith 1732da1bb401SStefano Zampini #undef __FUNCT__ 1733da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC" 17348cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1735da1bb401SStefano Zampini { 1736da1bb401SStefano Zampini PetscErrorCode ierr; 1737da1bb401SStefano Zampini PC_BDDC *pcbddc; 1738da1bb401SStefano Zampini 1739da1bb401SStefano Zampini PetscFunctionBegin; 1740da1bb401SStefano Zampini /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1741b00a9115SJed Brown ierr = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr); 1742da1bb401SStefano Zampini pc->data = (void*)pcbddc; 1743da1bb401SStefano Zampini 1744da1bb401SStefano Zampini /* create PCIS data structure */ 1745da1bb401SStefano Zampini ierr = PCISCreate(pc);CHKERRQ(ierr); 1746da1bb401SStefano Zampini 174747d04d0dSStefano Zampini /* BDDC customization */ 174808a5cf49SStefano Zampini pcbddc->use_local_adj = PETSC_TRUE; 174947d04d0dSStefano Zampini pcbddc->use_vertices = PETSC_TRUE; 175047d04d0dSStefano Zampini pcbddc->use_edges = PETSC_TRUE; 175147d04d0dSStefano Zampini pcbddc->use_faces = PETSC_FALSE; 175247d04d0dSStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 175347d04d0dSStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 175447d04d0dSStefano Zampini pcbddc->switch_static = PETSC_FALSE; 175547d04d0dSStefano Zampini pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 175647d04d0dSStefano Zampini pcbddc->dbg_flag = 0; 175747d04d0dSStefano Zampini 175839e2fb2aSStefano Zampini pcbddc->BtoNmap = 0; 1759727cdba6SStefano Zampini pcbddc->local_primal_size = 0; 1760e9189074SStefano Zampini pcbddc->n_vertices = 0; 1761e9189074SStefano Zampini pcbddc->n_actual_vertices = 0; 1762e9189074SStefano Zampini pcbddc->n_constraints = 0; 1763727cdba6SStefano Zampini pcbddc->primal_indices_local_idxs = 0; 1764fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_FALSE; 176568457ee5SStefano Zampini pcbddc->coarse_size = -1; 1766f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1767727cdba6SStefano Zampini pcbddc->new_primal_space_local = PETSC_FALSE; 1768f4ddd8eeSStefano Zampini pcbddc->global_primal_indices = 0; 1769f4ddd8eeSStefano Zampini pcbddc->onearnullspace = 0; 1770f4ddd8eeSStefano Zampini pcbddc->onearnullvecs_state = 0; 1771674ae819SStefano Zampini pcbddc->user_primal_vertices = 0; 17720bdf917eSStefano Zampini pcbddc->NullSpace = 0; 17733972b0daSStefano Zampini pcbddc->temp_solution = 0; 1774534831adSStefano Zampini pcbddc->original_rhs = 0; 1775534831adSStefano Zampini pcbddc->local_mat = 0; 1776534831adSStefano Zampini pcbddc->ChangeOfBasisMatrix = 0; 1777*b9b85e73SStefano Zampini pcbddc->user_ChangeOfBasisMatrix = 0; 1778da1bb401SStefano Zampini pcbddc->coarse_vec = 0; 1779da1bb401SStefano Zampini pcbddc->coarse_rhs = 0; 1780da1bb401SStefano Zampini pcbddc->coarse_ksp = 0; 1781da1bb401SStefano Zampini pcbddc->coarse_phi_B = 0; 1782da1bb401SStefano Zampini pcbddc->coarse_phi_D = 0; 178315aaf578SStefano Zampini pcbddc->coarse_psi_B = 0; 178415aaf578SStefano Zampini pcbddc->coarse_psi_D = 0; 1785da1bb401SStefano Zampini pcbddc->vec1_P = 0; 1786da1bb401SStefano Zampini pcbddc->vec1_R = 0; 1787da1bb401SStefano Zampini pcbddc->vec2_R = 0; 1788da1bb401SStefano Zampini pcbddc->local_auxmat1 = 0; 1789da1bb401SStefano Zampini pcbddc->local_auxmat2 = 0; 1790da1bb401SStefano Zampini pcbddc->R_to_B = 0; 1791da1bb401SStefano Zampini pcbddc->R_to_D = 0; 1792da1bb401SStefano Zampini pcbddc->ksp_D = 0; 1793da1bb401SStefano Zampini pcbddc->ksp_R = 0; 1794da1bb401SStefano Zampini pcbddc->NeumannBoundaries = 0; 1795785d1243SStefano Zampini pcbddc->NeumannBoundariesLocal = 0; 1796785d1243SStefano Zampini pcbddc->DirichletBoundaries = 0; 1797785d1243SStefano Zampini pcbddc->DirichletBoundariesLocal = 0; 179860d44989SStefano Zampini pcbddc->user_provided_isfordofs = PETSC_FALSE; 179960d44989SStefano Zampini pcbddc->n_ISForDofs = 0; 180063602bcaSStefano Zampini pcbddc->n_ISForDofsLocal = 0; 1801da1bb401SStefano Zampini pcbddc->ISForDofs = 0; 180263602bcaSStefano Zampini pcbddc->ISForDofsLocal = 0; 1803da1bb401SStefano Zampini pcbddc->ConstraintMatrix = 0; 180485c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 180547d04d0dSStefano Zampini pcbddc->coarse_loc_to_glob = 0; 180647d04d0dSStefano Zampini pcbddc->coarsening_ratio = 8; 18074fad6a16SStefano Zampini pcbddc->current_level = 0; 18082b510759SStefano Zampini pcbddc->max_levels = 0; 1809323d395dSStefano Zampini pcbddc->use_coarse_estimates = PETSC_FALSE; 1810f3bde8b3SStefano Zampini pcbddc->coarse_subassembling = 0; 1811323d395dSStefano Zampini pcbddc->coarse_subassembling_init = 0; 1812da1bb401SStefano Zampini 1813674ae819SStefano Zampini /* create local graph structure */ 1814674ae819SStefano Zampini ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1815674ae819SStefano Zampini 1816674ae819SStefano Zampini /* scaling */ 1817674ae819SStefano Zampini pcbddc->use_deluxe_scaling = PETSC_FALSE; 181808a5cf49SStefano Zampini pcbddc->deluxe_ctx = 0; 1819674ae819SStefano Zampini pcbddc->work_scaling = 0; 1820da1bb401SStefano Zampini 1821da1bb401SStefano Zampini /* function pointers */ 1822da1bb401SStefano Zampini pc->ops->apply = PCApply_BDDC; 182393bd9ae7SStefano Zampini pc->ops->applytranspose = PCApplyTranspose_BDDC; 1824da1bb401SStefano Zampini pc->ops->setup = PCSetUp_BDDC; 1825da1bb401SStefano Zampini pc->ops->destroy = PCDestroy_BDDC; 1826da1bb401SStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1827da1bb401SStefano Zampini pc->ops->view = 0; 1828da1bb401SStefano Zampini pc->ops->applyrichardson = 0; 1829da1bb401SStefano Zampini pc->ops->applysymmetricleft = 0; 1830da1bb401SStefano Zampini pc->ops->applysymmetricright = 0; 1831534831adSStefano Zampini pc->ops->presolve = PCPreSolve_BDDC; 1832534831adSStefano Zampini pc->ops->postsolve = PCPostSolve_BDDC; 1833da1bb401SStefano Zampini 1834da1bb401SStefano Zampini /* composing function */ 1835*b9b85e73SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisLocalMat_C",PCBDDCSetChangeOfBasisLocalMat_BDDC);CHKERRQ(ierr); 1836674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1837bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 18382b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1839b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 18402b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1841bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1842bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 184382ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1844bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 184582ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1846bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 184782ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1848bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 184982ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1850bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 185163602bcaSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr); 1852bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1853bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1854bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1855bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1856da1bb401SStefano Zampini PetscFunctionReturn(0); 1857da1bb401SStefano Zampini } 18583425bc38SStefano Zampini 1859