153cdbc3dSStefano Zampini /* TODOLIST 2eb97c9d2SStefano Zampini 3eb97c9d2SStefano Zampini ConstraintsSetup 4eb97c9d2SStefano Zampini - tolerances for constraints as an option (take care of single precision!) 5eb97c9d2SStefano Zampini 6eb97c9d2SStefano Zampini Solvers 7a0d3c3abSStefano Zampini - Add support for cholesky for coarse solver (similar to local solvers) 8eb97c9d2SStefano Zampini - Propagate ksp prefixes for solvers to mat objects? 9eb97c9d2SStefano Zampini - Propagate nearnullspace info among levels 10eb97c9d2SStefano Zampini 11eb97c9d2SStefano Zampini User interface 12b9b85e73SStefano Zampini - ** DofSplitting and DM attached to pc? 13eb97c9d2SStefano Zampini 14eb97c9d2SStefano Zampini Debugging output 15b9b85e73SStefano Zampini - * Better management of verbosity levels of debugging output 16eb97c9d2SStefano Zampini 17eb97c9d2SStefano Zampini Build 18eb97c9d2SStefano Zampini - make runexe59 19eb97c9d2SStefano Zampini 20eb97c9d2SStefano Zampini Extra 21b9b85e73SStefano Zampini - ** GetRid of PCBDDCApplySchur, use MatSchur instead 22b9b85e73SStefano Zampini - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)? 23eb97c9d2SStefano Zampini - add support for computing h,H and related using coordinates? 24c0b83709SStefano Zampini - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap) 25eb97c9d2SStefano Zampini - Better management in PCIS code 26eb97c9d2SStefano Zampini - BDDC with MG framework? 27eb97c9d2SStefano Zampini 28eb97c9d2SStefano Zampini FETIDP 29eb97c9d2SStefano Zampini - Move FETIDP code to its own classes 30eb97c9d2SStefano Zampini 31eb97c9d2SStefano Zampini MATIS related operations contained in BDDC code 32eb97c9d2SStefano Zampini - Provide general case for subassembling 33b9b85e73SStefano Zampini - *** Preallocation routines in MatISGetMPIAXAIJ 34eb97c9d2SStefano Zampini 3553cdbc3dSStefano Zampini */ 360c7d97c5SJed Brown 3753cdbc3dSStefano Zampini /* ---------------------------------------------------------------------------------------------------------------------------------------------- 380c7d97c5SJed Brown Implementation of BDDC preconditioner based on: 390c7d97c5SJed Brown C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 4053cdbc3dSStefano Zampini ---------------------------------------------------------------------------------------------------------------------------------------------- */ 4153cdbc3dSStefano Zampini 42ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 43ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 443b03a366Sstefano_zampini #include <petscblaslapack.h> 45674ae819SStefano Zampini 460c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 470c7d97c5SJed Brown #undef __FUNCT__ 480c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC" 490c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc) 500c7d97c5SJed Brown { 510c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 520c7d97c5SJed Brown PetscErrorCode ierr; 530c7d97c5SJed Brown 540c7d97c5SJed Brown PetscFunctionBegin; 550c7d97c5SJed Brown ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 568eeda7d8SStefano Zampini /* Verbose debugging */ 578eeda7d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 588eeda7d8SStefano Zampini /* Primal space cumstomization */ 5908a5cf49SStefano 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); 608eeda7d8SStefano 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); 618eeda7d8SStefano 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); 628eeda7d8SStefano 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); 638eeda7d8SStefano Zampini /* Change of basis */ 64b9b85e73SStefano 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); 65b9b85e73SStefano 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); 66674ae819SStefano Zampini if (!pcbddc->use_change_of_basis) { 67674ae819SStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 68674ae819SStefano Zampini } 698eeda7d8SStefano Zampini /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 708eeda7d8SStefano 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); 710298fd71SBarry Smith ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 722b510759SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 73323d395dSStefano 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); 74674ae819SStefano 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); 750c7d97c5SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 760c7d97c5SJed Brown PetscFunctionReturn(0); 770c7d97c5SJed Brown } 780c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 79674ae819SStefano Zampini #undef __FUNCT__ 80b9b85e73SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisLocalMat_BDDC" 81b9b85e73SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisLocalMat_BDDC(PC pc, Mat change) 82b9b85e73SStefano Zampini { 83b9b85e73SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 84b9b85e73SStefano Zampini PetscErrorCode ierr; 85b9b85e73SStefano Zampini 86b9b85e73SStefano Zampini PetscFunctionBegin; 87b9b85e73SStefano Zampini ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr); 88b9b85e73SStefano Zampini ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr); 89b9b85e73SStefano Zampini pcbddc->user_ChangeOfBasisMatrix = change; 90b9b85e73SStefano Zampini PetscFunctionReturn(0); 91b9b85e73SStefano Zampini } 92b9b85e73SStefano Zampini #undef __FUNCT__ 93b9b85e73SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisLocalMat" 94b9b85e73SStefano Zampini /*@ 95b9b85e73SStefano Zampini PCBDDCSetChangeOfBasisLocalMat - Set user defined change of basis for local boundary dofs 96b9b85e73SStefano Zampini 97b9b85e73SStefano Zampini Collective on PC 98b9b85e73SStefano Zampini 99b9b85e73SStefano Zampini Input Parameters: 100b9b85e73SStefano Zampini + pc - the preconditioning context 101b9b85e73SStefano Zampini - change - the local change of basis matrix, either in local (internal + boundary) or in local boundary numbering 102b9b85e73SStefano Zampini 103b9b85e73SStefano Zampini Level: intermediate 104b9b85e73SStefano Zampini 105b9b85e73SStefano Zampini Notes: 106b9b85e73SStefano Zampini 107b9b85e73SStefano Zampini .seealso: PCBDDC 108b9b85e73SStefano Zampini @*/ 109b9b85e73SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisLocalMat(PC pc, Mat change) 110b9b85e73SStefano Zampini { 111b9b85e73SStefano Zampini PetscErrorCode ierr; 112b9b85e73SStefano Zampini 113b9b85e73SStefano Zampini PetscFunctionBegin; 114b9b85e73SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 115b9b85e73SStefano Zampini PetscValidHeaderSpecific(change,MAT_CLASSID,2); 116b9b85e73SStefano Zampini PetscValidLogicalCollectiveBool(pc,change,2); 117b9b85e73SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisLocalMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr); 118b9b85e73SStefano Zampini PetscFunctionReturn(0); 119b9b85e73SStefano Zampini } 120b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */ 121b9b85e73SStefano Zampini #undef __FUNCT__ 122674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 123674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 124674ae819SStefano Zampini { 125674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 126674ae819SStefano Zampini PetscErrorCode ierr; 1271e6b0712SBarry Smith 128674ae819SStefano Zampini PetscFunctionBegin; 129674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 130674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 131674ae819SStefano Zampini pcbddc->user_primal_vertices = PrimalVertices; 132674ae819SStefano Zampini PetscFunctionReturn(0); 133674ae819SStefano Zampini } 134674ae819SStefano Zampini #undef __FUNCT__ 135674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 136674ae819SStefano Zampini /*@ 13728509bceSStefano Zampini PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 138674ae819SStefano Zampini 139674ae819SStefano Zampini Not collective 140674ae819SStefano Zampini 141674ae819SStefano Zampini Input Parameters: 142674ae819SStefano Zampini + pc - the preconditioning context 14328509bceSStefano Zampini - PrimalVertices - index set of primal vertices in local numbering 144674ae819SStefano Zampini 145674ae819SStefano Zampini Level: intermediate 146674ae819SStefano Zampini 147674ae819SStefano Zampini Notes: 148674ae819SStefano Zampini 149674ae819SStefano Zampini .seealso: PCBDDC 150674ae819SStefano Zampini @*/ 151674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 152674ae819SStefano Zampini { 153674ae819SStefano Zampini PetscErrorCode ierr; 154674ae819SStefano Zampini 155674ae819SStefano Zampini PetscFunctionBegin; 156674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 157674ae819SStefano Zampini PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 158674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 159674ae819SStefano Zampini PetscFunctionReturn(0); 160674ae819SStefano Zampini } 161674ae819SStefano Zampini /* -------------------------------------------------------------------------- */ 1620c7d97c5SJed Brown #undef __FUNCT__ 1634fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 1644fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 1654fad6a16SStefano Zampini { 1664fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1674fad6a16SStefano Zampini 1684fad6a16SStefano Zampini PetscFunctionBegin; 1694fad6a16SStefano Zampini pcbddc->coarsening_ratio = k; 1704fad6a16SStefano Zampini PetscFunctionReturn(0); 1714fad6a16SStefano Zampini } 1721e6b0712SBarry Smith 1734fad6a16SStefano Zampini #undef __FUNCT__ 1744fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio" 1754fad6a16SStefano Zampini /*@ 17628509bceSStefano Zampini PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 1774fad6a16SStefano Zampini 1784fad6a16SStefano Zampini Logically collective on PC 1794fad6a16SStefano Zampini 1804fad6a16SStefano Zampini Input Parameters: 1814fad6a16SStefano Zampini + pc - the preconditioning context 18228509bceSStefano Zampini - k - coarsening ratio (H/h at the coarser level) 1834fad6a16SStefano Zampini 18428509bceSStefano Zampini Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 1854fad6a16SStefano Zampini 1864fad6a16SStefano Zampini Level: intermediate 1874fad6a16SStefano Zampini 1884fad6a16SStefano Zampini Notes: 1894fad6a16SStefano Zampini 1904fad6a16SStefano Zampini .seealso: PCBDDC 1914fad6a16SStefano Zampini @*/ 1924fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 1934fad6a16SStefano Zampini { 1944fad6a16SStefano Zampini PetscErrorCode ierr; 1954fad6a16SStefano Zampini 1964fad6a16SStefano Zampini PetscFunctionBegin; 1974fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1982b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,k,2); 1994fad6a16SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 2004fad6a16SStefano Zampini PetscFunctionReturn(0); 2014fad6a16SStefano Zampini } 2022b510759SStefano Zampini 203b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 2042b510759SStefano Zampini #undef __FUNCT__ 205b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 206b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 207b8ffe317SStefano Zampini { 208b8ffe317SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 209b8ffe317SStefano Zampini 210b8ffe317SStefano Zampini PetscFunctionBegin; 21185c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = flg; 212b8ffe317SStefano Zampini PetscFunctionReturn(0); 213b8ffe317SStefano Zampini } 214b8ffe317SStefano Zampini 215b8ffe317SStefano Zampini #undef __FUNCT__ 216b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 217b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 2182b510759SStefano Zampini { 2192b510759SStefano Zampini PetscErrorCode ierr; 2202b510759SStefano Zampini 2212b510759SStefano Zampini PetscFunctionBegin; 2222b510759SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 223b8ffe317SStefano Zampini PetscValidLogicalCollectiveBool(pc,flg,2); 224b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 2252b510759SStefano Zampini PetscFunctionReturn(0); 2262b510759SStefano Zampini } 2271e6b0712SBarry Smith 2284fad6a16SStefano Zampini #undef __FUNCT__ 2292b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC" 2302b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 2314fad6a16SStefano Zampini { 2324fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2334fad6a16SStefano Zampini 2344fad6a16SStefano Zampini PetscFunctionBegin; 2352b510759SStefano Zampini pcbddc->current_level = level; 2364fad6a16SStefano Zampini PetscFunctionReturn(0); 2374fad6a16SStefano Zampini } 2381e6b0712SBarry Smith 2394fad6a16SStefano Zampini #undef __FUNCT__ 240b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel" 241b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 242b8ffe317SStefano Zampini { 243b8ffe317SStefano Zampini PetscErrorCode ierr; 244b8ffe317SStefano Zampini 245b8ffe317SStefano Zampini PetscFunctionBegin; 246b8ffe317SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 247b8ffe317SStefano Zampini PetscValidLogicalCollectiveInt(pc,level,2); 248b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 249b8ffe317SStefano Zampini PetscFunctionReturn(0); 250b8ffe317SStefano Zampini } 251b8ffe317SStefano Zampini 252b8ffe317SStefano Zampini #undef __FUNCT__ 2532b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC" 2542b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 2552b510759SStefano Zampini { 2562b510759SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2572b510759SStefano Zampini 2582b510759SStefano Zampini PetscFunctionBegin; 2592b510759SStefano Zampini pcbddc->max_levels = levels; 2602b510759SStefano Zampini PetscFunctionReturn(0); 2612b510759SStefano Zampini } 2622b510759SStefano Zampini 2632b510759SStefano Zampini #undef __FUNCT__ 2642b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels" 2654fad6a16SStefano Zampini /*@ 26628509bceSStefano Zampini PCBDDCSetLevels - Sets the maximum number of levels for multilevel 2674fad6a16SStefano Zampini 2684fad6a16SStefano Zampini Logically collective on PC 2694fad6a16SStefano Zampini 2704fad6a16SStefano Zampini Input Parameters: 2714fad6a16SStefano Zampini + pc - the preconditioning context 27228509bceSStefano Zampini - levels - the maximum number of levels (max 9) 2734fad6a16SStefano Zampini 27428509bceSStefano Zampini Default value is 0, i.e. traditional one-level BDDC 2754fad6a16SStefano Zampini 2764fad6a16SStefano Zampini Level: intermediate 2774fad6a16SStefano Zampini 2784fad6a16SStefano Zampini Notes: 2794fad6a16SStefano Zampini 2804fad6a16SStefano Zampini .seealso: PCBDDC 2814fad6a16SStefano Zampini @*/ 2822b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 2834fad6a16SStefano Zampini { 2844fad6a16SStefano Zampini PetscErrorCode ierr; 2854fad6a16SStefano Zampini 2864fad6a16SStefano Zampini PetscFunctionBegin; 2874fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2882b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,levels,2); 289312be037SStefano Zampini if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n"); 2902b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 2914fad6a16SStefano Zampini PetscFunctionReturn(0); 2924fad6a16SStefano Zampini } 2934fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */ 2941e6b0712SBarry Smith 2954fad6a16SStefano Zampini #undef __FUNCT__ 2960bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 2970bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 2980bdf917eSStefano Zampini { 2990bdf917eSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3000bdf917eSStefano Zampini PetscErrorCode ierr; 3010bdf917eSStefano Zampini 3020bdf917eSStefano Zampini PetscFunctionBegin; 3030bdf917eSStefano Zampini ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 3040bdf917eSStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 3050bdf917eSStefano Zampini pcbddc->NullSpace = NullSpace; 3060bdf917eSStefano Zampini PetscFunctionReturn(0); 3070bdf917eSStefano Zampini } 3081e6b0712SBarry Smith 3090bdf917eSStefano Zampini #undef __FUNCT__ 3100bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace" 3110bdf917eSStefano Zampini /*@ 31228509bceSStefano Zampini PCBDDCSetNullSpace - Set nullspace for BDDC operator 3130bdf917eSStefano Zampini 3140bdf917eSStefano Zampini Logically collective on PC and MatNullSpace 3150bdf917eSStefano Zampini 3160bdf917eSStefano Zampini Input Parameters: 3170bdf917eSStefano Zampini + pc - the preconditioning context 31828509bceSStefano Zampini - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 3190bdf917eSStefano Zampini 3200bdf917eSStefano Zampini Level: intermediate 3210bdf917eSStefano Zampini 3220bdf917eSStefano Zampini Notes: 3230bdf917eSStefano Zampini 3240bdf917eSStefano Zampini .seealso: PCBDDC 3250bdf917eSStefano Zampini @*/ 3260bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 3270bdf917eSStefano Zampini { 3280bdf917eSStefano Zampini PetscErrorCode ierr; 3290bdf917eSStefano Zampini 3300bdf917eSStefano Zampini PetscFunctionBegin; 3310bdf917eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 332674ae819SStefano Zampini PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 3332b510759SStefano Zampini PetscCheckSameComm(pc,1,NullSpace,2); 3340bdf917eSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 3350bdf917eSStefano Zampini PetscFunctionReturn(0); 3360bdf917eSStefano Zampini } 3370bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */ 3381e6b0712SBarry Smith 3390bdf917eSStefano Zampini #undef __FUNCT__ 3403b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 3413b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 3423b03a366Sstefano_zampini { 3433b03a366Sstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3443b03a366Sstefano_zampini PetscErrorCode ierr; 3453b03a366Sstefano_zampini 3463b03a366Sstefano_zampini PetscFunctionBegin; 347785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 348785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 3493b03a366Sstefano_zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 35036e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 35136e030ebSStefano Zampini pcbddc->DirichletBoundaries = DirichletBoundaries; 352fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 3533b03a366Sstefano_zampini PetscFunctionReturn(0); 3543b03a366Sstefano_zampini } 3551e6b0712SBarry Smith 3563b03a366Sstefano_zampini #undef __FUNCT__ 3573b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 3583b03a366Sstefano_zampini /*@ 35928509bceSStefano Zampini PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 3603b03a366Sstefano_zampini 361785d1243SStefano Zampini Collective 3623b03a366Sstefano_zampini 3633b03a366Sstefano_zampini Input Parameters: 3643b03a366Sstefano_zampini + pc - the preconditioning context 365785d1243SStefano Zampini - DirichletBoundaries - parallel IS defining the Dirichlet boundaries 3663b03a366Sstefano_zampini 3673b03a366Sstefano_zampini Level: intermediate 3683b03a366Sstefano_zampini 369785d1243SStefano Zampini Notes: Any process can list any global node 3703b03a366Sstefano_zampini 3713b03a366Sstefano_zampini .seealso: PCBDDC 3723b03a366Sstefano_zampini @*/ 3733b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 3743b03a366Sstefano_zampini { 3753b03a366Sstefano_zampini PetscErrorCode ierr; 3763b03a366Sstefano_zampini 3773b03a366Sstefano_zampini PetscFunctionBegin; 3783b03a366Sstefano_zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 379674ae819SStefano Zampini PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 380785d1243SStefano Zampini PetscCheckSameComm(pc,1,DirichletBoundaries,2); 3813b03a366Sstefano_zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 3823b03a366Sstefano_zampini PetscFunctionReturn(0); 3833b03a366Sstefano_zampini } 3843b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */ 3851e6b0712SBarry Smith 3863b03a366Sstefano_zampini #undef __FUNCT__ 38782ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC" 38882ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries) 3890c7d97c5SJed Brown { 3900c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3910c7d97c5SJed Brown PetscErrorCode ierr; 3920c7d97c5SJed Brown 3930c7d97c5SJed Brown PetscFunctionBegin; 394785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 395785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 3960c7d97c5SJed Brown ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 3970c7d97c5SJed Brown ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 398785d1243SStefano Zampini pcbddc->DirichletBoundariesLocal = DirichletBoundaries; 3997d24d155SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4000c7d97c5SJed Brown PetscFunctionReturn(0); 4010c7d97c5SJed Brown } 4020c7d97c5SJed Brown 4030c7d97c5SJed Brown #undef __FUNCT__ 40482ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal" 4059c0446d6SStefano Zampini /*@ 40682ba6b80SStefano Zampini PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering. 4079c0446d6SStefano Zampini 408785d1243SStefano Zampini Collective 40953cdbc3dSStefano Zampini 41053cdbc3dSStefano Zampini Input Parameters: 41153cdbc3dSStefano Zampini + pc - the preconditioning context 41282ba6b80SStefano Zampini - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering) 41353cdbc3dSStefano Zampini 41453cdbc3dSStefano Zampini Level: intermediate 41553cdbc3dSStefano Zampini 4169c0446d6SStefano Zampini Notes: 41753cdbc3dSStefano Zampini 41853cdbc3dSStefano Zampini .seealso: PCBDDC 41953cdbc3dSStefano Zampini @*/ 42082ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries) 4210c7d97c5SJed Brown { 4220c7d97c5SJed Brown PetscErrorCode ierr; 4230c7d97c5SJed Brown 4240c7d97c5SJed Brown PetscFunctionBegin; 4250c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4260c7d97c5SJed Brown PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 42782ba6b80SStefano Zampini PetscCheckSameComm(pc,1,DirichletBoundaries,2); 42882ba6b80SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 4290c7d97c5SJed Brown PetscFunctionReturn(0); 4300c7d97c5SJed Brown } 4310c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 4320c7d97c5SJed Brown 4330c7d97c5SJed Brown #undef __FUNCT__ 4340c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 43553cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 4360c7d97c5SJed Brown { 4370c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 43853cdbc3dSStefano Zampini PetscErrorCode ierr; 4390c7d97c5SJed Brown 4400c7d97c5SJed Brown PetscFunctionBegin; 441785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 442785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 44353cdbc3dSStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 44436e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 44536e030ebSStefano Zampini pcbddc->NeumannBoundaries = NeumannBoundaries; 446fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4470c7d97c5SJed Brown PetscFunctionReturn(0); 4480c7d97c5SJed Brown } 4491e6b0712SBarry Smith 4500c7d97c5SJed Brown #undef __FUNCT__ 4510c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 45257527edcSJed Brown /*@ 45328509bceSStefano Zampini PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 45457527edcSJed Brown 455785d1243SStefano Zampini Collective 45657527edcSJed Brown 45757527edcSJed Brown Input Parameters: 45857527edcSJed Brown + pc - the preconditioning context 459785d1243SStefano Zampini - NeumannBoundaries - parallel IS defining the Neumann boundaries 46057527edcSJed Brown 46157527edcSJed Brown Level: intermediate 46257527edcSJed Brown 463785d1243SStefano Zampini Notes: Any process can list any global node 46457527edcSJed Brown 46557527edcSJed Brown .seealso: PCBDDC 46657527edcSJed Brown @*/ 46753cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 4680c7d97c5SJed Brown { 4690c7d97c5SJed Brown PetscErrorCode ierr; 4700c7d97c5SJed Brown 4710c7d97c5SJed Brown PetscFunctionBegin; 4720c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 473674ae819SStefano Zampini PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 474785d1243SStefano Zampini PetscCheckSameComm(pc,1,NeumannBoundaries,2); 47553cdbc3dSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 47653cdbc3dSStefano Zampini PetscFunctionReturn(0); 47753cdbc3dSStefano Zampini } 47853cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 4791e6b0712SBarry Smith 48053cdbc3dSStefano Zampini #undef __FUNCT__ 48182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC" 48282ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries) 4830c7d97c5SJed Brown { 4840c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 4850c7d97c5SJed Brown PetscErrorCode ierr; 4860c7d97c5SJed Brown 4870c7d97c5SJed Brown PetscFunctionBegin; 488785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 489785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 4900c7d97c5SJed Brown ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 4910c7d97c5SJed Brown ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 492785d1243SStefano Zampini pcbddc->NeumannBoundariesLocal = NeumannBoundaries; 4937d24d155SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4940c7d97c5SJed Brown PetscFunctionReturn(0); 4950c7d97c5SJed Brown } 4960c7d97c5SJed Brown 4970c7d97c5SJed Brown #undef __FUNCT__ 49882ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal" 4990c7d97c5SJed Brown /*@ 50082ba6b80SStefano Zampini PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering. 5010c7d97c5SJed Brown 502785d1243SStefano Zampini Collective 5030c7d97c5SJed Brown 5040c7d97c5SJed Brown Input Parameters: 5050c7d97c5SJed Brown + pc - the preconditioning context 50682ba6b80SStefano Zampini - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering) 5070c7d97c5SJed Brown 5080c7d97c5SJed Brown Level: intermediate 5090c7d97c5SJed Brown 5100c7d97c5SJed Brown Notes: 5110c7d97c5SJed Brown 5120c7d97c5SJed Brown .seealso: PCBDDC 5130c7d97c5SJed Brown @*/ 51482ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries) 5150c7d97c5SJed Brown { 5160c7d97c5SJed Brown PetscErrorCode ierr; 5170c7d97c5SJed Brown 5180c7d97c5SJed Brown PetscFunctionBegin; 5190c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5200c7d97c5SJed Brown PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 52182ba6b80SStefano Zampini PetscCheckSameComm(pc,1,NeumannBoundaries,2); 52282ba6b80SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 52353cdbc3dSStefano Zampini PetscFunctionReturn(0); 52453cdbc3dSStefano Zampini } 52553cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 52653cdbc3dSStefano Zampini 52753cdbc3dSStefano Zampini #undef __FUNCT__ 528da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 529da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 530da1bb401SStefano Zampini { 531da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 532da1bb401SStefano Zampini 533da1bb401SStefano Zampini PetscFunctionBegin; 534da1bb401SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundaries; 535da1bb401SStefano Zampini PetscFunctionReturn(0); 536da1bb401SStefano Zampini } 5371e6b0712SBarry Smith 538da1bb401SStefano Zampini #undef __FUNCT__ 539da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 540da1bb401SStefano Zampini /*@ 541785d1243SStefano Zampini PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries 542da1bb401SStefano Zampini 543785d1243SStefano Zampini Collective 544785d1243SStefano Zampini 545785d1243SStefano Zampini Input Parameters: 546785d1243SStefano Zampini . pc - the preconditioning context 547785d1243SStefano Zampini 548785d1243SStefano Zampini Output Parameters: 549785d1243SStefano Zampini . DirichletBoundaries - index set defining the Dirichlet boundaries 550785d1243SStefano Zampini 551785d1243SStefano Zampini Level: intermediate 552785d1243SStefano Zampini 553785d1243SStefano Zampini Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries 554785d1243SStefano Zampini 555785d1243SStefano Zampini .seealso: PCBDDC 556785d1243SStefano Zampini @*/ 557785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 558785d1243SStefano Zampini { 559785d1243SStefano Zampini PetscErrorCode ierr; 560785d1243SStefano Zampini 561785d1243SStefano Zampini PetscFunctionBegin; 562785d1243SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 563785d1243SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 564785d1243SStefano Zampini PetscFunctionReturn(0); 565785d1243SStefano Zampini } 566785d1243SStefano Zampini /* -------------------------------------------------------------------------- */ 567785d1243SStefano Zampini 568785d1243SStefano Zampini #undef __FUNCT__ 569785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC" 570785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries) 571785d1243SStefano Zampini { 572785d1243SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 573785d1243SStefano Zampini 574785d1243SStefano Zampini PetscFunctionBegin; 575785d1243SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundariesLocal; 576785d1243SStefano Zampini PetscFunctionReturn(0); 577785d1243SStefano Zampini } 578785d1243SStefano Zampini 579785d1243SStefano Zampini #undef __FUNCT__ 58082ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal" 581da1bb401SStefano Zampini /*@ 58282ba6b80SStefano Zampini PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering) 583da1bb401SStefano Zampini 584785d1243SStefano Zampini Collective 585da1bb401SStefano Zampini 586da1bb401SStefano Zampini Input Parameters: 58728509bceSStefano Zampini . pc - the preconditioning context 588da1bb401SStefano Zampini 589da1bb401SStefano Zampini Output Parameters: 59028509bceSStefano Zampini . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 591da1bb401SStefano Zampini 592da1bb401SStefano Zampini Level: intermediate 593da1bb401SStefano Zampini 594da1bb401SStefano Zampini Notes: 595da1bb401SStefano Zampini 596da1bb401SStefano Zampini .seealso: PCBDDC 597da1bb401SStefano Zampini @*/ 59882ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries) 599da1bb401SStefano Zampini { 600da1bb401SStefano Zampini PetscErrorCode ierr; 601da1bb401SStefano Zampini 602da1bb401SStefano Zampini PetscFunctionBegin; 603da1bb401SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 60482ba6b80SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 605da1bb401SStefano Zampini PetscFunctionReturn(0); 606da1bb401SStefano Zampini } 607da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 6081e6b0712SBarry Smith 609da1bb401SStefano Zampini #undef __FUNCT__ 61053cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 61153cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 61253cdbc3dSStefano Zampini { 61353cdbc3dSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 61453cdbc3dSStefano Zampini 61553cdbc3dSStefano Zampini PetscFunctionBegin; 61653cdbc3dSStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundaries; 61753cdbc3dSStefano Zampini PetscFunctionReturn(0); 61853cdbc3dSStefano Zampini } 6191e6b0712SBarry Smith 62053cdbc3dSStefano Zampini #undef __FUNCT__ 62153cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 62253cdbc3dSStefano Zampini /*@ 623785d1243SStefano Zampini PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries 62453cdbc3dSStefano Zampini 625785d1243SStefano Zampini Collective 626785d1243SStefano Zampini 627785d1243SStefano Zampini Input Parameters: 628785d1243SStefano Zampini . pc - the preconditioning context 629785d1243SStefano Zampini 630785d1243SStefano Zampini Output Parameters: 631785d1243SStefano Zampini . NeumannBoundaries - index set defining the Neumann boundaries 632785d1243SStefano Zampini 633785d1243SStefano Zampini Level: intermediate 634785d1243SStefano Zampini 635785d1243SStefano Zampini Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries 636785d1243SStefano Zampini 637785d1243SStefano Zampini .seealso: PCBDDC 638785d1243SStefano Zampini @*/ 639785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 640785d1243SStefano Zampini { 641785d1243SStefano Zampini PetscErrorCode ierr; 642785d1243SStefano Zampini 643785d1243SStefano Zampini PetscFunctionBegin; 644785d1243SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 645785d1243SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 646785d1243SStefano Zampini PetscFunctionReturn(0); 647785d1243SStefano Zampini } 648785d1243SStefano Zampini /* -------------------------------------------------------------------------- */ 649785d1243SStefano Zampini 650785d1243SStefano Zampini #undef __FUNCT__ 651785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC" 652785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries) 653785d1243SStefano Zampini { 654785d1243SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 655785d1243SStefano Zampini 656785d1243SStefano Zampini PetscFunctionBegin; 657785d1243SStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundariesLocal; 658785d1243SStefano Zampini PetscFunctionReturn(0); 659785d1243SStefano Zampini } 660785d1243SStefano Zampini 661785d1243SStefano Zampini #undef __FUNCT__ 66282ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal" 66353cdbc3dSStefano Zampini /*@ 66482ba6b80SStefano Zampini PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering) 66553cdbc3dSStefano Zampini 666785d1243SStefano Zampini Collective 66753cdbc3dSStefano Zampini 66853cdbc3dSStefano Zampini Input Parameters: 66928509bceSStefano Zampini . pc - the preconditioning context 67053cdbc3dSStefano Zampini 67153cdbc3dSStefano Zampini Output Parameters: 67228509bceSStefano Zampini . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 67353cdbc3dSStefano Zampini 67453cdbc3dSStefano Zampini Level: intermediate 67553cdbc3dSStefano Zampini 67653cdbc3dSStefano Zampini Notes: 67753cdbc3dSStefano Zampini 67853cdbc3dSStefano Zampini .seealso: PCBDDC 67953cdbc3dSStefano Zampini @*/ 68082ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries) 68153cdbc3dSStefano Zampini { 68253cdbc3dSStefano Zampini PetscErrorCode ierr; 68353cdbc3dSStefano Zampini 68453cdbc3dSStefano Zampini PetscFunctionBegin; 68553cdbc3dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 68682ba6b80SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 6870c7d97c5SJed Brown PetscFunctionReturn(0); 6880c7d97c5SJed Brown } 68936e030ebSStefano Zampini /* -------------------------------------------------------------------------- */ 6901e6b0712SBarry Smith 69136e030ebSStefano Zampini #undef __FUNCT__ 692da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 6931a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 69436e030ebSStefano Zampini { 69536e030ebSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 696da1bb401SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 697da1bb401SStefano Zampini PetscErrorCode ierr; 69836e030ebSStefano Zampini 69936e030ebSStefano Zampini PetscFunctionBegin; 700674ae819SStefano Zampini /* free old CSR */ 701674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 702fb223d50SStefano Zampini /* TODO: PCBDDCGraphSetAdjacency */ 703674ae819SStefano Zampini /* get CSR into graph structure */ 704da1bb401SStefano Zampini if (copymode == PETSC_COPY_VALUES) { 705785e854fSJed Brown ierr = PetscMalloc1((nvtxs+1),&mat_graph->xadj);CHKERRQ(ierr); 706785e854fSJed Brown ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr); 707674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 708674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 709da1bb401SStefano Zampini } else if (copymode == PETSC_OWN_POINTER) { 7101a83f524SJed Brown mat_graph->xadj = (PetscInt*)xadj; 7111a83f524SJed Brown mat_graph->adjncy = (PetscInt*)adjncy; 712674ae819SStefano Zampini } 713575ad6abSStefano Zampini mat_graph->nvtxs_csr = nvtxs; 71436e030ebSStefano Zampini PetscFunctionReturn(0); 71536e030ebSStefano Zampini } 7161e6b0712SBarry Smith 71736e030ebSStefano Zampini #undef __FUNCT__ 718da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 71936e030ebSStefano Zampini /*@ 72028509bceSStefano Zampini PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix 72136e030ebSStefano Zampini 72236e030ebSStefano Zampini Not collective 72336e030ebSStefano Zampini 72436e030ebSStefano Zampini Input Parameters: 72536e030ebSStefano Zampini + pc - the preconditioning context 72628509bceSStefano Zampini . nvtxs - number of local vertices of the graph (i.e., the local size of your problem) 72728509bceSStefano Zampini . xadj, adjncy - the CSR graph 72828509bceSStefano Zampini - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 72936e030ebSStefano Zampini 73036e030ebSStefano Zampini Level: intermediate 73136e030ebSStefano Zampini 73236e030ebSStefano Zampini Notes: 73336e030ebSStefano Zampini 73428509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode 73536e030ebSStefano Zampini @*/ 7361a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 73736e030ebSStefano Zampini { 738575ad6abSStefano Zampini void (*f)(void) = 0; 73936e030ebSStefano Zampini PetscErrorCode ierr; 74036e030ebSStefano Zampini 74136e030ebSStefano Zampini PetscFunctionBegin; 74236e030ebSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 743674ae819SStefano Zampini PetscValidIntPointer(xadj,3); 744674ae819SStefano Zampini PetscValidIntPointer(xadj,4); 745674ae819SStefano Zampini if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 746674ae819SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 747da1bb401SStefano Zampini } 74836e030ebSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 749575ad6abSStefano Zampini /* free arrays if PCBDDC is not the PC type */ 750575ad6abSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 751575ad6abSStefano Zampini if (!f && copymode == PETSC_OWN_POINTER) { 752575ad6abSStefano Zampini ierr = PetscFree(xadj);CHKERRQ(ierr); 753575ad6abSStefano Zampini ierr = PetscFree(adjncy);CHKERRQ(ierr); 75436e030ebSStefano Zampini } 75536e030ebSStefano Zampini PetscFunctionReturn(0); 75636e030ebSStefano Zampini } 7579c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */ 7581e6b0712SBarry Smith 7599c0446d6SStefano Zampini #undef __FUNCT__ 76063602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC" 76163602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 76263602bcaSStefano Zampini { 76363602bcaSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 76463602bcaSStefano Zampini PetscInt i; 76563602bcaSStefano Zampini PetscErrorCode ierr; 76663602bcaSStefano Zampini 76763602bcaSStefano Zampini PetscFunctionBegin; 76863602bcaSStefano Zampini /* Destroy ISes if they were already set */ 76963602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 77063602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 77163602bcaSStefano Zampini } 77263602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 77363602bcaSStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 77463602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 77563602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 77663602bcaSStefano Zampini } 77763602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 77863602bcaSStefano Zampini pcbddc->n_ISForDofs = 0; 77963602bcaSStefano Zampini /* allocate space then set */ 78063602bcaSStefano Zampini ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr); 78163602bcaSStefano Zampini for (i=0;i<n_is;i++) { 78263602bcaSStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 78363602bcaSStefano Zampini pcbddc->ISForDofsLocal[i]=ISForDofs[i]; 78463602bcaSStefano Zampini } 78563602bcaSStefano Zampini pcbddc->n_ISForDofsLocal=n_is; 78663602bcaSStefano Zampini if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 7871a8a16d1SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 78863602bcaSStefano Zampini PetscFunctionReturn(0); 78963602bcaSStefano Zampini } 79063602bcaSStefano Zampini 79163602bcaSStefano Zampini #undef __FUNCT__ 79263602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal" 79363602bcaSStefano Zampini /*@ 79463602bcaSStefano Zampini PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix 79563602bcaSStefano Zampini 79663602bcaSStefano Zampini Collective 79763602bcaSStefano Zampini 79863602bcaSStefano Zampini Input Parameters: 79963602bcaSStefano Zampini + pc - the preconditioning context 80063602bcaSStefano Zampini - n_is - number of index sets defining the fields 80163602bcaSStefano Zampini . ISForDofs - array of IS describing the fields in local ordering 80263602bcaSStefano Zampini 80363602bcaSStefano Zampini Level: intermediate 80463602bcaSStefano Zampini 80563602bcaSStefano 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. 80663602bcaSStefano Zampini 80763602bcaSStefano Zampini .seealso: PCBDDC 80863602bcaSStefano Zampini @*/ 80963602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[]) 81063602bcaSStefano Zampini { 81163602bcaSStefano Zampini PetscInt i; 81263602bcaSStefano Zampini PetscErrorCode ierr; 81363602bcaSStefano Zampini 81463602bcaSStefano Zampini PetscFunctionBegin; 81563602bcaSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 81663602bcaSStefano Zampini PetscValidLogicalCollectiveInt(pc,n_is,2); 81763602bcaSStefano Zampini for (i=0;i<n_is;i++) { 81863602bcaSStefano Zampini PetscCheckSameComm(pc,1,ISForDofs[i],3); 81963602bcaSStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 82063602bcaSStefano Zampini } 82163602bcaSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 82263602bcaSStefano Zampini PetscFunctionReturn(0); 82363602bcaSStefano Zampini } 82463602bcaSStefano Zampini /* -------------------------------------------------------------------------- */ 82563602bcaSStefano Zampini 82663602bcaSStefano Zampini #undef __FUNCT__ 8279c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 8289c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 8299c0446d6SStefano Zampini { 8309c0446d6SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 8319c0446d6SStefano Zampini PetscInt i; 8329c0446d6SStefano Zampini PetscErrorCode ierr; 8339c0446d6SStefano Zampini 8349c0446d6SStefano Zampini PetscFunctionBegin; 835da1bb401SStefano Zampini /* Destroy ISes if they were already set */ 8369c0446d6SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 8379c0446d6SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 8389c0446d6SStefano Zampini } 839d11ae9bbSstefano_zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 84063602bcaSStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 84163602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 84263602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 84363602bcaSStefano Zampini } 84463602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 84563602bcaSStefano Zampini pcbddc->n_ISForDofsLocal = 0; 846da1bb401SStefano Zampini /* allocate space then set */ 847785e854fSJed Brown ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr); 8489c0446d6SStefano Zampini for (i=0;i<n_is;i++) { 849da1bb401SStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 850da1bb401SStefano Zampini pcbddc->ISForDofs[i]=ISForDofs[i]; 8519c0446d6SStefano Zampini } 8529c0446d6SStefano Zampini pcbddc->n_ISForDofs=n_is; 85363602bcaSStefano Zampini if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 8541a8a16d1SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 8559c0446d6SStefano Zampini PetscFunctionReturn(0); 8569c0446d6SStefano Zampini } 8571e6b0712SBarry Smith 8589c0446d6SStefano Zampini #undef __FUNCT__ 8599c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting" 8609c0446d6SStefano Zampini /*@ 86163602bcaSStefano Zampini PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix 8629c0446d6SStefano Zampini 86363602bcaSStefano Zampini Collective 8649c0446d6SStefano Zampini 8659c0446d6SStefano Zampini Input Parameters: 8669c0446d6SStefano Zampini + pc - the preconditioning context 86728509bceSStefano Zampini - n_is - number of index sets defining the fields 86863602bcaSStefano Zampini . ISForDofs - array of IS describing the fields in global ordering 8699c0446d6SStefano Zampini 8709c0446d6SStefano Zampini Level: intermediate 8719c0446d6SStefano Zampini 87263602bcaSStefano Zampini Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field. 8739c0446d6SStefano Zampini 8749c0446d6SStefano Zampini .seealso: PCBDDC 8759c0446d6SStefano Zampini @*/ 8769c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 8779c0446d6SStefano Zampini { 8782b510759SStefano Zampini PetscInt i; 8799c0446d6SStefano Zampini PetscErrorCode ierr; 8809c0446d6SStefano Zampini 8819c0446d6SStefano Zampini PetscFunctionBegin; 8829c0446d6SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 88363602bcaSStefano Zampini PetscValidLogicalCollectiveInt(pc,n_is,2); 8842b510759SStefano Zampini for (i=0;i<n_is;i++) { 88563602bcaSStefano Zampini PetscCheckSameComm(pc,1,ISForDofs[i],3); 88663602bcaSStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 8872b510759SStefano Zampini } 8889c0446d6SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 8899c0446d6SStefano Zampini PetscFunctionReturn(0); 8909c0446d6SStefano Zampini } 891da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 892534831adSStefano Zampini #undef __FUNCT__ 893534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC" 894534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 895534831adSStefano Zampini /* 896534831adSStefano Zampini PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 897534831adSStefano Zampini guess if a transformation of basis approach has been selected. 8989c0446d6SStefano Zampini 899534831adSStefano Zampini Input Parameter: 900534831adSStefano Zampini + pc - the preconditioner contex 901534831adSStefano Zampini 902534831adSStefano Zampini Application Interface Routine: PCPreSolve() 903534831adSStefano Zampini 904534831adSStefano Zampini Notes: 905534831adSStefano Zampini The interface routine PCPreSolve() is not usually called directly by 906534831adSStefano Zampini the user, but instead is called by KSPSolve(). 907534831adSStefano Zampini */ 908534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 909534831adSStefano Zampini { 910534831adSStefano Zampini PetscErrorCode ierr; 911534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 912534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 913534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 914534831adSStefano Zampini Mat temp_mat; 9153972b0daSStefano Zampini IS dirIS; 9163972b0daSStefano Zampini Vec used_vec; 91782ba6b80SStefano Zampini PetscBool guess_nonzero; 918534831adSStefano Zampini 919534831adSStefano Zampini PetscFunctionBegin; 92085c4d303SStefano Zampini /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 92185c4d303SStefano Zampini if (ksp) { 92285c4d303SStefano Zampini PetscBool iscg; 92385c4d303SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 92485c4d303SStefano Zampini if (!iscg) { 92585c4d303SStefano Zampini ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 92685c4d303SStefano Zampini } 92785c4d303SStefano Zampini } 92885c4d303SStefano Zampini /* Creates parallel work vectors used in presolve */ 92962a6ff1dSStefano Zampini if (!pcbddc->original_rhs) { 93062a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 93162a6ff1dSStefano Zampini } 93262a6ff1dSStefano Zampini if (!pcbddc->temp_solution) { 93362a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 93462a6ff1dSStefano Zampini } 9353972b0daSStefano Zampini if (x) { 9363972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 9373972b0daSStefano Zampini used_vec = x; 9383972b0daSStefano Zampini } else { 9393972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 9403972b0daSStefano Zampini used_vec = pcbddc->temp_solution; 9413972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 9423972b0daSStefano Zampini } 9433972b0daSStefano Zampini /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 9443972b0daSStefano Zampini if (ksp) { 9453972b0daSStefano Zampini ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 9463972b0daSStefano Zampini if (!guess_nonzero) { 9473972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 9483972b0daSStefano Zampini } 9493972b0daSStefano Zampini } 9503308cffdSStefano Zampini 9513972b0daSStefano Zampini /* store the original rhs */ 9523972b0daSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 9533972b0daSStefano Zampini 9543972b0daSStefano Zampini /* Take into account zeroed rows -> change rhs and store solution removed */ 955785d1243SStefano Zampini /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */ 95682ba6b80SStefano Zampini ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr); 95782ba6b80SStefano Zampini if (rhs && dirIS) { 958785d1243SStefano Zampini PetscInt dirsize,i,*is_indices; 959785d1243SStefano Zampini PetscScalar *array_x,*array_diagonal; 960785d1243SStefano Zampini 9613972b0daSStefano Zampini ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 9623972b0daSStefano Zampini ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 9633972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9643972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9653972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9663972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 96782ba6b80SStefano Zampini ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr); 9683972b0daSStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 9693972b0daSStefano Zampini ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 9703972b0daSStefano Zampini ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 9712fa5cd67SKarl Rupp for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 9723972b0daSStefano Zampini ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 9733972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 9743972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 9753972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9763972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 977b76ba322SStefano Zampini 9783972b0daSStefano Zampini /* remove the computed solution from the rhs */ 9793972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 9803972b0daSStefano Zampini ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 9813972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 9823308cffdSStefano Zampini } 983b76ba322SStefano Zampini 984b76ba322SStefano Zampini /* store partially computed solution and set initial guess */ 9853972b0daSStefano Zampini if (x) { 9863972b0daSStefano Zampini ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 9873972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 98885c4d303SStefano Zampini if (pcbddc->use_exact_dirichlet_trick) { 989b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 990b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 991b76ba322SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 992b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 993b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 994b76ba322SStefano Zampini if (ksp) { 995b76ba322SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 996b76ba322SStefano Zampini } 997b76ba322SStefano Zampini } 9983972b0daSStefano Zampini } 999b76ba322SStefano Zampini 100092e3dcfbSStefano Zampini /* prepare MatMult and rhs for solver */ 1001b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix) { 1002b76ba322SStefano Zampini /* swap pointers for local matrices */ 1003b76ba322SStefano Zampini temp_mat = matis->A; 1004b76ba322SStefano Zampini matis->A = pcbddc->local_mat; 1005b76ba322SStefano Zampini pcbddc->local_mat = temp_mat; 100692e3dcfbSStefano Zampini if (rhs) { 1007b76ba322SStefano Zampini /* Get local rhs and apply transformation of basis */ 1008b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1009b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1010b76ba322SStefano Zampini /* from original basis to modified basis */ 1011b76ba322SStefano Zampini ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 1012b76ba322SStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 1013b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1014b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1015674ae819SStefano Zampini } 101692e3dcfbSStefano Zampini } 101792e3dcfbSStefano Zampini 101892e3dcfbSStefano Zampini /* remove nullspace if present */ 10190bdf917eSStefano Zampini if (ksp && pcbddc->NullSpace) { 1020d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 1021d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 1022b76ba322SStefano Zampini } 10230bdf917eSStefano Zampini ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 1024534831adSStefano Zampini PetscFunctionReturn(0); 1025534831adSStefano Zampini } 1026534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 1027534831adSStefano Zampini #undef __FUNCT__ 1028534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC" 1029534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 1030534831adSStefano Zampini /* 1031534831adSStefano Zampini PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 1032534831adSStefano Zampini approach has been selected. Also, restores rhs to its original state. 1033534831adSStefano Zampini 1034534831adSStefano Zampini Input Parameter: 1035534831adSStefano Zampini + pc - the preconditioner contex 1036534831adSStefano Zampini 1037534831adSStefano Zampini Application Interface Routine: PCPostSolve() 1038534831adSStefano Zampini 1039534831adSStefano Zampini Notes: 1040534831adSStefano Zampini The interface routine PCPostSolve() is not usually called directly by 1041534831adSStefano Zampini the user, but instead is called by KSPSolve(). 1042534831adSStefano Zampini */ 1043534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1044534831adSStefano Zampini { 1045534831adSStefano Zampini PetscErrorCode ierr; 1046534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1047534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 1048534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 1049534831adSStefano Zampini Mat temp_mat; 1050534831adSStefano Zampini 1051534831adSStefano Zampini PetscFunctionBegin; 1052b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix) { 1053534831adSStefano Zampini /* swap pointers for local matrices */ 1054534831adSStefano Zampini temp_mat = matis->A; 1055534831adSStefano Zampini matis->A = pcbddc->local_mat; 1056534831adSStefano Zampini pcbddc->local_mat = temp_mat; 10573425bc38SStefano Zampini } 1058b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix && x) { 1059534831adSStefano Zampini /* Get Local boundary and apply transformation of basis to solution vector */ 1060534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1061534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1062534831adSStefano Zampini /* from modified basis to original basis */ 1063534831adSStefano Zampini ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 1064534831adSStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 1065534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1066534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1067534831adSStefano Zampini } 10683972b0daSStefano Zampini /* add solution removed in presolve */ 10693425bc38SStefano Zampini if (x) { 10703425bc38SStefano Zampini ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 10713425bc38SStefano Zampini } 1072fb223d50SStefano Zampini /* restore rhs to its original state */ 1073fb223d50SStefano Zampini if (rhs) { 1074fb223d50SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 1075fb223d50SStefano Zampini } 1076534831adSStefano Zampini PetscFunctionReturn(0); 1077534831adSStefano Zampini } 1078534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 107953cdbc3dSStefano Zampini #undef __FUNCT__ 108053cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC" 10810c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 10820c7d97c5SJed Brown /* 10830c7d97c5SJed Brown PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 10840c7d97c5SJed Brown by setting data structures and options. 10850c7d97c5SJed Brown 10860c7d97c5SJed Brown Input Parameter: 108753cdbc3dSStefano Zampini + pc - the preconditioner context 10880c7d97c5SJed Brown 10890c7d97c5SJed Brown Application Interface Routine: PCSetUp() 10900c7d97c5SJed Brown 10910c7d97c5SJed Brown Notes: 10920c7d97c5SJed Brown The interface routine PCSetUp() is not usually called directly by 10930c7d97c5SJed Brown the user, but instead is called by PCApply() if necessary. 10940c7d97c5SJed Brown */ 109553cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc) 10960c7d97c5SJed Brown { 10970c7d97c5SJed Brown PetscErrorCode ierr; 10980c7d97c5SJed Brown PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 1099f4ddd8eeSStefano Zampini MatNullSpace nearnullspace; 1100674ae819SStefano Zampini PetscBool computeis,computetopography,computesolvers; 1101165b64e2SStefano Zampini PetscBool new_nearnullspace_provided; 11020c7d97c5SJed Brown 11030c7d97c5SJed Brown PetscFunctionBegin; 1104f4ddd8eeSStefano Zampini /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 11053b03a366Sstefano_zampini /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1106fcd409f5SStefano Zampini Also, BDDC directly build the Dirichlet problem */ 1107f4ddd8eeSStefano Zampini /* split work */ 1108674ae819SStefano Zampini if (pc->setupcalled) { 1109674ae819SStefano Zampini computeis = PETSC_FALSE; 1110d1e9a80fSBarry Smith if (pc->flag == SAME_NONZERO_PATTERN) { 1111674ae819SStefano Zampini computetopography = PETSC_FALSE; 1112674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1113674ae819SStefano Zampini } else { /* DIFFERENT_NONZERO_PATTERN */ 1114674ae819SStefano Zampini computetopography = PETSC_TRUE; 1115674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1116674ae819SStefano Zampini } 1117674ae819SStefano Zampini } else { 1118674ae819SStefano Zampini computeis = PETSC_TRUE; 1119674ae819SStefano Zampini computetopography = PETSC_TRUE; 1120674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1121674ae819SStefano Zampini } 1122fb180af4SStefano Zampini if (pcbddc->recompute_topography) { 1123fb180af4SStefano Zampini computetopography = PETSC_TRUE; 1124fb180af4SStefano Zampini } 1125f4ddd8eeSStefano Zampini 1126f4ddd8eeSStefano Zampini /* Get stdout for dbg */ 112770cf5478SStefano Zampini if (pcbddc->dbg_flag) { 112870cf5478SStefano Zampini if (!pcbddc->dbg_viewer) { 112958a03d70SStefano Zampini pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)); 1130f4ddd8eeSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 113170cf5478SStefano Zampini } 113258a03d70SStefano Zampini ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1133f4ddd8eeSStefano Zampini } 1134f4ddd8eeSStefano Zampini 1135fcd409f5SStefano Zampini /* Set up all the "iterative substructuring" common block without computing solvers */ 1136674ae819SStefano Zampini if (computeis) { 1137fcd409f5SStefano Zampini /* HACK INTO PCIS */ 1138fcd409f5SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 1139fcd409f5SStefano Zampini pcis->computesolvers = PETSC_FALSE; 1140674ae819SStefano Zampini ierr = PCISSetUp(pc);CHKERRQ(ierr); 114139e2fb2aSStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr); 1142674ae819SStefano Zampini } 1143f4ddd8eeSStefano Zampini 1144b9b85e73SStefano Zampini /* check user defined change of basis (if any) */ 1145b9b85e73SStefano Zampini if (pcbddc->user_ChangeOfBasisMatrix) { 1146b9b85e73SStefano Zampini PC_IS* pcis= (PC_IS*)pc->data; 1147b9b85e73SStefano Zampini PetscInt n,m; 1148b9b85e73SStefano Zampini ierr = MatGetSize(pcbddc->user_ChangeOfBasisMatrix,&n,&m);CHKERRQ(ierr); 1149b9b85e73SStefano Zampini if (n != m) { 1150b9b85e73SStefano Zampini SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Change of basis matrix should be square %d != %d\n",n,m); 1151b9b85e73SStefano Zampini } else if (n != pcis->n_B && n != pcis->n) { 1152b9b85e73SStefano 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); 1153b9b85e73SStefano Zampini } 1154b9b85e73SStefano Zampini /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */ 1155b9b85e73SStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 1156b9b85e73SStefano Zampini } 1157f4ddd8eeSStefano Zampini /* Analyze interface */ 1158674ae819SStefano Zampini if (computetopography) { 1159674ae819SStefano Zampini ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 1160fb8d54d4SStefano Zampini } 1161fb8d54d4SStefano Zampini 1162f4ddd8eeSStefano Zampini /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1163fb8d54d4SStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 1164f4ddd8eeSStefano Zampini ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr); 1165f4ddd8eeSStefano Zampini if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1166f4ddd8eeSStefano Zampini if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1167f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1168f4ddd8eeSStefano Zampini } else { 1169f4ddd8eeSStefano Zampini /* determine if the two nullspaces are different (should be lightweight) */ 1170f4ddd8eeSStefano Zampini if (nearnullspace != pcbddc->onearnullspace) { 1171f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1172165b64e2SStefano Zampini } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1173f4ddd8eeSStefano Zampini PetscInt i; 1174165b64e2SStefano Zampini const Vec *nearnullvecs; 1175165b64e2SStefano Zampini PetscObjectState state; 1176165b64e2SStefano Zampini PetscInt nnsp_size; 1177165b64e2SStefano Zampini ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 1178f4ddd8eeSStefano Zampini for (i=0;i<nnsp_size;i++) { 1179f4ddd8eeSStefano Zampini ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr); 1180165b64e2SStefano Zampini if (pcbddc->onearnullvecs_state[i] != state) { 1181f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1182f4ddd8eeSStefano Zampini break; 1183f4ddd8eeSStefano Zampini } 1184f4ddd8eeSStefano Zampini } 1185f4ddd8eeSStefano Zampini } 1186f4ddd8eeSStefano Zampini } 1187f4ddd8eeSStefano Zampini } else { 1188f4ddd8eeSStefano Zampini if (!nearnullspace) { /* both nearnullspaces are null */ 1189f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 1190f4ddd8eeSStefano Zampini } else { /* nearnullspace attached later */ 1191f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1192f4ddd8eeSStefano Zampini } 1193f4ddd8eeSStefano Zampini } 1194f4ddd8eeSStefano Zampini 1195f4ddd8eeSStefano Zampini /* Setup constraints and related work vectors */ 1196727cdba6SStefano Zampini /* reset primal space flags */ 1197f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1198727cdba6SStefano Zampini pcbddc->new_primal_space_local = PETSC_FALSE; 1199fb8d54d4SStefano Zampini if (computetopography || new_nearnullspace_provided) { 1200727cdba6SStefano Zampini /* It also sets the primal space flags */ 1201674ae819SStefano Zampini ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 1202e7b262bdSStefano Zampini /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1203f4ddd8eeSStefano Zampini ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr); 1204674ae819SStefano Zampini } 1205fb8d54d4SStefano Zampini 1206f4ddd8eeSStefano Zampini if (computesolvers || pcbddc->new_primal_space) { 1207674ae819SStefano Zampini /* reset data */ 1208674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 120999cc7994SStefano Zampini /* Create coarse and local stuffs */ 121099cc7994SStefano Zampini ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 1211674ae819SStefano Zampini ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 12120c7d97c5SJed Brown } 121370cf5478SStefano Zampini 121458a03d70SStefano Zampini if (pcbddc->dbg_flag) { 121558a03d70SStefano Zampini ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 12162b510759SStefano Zampini } 12170c7d97c5SJed Brown PetscFunctionReturn(0); 12180c7d97c5SJed Brown } 12190c7d97c5SJed Brown 12200c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 12210c7d97c5SJed Brown /* 122250efa1b5SStefano Zampini PCApply_BDDC - Applies the BDDC operator to a vector. 12230c7d97c5SJed Brown 12240c7d97c5SJed Brown Input Parameters: 12250c7d97c5SJed Brown . pc - the preconditioner context 12260c7d97c5SJed Brown . r - input vector (global) 12270c7d97c5SJed Brown 12280c7d97c5SJed Brown Output Parameter: 12290c7d97c5SJed Brown . z - output vector (global) 12300c7d97c5SJed Brown 12310c7d97c5SJed Brown Application Interface Routine: PCApply() 12320c7d97c5SJed Brown */ 12330c7d97c5SJed Brown #undef __FUNCT__ 12340c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC" 123553cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 12360c7d97c5SJed Brown { 12370c7d97c5SJed Brown PC_IS *pcis = (PC_IS*)(pc->data); 12380c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 12390c7d97c5SJed Brown PetscErrorCode ierr; 12403b03a366Sstefano_zampini const PetscScalar one = 1.0; 12413b03a366Sstefano_zampini const PetscScalar m_one = -1.0; 12422617d88aSStefano Zampini const PetscScalar zero = 0.0; 12430c7d97c5SJed Brown 12440c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN 12450c7d97c5SJed Brown NN interface preconditioner changed to BDDC 12468eeda7d8SStefano Zampini Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 12470c7d97c5SJed Brown 12480c7d97c5SJed Brown PetscFunctionBegin; 124985c4d303SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 12500c7d97c5SJed Brown /* First Dirichlet solve */ 12510c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12520c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 125353cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 12540c7d97c5SJed Brown /* 12550c7d97c5SJed Brown Assembling right hand side for BDDC operator 1256674ae819SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 1257674ae819SStefano Zampini - pcis->vec1_B the interface part of the global vector z 12580c7d97c5SJed Brown */ 12590c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 12600c7d97c5SJed Brown ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 12618eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 12620c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 12630c7d97c5SJed Brown ierr = VecCopy(r,z);CHKERRQ(ierr); 12640c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12650c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1266674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1267b76ba322SStefano Zampini } else { 12680bdf917eSStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 1269b76ba322SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 1270674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1271b76ba322SStefano Zampini } 1272b76ba322SStefano Zampini 12732617d88aSStefano Zampini /* Apply interface preconditioner 12742617d88aSStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1275dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr); 12762617d88aSStefano Zampini 1277674ae819SStefano Zampini /* Apply transpose of partition of unity operator */ 1278674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 12790c7d97c5SJed Brown 12803b03a366Sstefano_zampini /* Second Dirichlet solve and assembling of output */ 12810c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12820c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12830c7d97c5SJed Brown ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 12848eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1285df187020SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1286df187020SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1287df187020SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1288df187020SStefano Zampini ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 12890c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12900c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12910c7d97c5SJed Brown PetscFunctionReturn(0); 12920c7d97c5SJed Brown } 129350efa1b5SStefano Zampini 129450efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */ 129550efa1b5SStefano Zampini /* 129650efa1b5SStefano Zampini PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector. 129750efa1b5SStefano Zampini 129850efa1b5SStefano Zampini Input Parameters: 129950efa1b5SStefano Zampini . pc - the preconditioner context 130050efa1b5SStefano Zampini . r - input vector (global) 130150efa1b5SStefano Zampini 130250efa1b5SStefano Zampini Output Parameter: 130350efa1b5SStefano Zampini . z - output vector (global) 130450efa1b5SStefano Zampini 130550efa1b5SStefano Zampini Application Interface Routine: PCApplyTranspose() 130650efa1b5SStefano Zampini */ 130750efa1b5SStefano Zampini #undef __FUNCT__ 130850efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC" 130950efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z) 131050efa1b5SStefano Zampini { 131150efa1b5SStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 131250efa1b5SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 131350efa1b5SStefano Zampini PetscErrorCode ierr; 131450efa1b5SStefano Zampini const PetscScalar one = 1.0; 131550efa1b5SStefano Zampini const PetscScalar m_one = -1.0; 131650efa1b5SStefano Zampini const PetscScalar zero = 0.0; 131750efa1b5SStefano Zampini 131850efa1b5SStefano Zampini PetscFunctionBegin; 131950efa1b5SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 132050efa1b5SStefano Zampini /* First Dirichlet solve */ 132150efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 132250efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 132350efa1b5SStefano Zampini ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 132450efa1b5SStefano Zampini /* 132550efa1b5SStefano Zampini Assembling right hand side for BDDC operator 132650efa1b5SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 132750efa1b5SStefano Zampini - pcis->vec1_B the interface part of the global vector z 132850efa1b5SStefano Zampini */ 132950efa1b5SStefano Zampini ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 133050efa1b5SStefano Zampini ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 133150efa1b5SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 133250efa1b5SStefano Zampini ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 133350efa1b5SStefano Zampini ierr = VecCopy(r,z);CHKERRQ(ierr); 133450efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 133550efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 133650efa1b5SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 133750efa1b5SStefano Zampini } else { 133850efa1b5SStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 133950efa1b5SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 134050efa1b5SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 134150efa1b5SStefano Zampini } 134250efa1b5SStefano Zampini 134350efa1b5SStefano Zampini /* Apply interface preconditioner 134450efa1b5SStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1345dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr); 134650efa1b5SStefano Zampini 134750efa1b5SStefano Zampini /* Apply transpose of partition of unity operator */ 134850efa1b5SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 134950efa1b5SStefano Zampini 135050efa1b5SStefano Zampini /* Second Dirichlet solve and assembling of output */ 135150efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 135250efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 135350efa1b5SStefano Zampini ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 135450efa1b5SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1355b0147a47SStefano Zampini ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1356b0147a47SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1357b0147a47SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1358b0147a47SStefano Zampini ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 135950efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 136050efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 136150efa1b5SStefano Zampini PetscFunctionReturn(0); 136250efa1b5SStefano Zampini } 1363da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 1364674ae819SStefano Zampini 1365da1bb401SStefano Zampini #undef __FUNCT__ 1366da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC" 1367da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc) 1368da1bb401SStefano Zampini { 1369da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1370da1bb401SStefano Zampini PetscErrorCode ierr; 1371da1bb401SStefano Zampini 1372da1bb401SStefano Zampini PetscFunctionBegin; 1373da1bb401SStefano Zampini /* free data created by PCIS */ 1374da1bb401SStefano Zampini ierr = PCISDestroy(pc);CHKERRQ(ierr); 1375674ae819SStefano Zampini /* free BDDC custom data */ 1376674ae819SStefano Zampini ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 1377674ae819SStefano Zampini /* destroy objects related to topography */ 1378674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 1379674ae819SStefano Zampini /* free allocated graph structure */ 1380da1bb401SStefano Zampini ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 1381674ae819SStefano Zampini /* free data for scaling operator */ 1382674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 1383674ae819SStefano Zampini /* free solvers stuff */ 1384674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 138562a6ff1dSStefano Zampini /* free global vectors needed in presolve */ 138662a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 138762a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 138839e2fb2aSStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr); 13893425bc38SStefano Zampini /* remove functions */ 1390b9b85e73SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisLocalMat_C",NULL);CHKERRQ(ierr); 1391674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 1392bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 13932b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 1394b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 13952b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 1396bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 1397bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 139882ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1399bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 140082ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1401bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 140282ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1403bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1404785d1243SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1405bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 140663602bcaSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr); 1407bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 1408bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 1409bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 1410bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 1411674ae819SStefano Zampini /* Free the private data structure */ 1412674ae819SStefano Zampini ierr = PetscFree(pc->data);CHKERRQ(ierr); 1413da1bb401SStefano Zampini PetscFunctionReturn(0); 1414da1bb401SStefano Zampini } 14153425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 14161e6b0712SBarry Smith 14173425bc38SStefano Zampini #undef __FUNCT__ 14183425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 14193425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 14203425bc38SStefano Zampini { 1421674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 14223425bc38SStefano Zampini PC_IS* pcis; 14233425bc38SStefano Zampini PC_BDDC* pcbddc; 14243425bc38SStefano Zampini PetscErrorCode ierr; 14250c7d97c5SJed Brown 14263425bc38SStefano Zampini PetscFunctionBegin; 14273425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 14283425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 14293425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 14303425bc38SStefano Zampini 14313425bc38SStefano Zampini /* change of basis for physical rhs if needed 14323425bc38SStefano Zampini It also changes the rhs in case of dirichlet boundaries */ 14333308cffdSStefano Zampini ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 14343425bc38SStefano Zampini /* store vectors for computation of fetidp final solution */ 14353425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14363425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1437fb223d50SStefano Zampini /* scale rhs since it should be unassembled */ 1438fb223d50SStefano Zampini /* TODO use counter scaling? (also below) */ 14393425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14403425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1441674ae819SStefano Zampini /* Apply partition of unity */ 14423425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1443674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 14448eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 14453425bc38SStefano Zampini /* compute partially subassembled Schur complement right-hand side */ 14463425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 14473425bc38SStefano Zampini ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 14483425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 14493425bc38SStefano Zampini ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 14503425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14513425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1452674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 14533425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14543425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14553425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 14563425bc38SStefano Zampini } 14573425bc38SStefano Zampini /* BDDC rhs */ 14583425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 14598eeda7d8SStefano Zampini if (pcbddc->switch_static) { 14603425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 14613425bc38SStefano Zampini } 14623425bc38SStefano Zampini /* apply BDDC */ 1463dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 14643425bc38SStefano Zampini /* Application of B_delta and assembling of rhs for fetidp fluxes */ 14653425bc38SStefano Zampini ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 14663425bc38SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 14673425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14683425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14693425bc38SStefano Zampini /* restore original rhs */ 14703425bc38SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 14713425bc38SStefano Zampini PetscFunctionReturn(0); 14723425bc38SStefano Zampini } 14731e6b0712SBarry Smith 14743425bc38SStefano Zampini #undef __FUNCT__ 14753425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 14763425bc38SStefano Zampini /*@ 147728509bceSStefano Zampini PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system 14783425bc38SStefano Zampini 14793425bc38SStefano Zampini Collective 14803425bc38SStefano Zampini 14813425bc38SStefano Zampini Input Parameters: 148228509bceSStefano Zampini + fetidp_mat - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators 148328509bceSStefano Zampini . standard_rhs - the right-hand side for your linear system 14843425bc38SStefano Zampini 14853425bc38SStefano Zampini Output Parameters: 148628509bceSStefano Zampini - fetidp_flux_rhs - the right-hand side for the FETIDP linear system 14873425bc38SStefano Zampini 14883425bc38SStefano Zampini Level: developer 14893425bc38SStefano Zampini 14903425bc38SStefano Zampini Notes: 14913425bc38SStefano Zampini 149228509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 14933425bc38SStefano Zampini @*/ 14943425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 14953425bc38SStefano Zampini { 1496674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 14973425bc38SStefano Zampini PetscErrorCode ierr; 14983425bc38SStefano Zampini 14993425bc38SStefano Zampini PetscFunctionBegin; 15003425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15013425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 15023425bc38SStefano Zampini PetscFunctionReturn(0); 15033425bc38SStefano Zampini } 15043425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 15051e6b0712SBarry Smith 15063425bc38SStefano Zampini #undef __FUNCT__ 15073425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 15083425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 15093425bc38SStefano Zampini { 1510674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 15113425bc38SStefano Zampini PC_IS* pcis; 15123425bc38SStefano Zampini PC_BDDC* pcbddc; 15133425bc38SStefano Zampini PetscErrorCode ierr; 15143425bc38SStefano Zampini 15153425bc38SStefano Zampini PetscFunctionBegin; 15163425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15173425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 15183425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 15193425bc38SStefano Zampini 15203425bc38SStefano Zampini /* apply B_delta^T */ 15213425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15223425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15233425bc38SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 15243425bc38SStefano Zampini /* compute rhs for BDDC application */ 15253425bc38SStefano Zampini ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 15268eeda7d8SStefano Zampini if (pcbddc->switch_static) { 15273425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 15283425bc38SStefano Zampini } 15293425bc38SStefano Zampini /* apply BDDC */ 1530dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 15313425bc38SStefano Zampini /* put values into standard global vector */ 15323425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15333425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15348eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 15353425bc38SStefano Zampini /* compute values into the interior if solved for the partially subassembled Schur complement */ 15363425bc38SStefano Zampini ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 15373425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 15383425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 15393425bc38SStefano Zampini } 15403425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15413425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15423425bc38SStefano Zampini /* final change of basis if needed 15433425bc38SStefano Zampini Is also sums the dirichlet part removed during RHS assembling */ 15443308cffdSStefano Zampini ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 15453425bc38SStefano Zampini PetscFunctionReturn(0); 15463425bc38SStefano Zampini } 15471e6b0712SBarry Smith 15483425bc38SStefano Zampini #undef __FUNCT__ 15493425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 15503425bc38SStefano Zampini /*@ 155128509bceSStefano Zampini PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system 15523425bc38SStefano Zampini 15533425bc38SStefano Zampini Collective 15543425bc38SStefano Zampini 15553425bc38SStefano Zampini Input Parameters: 155628509bceSStefano Zampini + fetidp_mat - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators 155728509bceSStefano Zampini . fetidp_flux_sol - the solution of the FETIDP linear system 15583425bc38SStefano Zampini 15593425bc38SStefano Zampini Output Parameters: 156028509bceSStefano Zampini - standard_sol - the solution defined on the physical domain 15613425bc38SStefano Zampini 15623425bc38SStefano Zampini Level: developer 15633425bc38SStefano Zampini 15643425bc38SStefano Zampini Notes: 15653425bc38SStefano Zampini 156628509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 15673425bc38SStefano Zampini @*/ 15683425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 15693425bc38SStefano Zampini { 1570674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 15713425bc38SStefano Zampini PetscErrorCode ierr; 15723425bc38SStefano Zampini 15733425bc38SStefano Zampini PetscFunctionBegin; 15743425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15753425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 15763425bc38SStefano Zampini PetscFunctionReturn(0); 15773425bc38SStefano Zampini } 15783425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 15791e6b0712SBarry Smith 1580f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1581edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec); 1582f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1583f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1584edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec); 1585f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1586674ae819SStefano Zampini 15873425bc38SStefano Zampini #undef __FUNCT__ 15883425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 15893425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 15903425bc38SStefano Zampini { 1591674ae819SStefano Zampini 1592674ae819SStefano Zampini FETIDPMat_ctx fetidpmat_ctx; 15933425bc38SStefano Zampini Mat newmat; 1594674ae819SStefano Zampini FETIDPPC_ctx fetidppc_ctx; 15953425bc38SStefano Zampini PC newpc; 1596ce94432eSBarry Smith MPI_Comm comm; 15973425bc38SStefano Zampini PetscErrorCode ierr; 15983425bc38SStefano Zampini 15993425bc38SStefano Zampini PetscFunctionBegin; 1600ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 16013425bc38SStefano Zampini /* FETIDP linear matrix */ 16023425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 16033425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 16043425bc38SStefano Zampini ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 16053425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 1606edf7251bSStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr); 16073425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 16083425bc38SStefano Zampini ierr = MatSetUp(newmat);CHKERRQ(ierr); 16093425bc38SStefano Zampini /* FETIDP preconditioner */ 16103425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 16113425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 16123425bc38SStefano Zampini ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 16133425bc38SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 16143425bc38SStefano Zampini ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 16153425bc38SStefano Zampini ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 1616edf7251bSStefano Zampini ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr); 16173425bc38SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 161823ee1639SBarry Smith ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr); 16193425bc38SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 16203425bc38SStefano Zampini /* return pointers for objects created */ 16213425bc38SStefano Zampini *fetidp_mat=newmat; 16223425bc38SStefano Zampini *fetidp_pc=newpc; 16233425bc38SStefano Zampini PetscFunctionReturn(0); 16243425bc38SStefano Zampini } 16251e6b0712SBarry Smith 16263425bc38SStefano Zampini #undef __FUNCT__ 16273425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 16283425bc38SStefano Zampini /*@ 162928509bceSStefano Zampini PCBDDCCreateFETIDPOperators - Create operators for FETIDP 16303425bc38SStefano Zampini 16313425bc38SStefano Zampini Collective 16323425bc38SStefano Zampini 16333425bc38SStefano Zampini Input Parameters: 163428509bceSStefano Zampini + pc - the BDDC preconditioning context already setup 163528509bceSStefano Zampini 163628509bceSStefano Zampini Output Parameters: 163728509bceSStefano Zampini . fetidp_mat - shell FETIDP matrix object 163828509bceSStefano Zampini . fetidp_pc - shell Dirichlet preconditioner for FETIDP matrix 163928509bceSStefano Zampini 164028509bceSStefano Zampini Options Database Keys: 164128509bceSStefano Zampini - -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers 16423425bc38SStefano Zampini 16433425bc38SStefano Zampini Level: developer 16443425bc38SStefano Zampini 16453425bc38SStefano Zampini Notes: 164628509bceSStefano Zampini Currently the only operation provided for FETIDP matrix is MatMult 16473425bc38SStefano Zampini 16483425bc38SStefano Zampini .seealso: PCBDDC 16493425bc38SStefano Zampini @*/ 16503425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 16513425bc38SStefano Zampini { 16523425bc38SStefano Zampini PetscErrorCode ierr; 16533425bc38SStefano Zampini 16543425bc38SStefano Zampini PetscFunctionBegin; 16553425bc38SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16563425bc38SStefano Zampini if (pc->setupcalled) { 1657516d51a7SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1658f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 16593425bc38SStefano Zampini PetscFunctionReturn(0); 16603425bc38SStefano Zampini } 16610c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 1662da1bb401SStefano Zampini /*MC 1663da1bb401SStefano Zampini PCBDDC - Balancing Domain Decomposition by Constraints. 16640c7d97c5SJed Brown 166528509bceSStefano Zampini An implementation of the BDDC preconditioner based on 166628509bceSStefano Zampini 166728509bceSStefano Zampini .vb 166828509bceSStefano Zampini [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 166928509bceSStefano 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 167028509bceSStefano Zampini [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 167128509bceSStefano Zampini .ve 167228509bceSStefano Zampini 167328509bceSStefano Zampini The matrix to be preconditioned (Pmat) must be of type MATIS. 167428509bceSStefano Zampini 1675b6fdb6dfSStefano Zampini Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 167628509bceSStefano Zampini 167728509bceSStefano Zampini It also works with unsymmetric and indefinite problems. 167828509bceSStefano Zampini 1679b6fdb6dfSStefano 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. 1680b6fdb6dfSStefano Zampini 168128509bceSStefano Zampini Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace 168228509bceSStefano Zampini 168328509bceSStefano 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 168428509bceSStefano Zampini 168528509bceSStefano Zampini Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace. 168628509bceSStefano Zampini 1687b6fdb6dfSStefano 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. 168828509bceSStefano Zampini 168928509bceSStefano Zampini The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object. 169028509bceSStefano Zampini 1691da1bb401SStefano Zampini Options Database Keys: 169228509bceSStefano Zampini 169328509bceSStefano Zampini . -pc_bddc_use_vertices <1> - use or not vertices in primal space 169428509bceSStefano Zampini . -pc_bddc_use_edges <1> - use or not edges in primal space 1695b6fdb6dfSStefano Zampini . -pc_bddc_use_faces <0> - use or not faces in primal space 169628509bceSStefano Zampini . -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only) 169728509bceSStefano Zampini . -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested 169828509bceSStefano Zampini . -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1]) 169928509bceSStefano Zampini . -pc_bddc_levels <0> - maximum number of levels for multilevel 170028509bceSStefano Zampini . -pc_bddc_coarsening_ratio - H/h ratio at the coarser level 170128509bceSStefano Zampini - -pc_bddc_check_level <0> - set verbosity level of debugging output 170228509bceSStefano Zampini 170328509bceSStefano Zampini Options for Dirichlet, Neumann or coarse solver can be set with 170428509bceSStefano Zampini .vb 170528509bceSStefano Zampini -pc_bddc_dirichlet_ 170628509bceSStefano Zampini -pc_bddc_neumann_ 170728509bceSStefano Zampini -pc_bddc_coarse_ 170828509bceSStefano Zampini .ve 170928509bceSStefano Zampini e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg 171028509bceSStefano Zampini 171128509bceSStefano Zampini When using a multilevel approach, solvers' options at the N-th level can be specified as 171228509bceSStefano Zampini .vb 1713312be037SStefano Zampini -pc_bddc_dirichlet_lN_ 1714312be037SStefano Zampini -pc_bddc_neumann_lN_ 1715312be037SStefano Zampini -pc_bddc_coarse_lN_ 171628509bceSStefano Zampini .ve 1717312be037SStefano Zampini Note that level number ranges from the finest 0 to the coarsest N. 1718da1bb401SStefano Zampini 1719da1bb401SStefano Zampini Level: intermediate 1720da1bb401SStefano Zampini 1721b6fdb6dfSStefano Zampini Developer notes: 1722da1bb401SStefano Zampini 172328509bceSStefano Zampini New deluxe scaling operator will be available soon. 1724da1bb401SStefano Zampini 1725da1bb401SStefano Zampini Contributed by Stefano Zampini 1726da1bb401SStefano Zampini 1727da1bb401SStefano Zampini .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1728da1bb401SStefano Zampini M*/ 1729b2573a8aSBarry Smith 1730da1bb401SStefano Zampini #undef __FUNCT__ 1731da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC" 17328cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1733da1bb401SStefano Zampini { 1734da1bb401SStefano Zampini PetscErrorCode ierr; 1735da1bb401SStefano Zampini PC_BDDC *pcbddc; 1736da1bb401SStefano Zampini 1737da1bb401SStefano Zampini PetscFunctionBegin; 1738da1bb401SStefano Zampini /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1739b00a9115SJed Brown ierr = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr); 1740da1bb401SStefano Zampini pc->data = (void*)pcbddc; 1741da1bb401SStefano Zampini 1742da1bb401SStefano Zampini /* create PCIS data structure */ 1743da1bb401SStefano Zampini ierr = PCISCreate(pc);CHKERRQ(ierr); 1744da1bb401SStefano Zampini 174547d04d0dSStefano Zampini /* BDDC customization */ 174608a5cf49SStefano Zampini pcbddc->use_local_adj = PETSC_TRUE; 174747d04d0dSStefano Zampini pcbddc->use_vertices = PETSC_TRUE; 174847d04d0dSStefano Zampini pcbddc->use_edges = PETSC_TRUE; 174947d04d0dSStefano Zampini pcbddc->use_faces = PETSC_FALSE; 175047d04d0dSStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 175147d04d0dSStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 175247d04d0dSStefano Zampini pcbddc->switch_static = PETSC_FALSE; 175347d04d0dSStefano Zampini pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 175447d04d0dSStefano Zampini pcbddc->dbg_flag = 0; 1755*b9d89cd5SStefano Zampini /* private */ 1756*b9d89cd5SStefano Zampini pcbddc->issym = PETSC_FALSE; 175739e2fb2aSStefano Zampini pcbddc->BtoNmap = 0; 1758727cdba6SStefano Zampini pcbddc->local_primal_size = 0; 1759e9189074SStefano Zampini pcbddc->n_vertices = 0; 1760e9189074SStefano Zampini pcbddc->n_actual_vertices = 0; 1761e9189074SStefano Zampini pcbddc->n_constraints = 0; 1762727cdba6SStefano Zampini pcbddc->primal_indices_local_idxs = 0; 1763fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_FALSE; 176468457ee5SStefano Zampini pcbddc->coarse_size = -1; 1765f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1766727cdba6SStefano Zampini pcbddc->new_primal_space_local = PETSC_FALSE; 1767f4ddd8eeSStefano Zampini pcbddc->global_primal_indices = 0; 1768f4ddd8eeSStefano Zampini pcbddc->onearnullspace = 0; 1769f4ddd8eeSStefano Zampini pcbddc->onearnullvecs_state = 0; 1770674ae819SStefano Zampini pcbddc->user_primal_vertices = 0; 17710bdf917eSStefano Zampini pcbddc->NullSpace = 0; 17723972b0daSStefano Zampini pcbddc->temp_solution = 0; 1773534831adSStefano Zampini pcbddc->original_rhs = 0; 1774534831adSStefano Zampini pcbddc->local_mat = 0; 1775534831adSStefano Zampini pcbddc->ChangeOfBasisMatrix = 0; 1776b9b85e73SStefano Zampini pcbddc->user_ChangeOfBasisMatrix = 0; 1777da1bb401SStefano Zampini pcbddc->coarse_vec = 0; 1778da1bb401SStefano Zampini pcbddc->coarse_rhs = 0; 1779da1bb401SStefano Zampini pcbddc->coarse_ksp = 0; 1780da1bb401SStefano Zampini pcbddc->coarse_phi_B = 0; 1781da1bb401SStefano Zampini pcbddc->coarse_phi_D = 0; 178215aaf578SStefano Zampini pcbddc->coarse_psi_B = 0; 178315aaf578SStefano Zampini pcbddc->coarse_psi_D = 0; 1784da1bb401SStefano Zampini pcbddc->vec1_P = 0; 1785da1bb401SStefano Zampini pcbddc->vec1_R = 0; 1786da1bb401SStefano Zampini pcbddc->vec2_R = 0; 1787da1bb401SStefano Zampini pcbddc->local_auxmat1 = 0; 1788da1bb401SStefano Zampini pcbddc->local_auxmat2 = 0; 1789da1bb401SStefano Zampini pcbddc->R_to_B = 0; 1790da1bb401SStefano Zampini pcbddc->R_to_D = 0; 1791da1bb401SStefano Zampini pcbddc->ksp_D = 0; 1792da1bb401SStefano Zampini pcbddc->ksp_R = 0; 1793da1bb401SStefano Zampini pcbddc->NeumannBoundaries = 0; 1794785d1243SStefano Zampini pcbddc->NeumannBoundariesLocal = 0; 1795785d1243SStefano Zampini pcbddc->DirichletBoundaries = 0; 1796785d1243SStefano Zampini pcbddc->DirichletBoundariesLocal = 0; 179760d44989SStefano Zampini pcbddc->user_provided_isfordofs = PETSC_FALSE; 179860d44989SStefano Zampini pcbddc->n_ISForDofs = 0; 179963602bcaSStefano Zampini pcbddc->n_ISForDofsLocal = 0; 1800da1bb401SStefano Zampini pcbddc->ISForDofs = 0; 180163602bcaSStefano Zampini pcbddc->ISForDofsLocal = 0; 1802da1bb401SStefano Zampini pcbddc->ConstraintMatrix = 0; 180385c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 180447d04d0dSStefano Zampini pcbddc->coarse_loc_to_glob = 0; 180547d04d0dSStefano Zampini pcbddc->coarsening_ratio = 8; 18064fad6a16SStefano Zampini pcbddc->current_level = 0; 18072b510759SStefano Zampini pcbddc->max_levels = 0; 1808323d395dSStefano Zampini pcbddc->use_coarse_estimates = PETSC_FALSE; 1809f3bde8b3SStefano Zampini pcbddc->coarse_subassembling = 0; 1810323d395dSStefano Zampini pcbddc->coarse_subassembling_init = 0; 1811da1bb401SStefano Zampini 1812674ae819SStefano Zampini /* create local graph structure */ 1813674ae819SStefano Zampini ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1814674ae819SStefano Zampini 1815674ae819SStefano Zampini /* scaling */ 1816674ae819SStefano Zampini pcbddc->use_deluxe_scaling = PETSC_FALSE; 181708a5cf49SStefano Zampini pcbddc->deluxe_ctx = 0; 1818674ae819SStefano Zampini pcbddc->work_scaling = 0; 1819da1bb401SStefano Zampini 1820da1bb401SStefano Zampini /* function pointers */ 1821da1bb401SStefano Zampini pc->ops->apply = PCApply_BDDC; 182293bd9ae7SStefano Zampini pc->ops->applytranspose = PCApplyTranspose_BDDC; 1823da1bb401SStefano Zampini pc->ops->setup = PCSetUp_BDDC; 1824da1bb401SStefano Zampini pc->ops->destroy = PCDestroy_BDDC; 1825da1bb401SStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1826da1bb401SStefano Zampini pc->ops->view = 0; 1827da1bb401SStefano Zampini pc->ops->applyrichardson = 0; 1828da1bb401SStefano Zampini pc->ops->applysymmetricleft = 0; 1829da1bb401SStefano Zampini pc->ops->applysymmetricright = 0; 1830534831adSStefano Zampini pc->ops->presolve = PCPreSolve_BDDC; 1831534831adSStefano Zampini pc->ops->postsolve = PCPostSolve_BDDC; 1832da1bb401SStefano Zampini 1833da1bb401SStefano Zampini /* composing function */ 1834b9b85e73SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisLocalMat_C",PCBDDCSetChangeOfBasisLocalMat_BDDC);CHKERRQ(ierr); 1835674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1836bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 18372b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1838b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 18392b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1840bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1841bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 184282ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1843bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 184482ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1845bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 184682ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1847bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 184882ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1849bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 185063602bcaSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr); 1851bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1852bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1853bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1854bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1855da1bb401SStefano Zampini PetscFunctionReturn(0); 1856da1bb401SStefano Zampini } 18573425bc38SStefano Zampini 1858