153cdbc3dSStefano Zampini /* TODOLIST 2eb97c9d2SStefano Zampini 3eb97c9d2SStefano Zampini ConstraintsSetup 4eb97c9d2SStefano Zampini - assure same constraints between neighbours by sorting vals by global index before SVD! 5eb97c9d2SStefano Zampini - tolerances for constraints as an option (take care of single precision!) 6eb97c9d2SStefano Zampini - Allow different constraints customizations among different linear solves (requires also reset/destroy of ksp_R and coarse_ksp) 7eb97c9d2SStefano Zampini - MAT_IGNORE_ZERO_ENTRIES for Constraints Matrix 8eb97c9d2SStefano Zampini 9eb97c9d2SStefano Zampini Solvers 10eb97c9d2SStefano Zampini - Try to reduce the work when reusing the solvers 11eb97c9d2SStefano Zampini - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers) 12eb97c9d2SStefano Zampini - reuse already allocated coarse matrix if possible 13eb97c9d2SStefano Zampini - Propagate ksp prefixes for solvers to mat objects? 14eb97c9d2SStefano Zampini - Propagate nearnullspace info among levels 15eb97c9d2SStefano Zampini 16eb97c9d2SStefano Zampini User interface 17eb97c9d2SStefano Zampini - Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet) 18eb97c9d2SStefano Zampini - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access) 19eb97c9d2SStefano Zampini - Provide PCApplyTranpose_BDDC 20eb97c9d2SStefano Zampini - DofSplitting and DM attached to pc? 21eb97c9d2SStefano Zampini 22eb97c9d2SStefano Zampini Debugging output 23eb97c9d2SStefano Zampini - Better management of verbosity levels of debugging output 24eb97c9d2SStefano Zampini - Crashes on some architecture -> call SynchronizedAllow before every SynchronizedPrintf 25eb97c9d2SStefano Zampini 26eb97c9d2SStefano Zampini Build 27eb97c9d2SStefano Zampini - make runexe59 28eb97c9d2SStefano Zampini 29eb97c9d2SStefano Zampini Extra 30eb97c9d2SStefano Zampini - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)? 31eb97c9d2SStefano Zampini - Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver? 32eb97c9d2SStefano Zampini - add support for computing h,H and related using coordinates? 33c0b83709SStefano Zampini - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap) 34eb97c9d2SStefano Zampini - Better management in PCIS code 35eb97c9d2SStefano Zampini - BDDC with MG framework? 36eb97c9d2SStefano Zampini 37eb97c9d2SStefano Zampini FETIDP 38eb97c9d2SStefano Zampini - Move FETIDP code to its own classes 39eb97c9d2SStefano Zampini 40eb97c9d2SStefano Zampini MATIS related operations contained in BDDC code 41eb97c9d2SStefano Zampini - Provide general case for subassembling 42eb97c9d2SStefano Zampini - Preallocation routines in MatConvert_IS_AIJ 43eb97c9d2SStefano Zampini 4453cdbc3dSStefano Zampini */ 450c7d97c5SJed Brown 4653cdbc3dSStefano Zampini /* ---------------------------------------------------------------------------------------------------------------------------------------------- 470c7d97c5SJed Brown Implementation of BDDC preconditioner based on: 480c7d97c5SJed Brown C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 4953cdbc3dSStefano Zampini ---------------------------------------------------------------------------------------------------------------------------------------------- */ 5053cdbc3dSStefano Zampini 51674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 52674ae819SStefano Zampini #include "bddcprivate.h" 533b03a366Sstefano_zampini #include <petscblaslapack.h> 54674ae819SStefano Zampini 550c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 560c7d97c5SJed Brown #undef __FUNCT__ 570c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC" 580c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc) 590c7d97c5SJed Brown { 600c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 610c7d97c5SJed Brown PetscErrorCode ierr; 620c7d97c5SJed Brown 630c7d97c5SJed Brown PetscFunctionBegin; 640c7d97c5SJed Brown ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 658eeda7d8SStefano Zampini /* Verbose debugging */ 668eeda7d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 678eeda7d8SStefano Zampini /* Primal space cumstomization */ 688eeda7d8SStefano 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); 698eeda7d8SStefano 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); 708eeda7d8SStefano 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); 718eeda7d8SStefano Zampini /* Change of basis */ 728eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr); 738eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr); 74674ae819SStefano Zampini if (!pcbddc->use_change_of_basis) { 75674ae819SStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 76674ae819SStefano Zampini } 778eeda7d8SStefano Zampini /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 788eeda7d8SStefano 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); 790298fd71SBarry Smith ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 802b510759SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 81674ae819SStefano 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); 820c7d97c5SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 830c7d97c5SJed Brown PetscFunctionReturn(0); 840c7d97c5SJed Brown } 850c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 86674ae819SStefano Zampini #undef __FUNCT__ 87674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 88674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 89674ae819SStefano Zampini { 90674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 91674ae819SStefano Zampini PetscErrorCode ierr; 921e6b0712SBarry Smith 93674ae819SStefano Zampini PetscFunctionBegin; 94674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 95674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 96674ae819SStefano Zampini pcbddc->user_primal_vertices = PrimalVertices; 97674ae819SStefano Zampini PetscFunctionReturn(0); 98674ae819SStefano Zampini } 99674ae819SStefano Zampini #undef __FUNCT__ 100674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 101674ae819SStefano Zampini /*@ 10228509bceSStefano Zampini PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 103674ae819SStefano Zampini 104674ae819SStefano Zampini Not collective 105674ae819SStefano Zampini 106674ae819SStefano Zampini Input Parameters: 107674ae819SStefano Zampini + pc - the preconditioning context 10828509bceSStefano Zampini - PrimalVertices - index set of primal vertices in local numbering 109674ae819SStefano Zampini 110674ae819SStefano Zampini Level: intermediate 111674ae819SStefano Zampini 112674ae819SStefano Zampini Notes: 113674ae819SStefano Zampini 114674ae819SStefano Zampini .seealso: PCBDDC 115674ae819SStefano Zampini @*/ 116674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 117674ae819SStefano Zampini { 118674ae819SStefano Zampini PetscErrorCode ierr; 119674ae819SStefano Zampini 120674ae819SStefano Zampini PetscFunctionBegin; 121674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 122674ae819SStefano Zampini PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 123674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 124674ae819SStefano Zampini PetscFunctionReturn(0); 125674ae819SStefano Zampini } 126674ae819SStefano Zampini /* -------------------------------------------------------------------------- */ 1270c7d97c5SJed Brown #undef __FUNCT__ 1284fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 1294fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 1304fad6a16SStefano Zampini { 1314fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1324fad6a16SStefano Zampini 1334fad6a16SStefano Zampini PetscFunctionBegin; 1344fad6a16SStefano Zampini pcbddc->coarsening_ratio = k; 1354fad6a16SStefano Zampini PetscFunctionReturn(0); 1364fad6a16SStefano Zampini } 1371e6b0712SBarry Smith 1384fad6a16SStefano Zampini #undef __FUNCT__ 1394fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio" 1404fad6a16SStefano Zampini /*@ 14128509bceSStefano Zampini PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 1424fad6a16SStefano Zampini 1434fad6a16SStefano Zampini Logically collective on PC 1444fad6a16SStefano Zampini 1454fad6a16SStefano Zampini Input Parameters: 1464fad6a16SStefano Zampini + pc - the preconditioning context 14728509bceSStefano Zampini - k - coarsening ratio (H/h at the coarser level) 1484fad6a16SStefano Zampini 14928509bceSStefano Zampini Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 1504fad6a16SStefano Zampini 1514fad6a16SStefano Zampini Level: intermediate 1524fad6a16SStefano Zampini 1534fad6a16SStefano Zampini Notes: 1544fad6a16SStefano Zampini 1554fad6a16SStefano Zampini .seealso: PCBDDC 1564fad6a16SStefano Zampini @*/ 1574fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 1584fad6a16SStefano Zampini { 1594fad6a16SStefano Zampini PetscErrorCode ierr; 1604fad6a16SStefano Zampini 1614fad6a16SStefano Zampini PetscFunctionBegin; 1624fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1632b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,k,2); 1644fad6a16SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 1654fad6a16SStefano Zampini PetscFunctionReturn(0); 1664fad6a16SStefano Zampini } 1672b510759SStefano Zampini 168b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 1692b510759SStefano Zampini #undef __FUNCT__ 170b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 171b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 172b8ffe317SStefano Zampini { 173b8ffe317SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 174b8ffe317SStefano Zampini 175b8ffe317SStefano Zampini PetscFunctionBegin; 17685c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = flg; 177b8ffe317SStefano Zampini PetscFunctionReturn(0); 178b8ffe317SStefano Zampini } 179b8ffe317SStefano Zampini 180b8ffe317SStefano Zampini #undef __FUNCT__ 181b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 182b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 1832b510759SStefano Zampini { 1842b510759SStefano Zampini PetscErrorCode ierr; 1852b510759SStefano Zampini 1862b510759SStefano Zampini PetscFunctionBegin; 1872b510759SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 188b8ffe317SStefano Zampini PetscValidLogicalCollectiveBool(pc,flg,2); 189b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1902b510759SStefano Zampini PetscFunctionReturn(0); 1912b510759SStefano Zampini } 1921e6b0712SBarry Smith 1934fad6a16SStefano Zampini #undef __FUNCT__ 1942b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC" 1952b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 1964fad6a16SStefano Zampini { 1974fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1984fad6a16SStefano Zampini 1994fad6a16SStefano Zampini PetscFunctionBegin; 2002b510759SStefano Zampini pcbddc->current_level = level; 2014fad6a16SStefano Zampini PetscFunctionReturn(0); 2024fad6a16SStefano Zampini } 2031e6b0712SBarry Smith 2044fad6a16SStefano Zampini #undef __FUNCT__ 205b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel" 206b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 207b8ffe317SStefano Zampini { 208b8ffe317SStefano Zampini PetscErrorCode ierr; 209b8ffe317SStefano Zampini 210b8ffe317SStefano Zampini PetscFunctionBegin; 211b8ffe317SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 212b8ffe317SStefano Zampini PetscValidLogicalCollectiveInt(pc,level,2); 213b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 214b8ffe317SStefano Zampini PetscFunctionReturn(0); 215b8ffe317SStefano Zampini } 216b8ffe317SStefano Zampini 217b8ffe317SStefano Zampini #undef __FUNCT__ 2182b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC" 2192b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 2202b510759SStefano Zampini { 2212b510759SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2222b510759SStefano Zampini 2232b510759SStefano Zampini PetscFunctionBegin; 2242b510759SStefano Zampini pcbddc->max_levels = levels; 2252b510759SStefano Zampini PetscFunctionReturn(0); 2262b510759SStefano Zampini } 2272b510759SStefano Zampini 2282b510759SStefano Zampini #undef __FUNCT__ 2292b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels" 2304fad6a16SStefano Zampini /*@ 23128509bceSStefano Zampini PCBDDCSetLevels - Sets the maximum number of levels for multilevel 2324fad6a16SStefano Zampini 2334fad6a16SStefano Zampini Logically collective on PC 2344fad6a16SStefano Zampini 2354fad6a16SStefano Zampini Input Parameters: 2364fad6a16SStefano Zampini + pc - the preconditioning context 23728509bceSStefano Zampini - levels - the maximum number of levels (max 9) 2384fad6a16SStefano Zampini 23928509bceSStefano Zampini Default value is 0, i.e. traditional one-level BDDC 2404fad6a16SStefano Zampini 2414fad6a16SStefano Zampini Level: intermediate 2424fad6a16SStefano Zampini 2434fad6a16SStefano Zampini Notes: 2444fad6a16SStefano Zampini 2454fad6a16SStefano Zampini .seealso: PCBDDC 2464fad6a16SStefano Zampini @*/ 2472b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 2484fad6a16SStefano Zampini { 2494fad6a16SStefano Zampini PetscErrorCode ierr; 2504fad6a16SStefano Zampini 2514fad6a16SStefano Zampini PetscFunctionBegin; 2524fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2532b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,levels,2); 2542b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 2554fad6a16SStefano Zampini PetscFunctionReturn(0); 2564fad6a16SStefano Zampini } 2574fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */ 2581e6b0712SBarry Smith 2594fad6a16SStefano Zampini #undef __FUNCT__ 2600bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 2610bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 2620bdf917eSStefano Zampini { 2630bdf917eSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2640bdf917eSStefano Zampini PetscErrorCode ierr; 2650bdf917eSStefano Zampini 2660bdf917eSStefano Zampini PetscFunctionBegin; 2670bdf917eSStefano Zampini ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 2680bdf917eSStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 2690bdf917eSStefano Zampini pcbddc->NullSpace = NullSpace; 2700bdf917eSStefano Zampini PetscFunctionReturn(0); 2710bdf917eSStefano Zampini } 2721e6b0712SBarry Smith 2730bdf917eSStefano Zampini #undef __FUNCT__ 2740bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace" 2750bdf917eSStefano Zampini /*@ 27628509bceSStefano Zampini PCBDDCSetNullSpace - Set nullspace for BDDC operator 2770bdf917eSStefano Zampini 2780bdf917eSStefano Zampini Logically collective on PC and MatNullSpace 2790bdf917eSStefano Zampini 2800bdf917eSStefano Zampini Input Parameters: 2810bdf917eSStefano Zampini + pc - the preconditioning context 28228509bceSStefano Zampini - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 2830bdf917eSStefano Zampini 2840bdf917eSStefano Zampini Level: intermediate 2850bdf917eSStefano Zampini 2860bdf917eSStefano Zampini Notes: 2870bdf917eSStefano Zampini 2880bdf917eSStefano Zampini .seealso: PCBDDC 2890bdf917eSStefano Zampini @*/ 2900bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 2910bdf917eSStefano Zampini { 2920bdf917eSStefano Zampini PetscErrorCode ierr; 2930bdf917eSStefano Zampini 2940bdf917eSStefano Zampini PetscFunctionBegin; 2950bdf917eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 296674ae819SStefano Zampini PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 2972b510759SStefano Zampini PetscCheckSameComm(pc,1,NullSpace,2); 2980bdf917eSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 2990bdf917eSStefano Zampini PetscFunctionReturn(0); 3000bdf917eSStefano Zampini } 3010bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */ 3021e6b0712SBarry Smith 3030bdf917eSStefano Zampini #undef __FUNCT__ 3043b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 3053b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 3063b03a366Sstefano_zampini { 3073b03a366Sstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3083b03a366Sstefano_zampini PetscErrorCode ierr; 3093b03a366Sstefano_zampini 3103b03a366Sstefano_zampini PetscFunctionBegin; 3113b03a366Sstefano_zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 31236e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 31336e030ebSStefano Zampini pcbddc->DirichletBoundaries=DirichletBoundaries; 3143b03a366Sstefano_zampini PetscFunctionReturn(0); 3153b03a366Sstefano_zampini } 3161e6b0712SBarry Smith 3173b03a366Sstefano_zampini #undef __FUNCT__ 3183b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 3193b03a366Sstefano_zampini /*@ 32028509bceSStefano Zampini PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 3213b03a366Sstefano_zampini 3223b03a366Sstefano_zampini Not collective 3233b03a366Sstefano_zampini 3243b03a366Sstefano_zampini Input Parameters: 3253b03a366Sstefano_zampini + pc - the preconditioning context 32628509bceSStefano Zampini - DirichletBoundaries - sequential IS defining the subdomain part of Dirichlet boundaries (in local ordering) 3273b03a366Sstefano_zampini 3283b03a366Sstefano_zampini Level: intermediate 3293b03a366Sstefano_zampini 3303b03a366Sstefano_zampini Notes: 3313b03a366Sstefano_zampini 3323b03a366Sstefano_zampini .seealso: PCBDDC 3333b03a366Sstefano_zampini @*/ 3343b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 3353b03a366Sstefano_zampini { 3363b03a366Sstefano_zampini PetscErrorCode ierr; 3373b03a366Sstefano_zampini 3383b03a366Sstefano_zampini PetscFunctionBegin; 3393b03a366Sstefano_zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 340674ae819SStefano Zampini PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 3413b03a366Sstefano_zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 3423b03a366Sstefano_zampini PetscFunctionReturn(0); 3433b03a366Sstefano_zampini } 3443b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */ 3451e6b0712SBarry Smith 3463b03a366Sstefano_zampini #undef __FUNCT__ 3470c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 34853cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 3490c7d97c5SJed Brown { 3500c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 35153cdbc3dSStefano Zampini PetscErrorCode ierr; 3520c7d97c5SJed Brown 3530c7d97c5SJed Brown PetscFunctionBegin; 35453cdbc3dSStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 35536e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 35636e030ebSStefano Zampini pcbddc->NeumannBoundaries=NeumannBoundaries; 3570c7d97c5SJed Brown PetscFunctionReturn(0); 3580c7d97c5SJed Brown } 3591e6b0712SBarry Smith 3600c7d97c5SJed Brown #undef __FUNCT__ 3610c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 36257527edcSJed Brown /*@ 36328509bceSStefano Zampini PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 36457527edcSJed Brown 3659c0446d6SStefano Zampini Not collective 36657527edcSJed Brown 36757527edcSJed Brown Input Parameters: 36857527edcSJed Brown + pc - the preconditioning context 36928509bceSStefano Zampini - NeumannBoundaries - sequential IS defining the subdomain part of Neumann boundaries (in local ordering) 37057527edcSJed Brown 37157527edcSJed Brown Level: intermediate 37257527edcSJed Brown 37357527edcSJed Brown Notes: 37457527edcSJed Brown 37557527edcSJed Brown .seealso: PCBDDC 37657527edcSJed Brown @*/ 37753cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 3780c7d97c5SJed Brown { 3790c7d97c5SJed Brown PetscErrorCode ierr; 3800c7d97c5SJed Brown 3810c7d97c5SJed Brown PetscFunctionBegin; 3820c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 383674ae819SStefano Zampini PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 38453cdbc3dSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 38553cdbc3dSStefano Zampini PetscFunctionReturn(0); 38653cdbc3dSStefano Zampini } 38753cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 3881e6b0712SBarry Smith 38953cdbc3dSStefano Zampini #undef __FUNCT__ 390da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 391da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 392da1bb401SStefano Zampini { 393da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 394da1bb401SStefano Zampini 395da1bb401SStefano Zampini PetscFunctionBegin; 396da1bb401SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundaries; 397da1bb401SStefano Zampini PetscFunctionReturn(0); 398da1bb401SStefano Zampini } 3991e6b0712SBarry Smith 400da1bb401SStefano Zampini #undef __FUNCT__ 401da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 402da1bb401SStefano Zampini /*@ 40328509bceSStefano Zampini PCBDDCGetDirichletBoundaries - Get IS for local Dirichlet boundaries 404da1bb401SStefano Zampini 405da1bb401SStefano Zampini Not collective 406da1bb401SStefano Zampini 407da1bb401SStefano Zampini Input Parameters: 40828509bceSStefano Zampini . pc - the preconditioning context 409da1bb401SStefano Zampini 410da1bb401SStefano Zampini Output Parameters: 41128509bceSStefano Zampini . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 412da1bb401SStefano Zampini 413da1bb401SStefano Zampini Level: intermediate 414da1bb401SStefano Zampini 415da1bb401SStefano Zampini Notes: 416da1bb401SStefano Zampini 417da1bb401SStefano Zampini .seealso: PCBDDC 418da1bb401SStefano Zampini @*/ 419da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 420da1bb401SStefano Zampini { 421da1bb401SStefano Zampini PetscErrorCode ierr; 422da1bb401SStefano Zampini 423da1bb401SStefano Zampini PetscFunctionBegin; 424da1bb401SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 425da1bb401SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 426da1bb401SStefano Zampini PetscFunctionReturn(0); 427da1bb401SStefano Zampini } 428da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 4291e6b0712SBarry Smith 430da1bb401SStefano Zampini #undef __FUNCT__ 43153cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 43253cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 43353cdbc3dSStefano Zampini { 43453cdbc3dSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 43553cdbc3dSStefano Zampini 43653cdbc3dSStefano Zampini PetscFunctionBegin; 43753cdbc3dSStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundaries; 43853cdbc3dSStefano Zampini PetscFunctionReturn(0); 43953cdbc3dSStefano Zampini } 4401e6b0712SBarry Smith 44153cdbc3dSStefano Zampini #undef __FUNCT__ 44253cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 44353cdbc3dSStefano Zampini /*@ 44428509bceSStefano Zampini PCBDDCGetNeumannBoundaries - Get IS for local Neumann boundaries 44553cdbc3dSStefano Zampini 4469c0446d6SStefano Zampini Not collective 44753cdbc3dSStefano Zampini 44853cdbc3dSStefano Zampini Input Parameters: 44928509bceSStefano Zampini . pc - the preconditioning context 45053cdbc3dSStefano Zampini 45153cdbc3dSStefano Zampini Output Parameters: 45228509bceSStefano Zampini . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 45353cdbc3dSStefano Zampini 45453cdbc3dSStefano Zampini Level: intermediate 45553cdbc3dSStefano Zampini 45653cdbc3dSStefano Zampini Notes: 45753cdbc3dSStefano Zampini 45853cdbc3dSStefano Zampini .seealso: PCBDDC 45953cdbc3dSStefano Zampini @*/ 46053cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 46153cdbc3dSStefano Zampini { 46253cdbc3dSStefano Zampini PetscErrorCode ierr; 46353cdbc3dSStefano Zampini 46453cdbc3dSStefano Zampini PetscFunctionBegin; 46553cdbc3dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 46653cdbc3dSStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 4670c7d97c5SJed Brown PetscFunctionReturn(0); 4680c7d97c5SJed Brown } 46936e030ebSStefano Zampini /* -------------------------------------------------------------------------- */ 4701e6b0712SBarry Smith 47136e030ebSStefano Zampini #undef __FUNCT__ 472da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 4731a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 47436e030ebSStefano Zampini { 47536e030ebSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 476da1bb401SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 477da1bb401SStefano Zampini PetscErrorCode ierr; 47836e030ebSStefano Zampini 47936e030ebSStefano Zampini PetscFunctionBegin; 480674ae819SStefano Zampini /* free old CSR */ 481674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 482fb223d50SStefano Zampini /* TODO: PCBDDCGraphSetAdjacency */ 483674ae819SStefano Zampini /* get CSR into graph structure */ 484da1bb401SStefano Zampini if (copymode == PETSC_COPY_VALUES) { 485674ae819SStefano Zampini ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr); 486674ae819SStefano Zampini ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr); 487674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 488674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 489da1bb401SStefano Zampini } else if (copymode == PETSC_OWN_POINTER) { 4901a83f524SJed Brown mat_graph->xadj = (PetscInt*)xadj; 4911a83f524SJed Brown mat_graph->adjncy = (PetscInt*)adjncy; 492674ae819SStefano Zampini } 493575ad6abSStefano Zampini mat_graph->nvtxs_csr = nvtxs; 49436e030ebSStefano Zampini PetscFunctionReturn(0); 49536e030ebSStefano Zampini } 4961e6b0712SBarry Smith 49736e030ebSStefano Zampini #undef __FUNCT__ 498da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 49936e030ebSStefano Zampini /*@ 50028509bceSStefano Zampini PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix 50136e030ebSStefano Zampini 50236e030ebSStefano Zampini Not collective 50336e030ebSStefano Zampini 50436e030ebSStefano Zampini Input Parameters: 50536e030ebSStefano Zampini + pc - the preconditioning context 50628509bceSStefano Zampini . nvtxs - number of local vertices of the graph (i.e., the local size of your problem) 50728509bceSStefano Zampini . xadj, adjncy - the CSR graph 50828509bceSStefano Zampini - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 50936e030ebSStefano Zampini 51036e030ebSStefano Zampini Level: intermediate 51136e030ebSStefano Zampini 51236e030ebSStefano Zampini Notes: 51336e030ebSStefano Zampini 51428509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode 51536e030ebSStefano Zampini @*/ 5161a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 51736e030ebSStefano Zampini { 518575ad6abSStefano Zampini void (*f)(void) = 0; 51936e030ebSStefano Zampini PetscErrorCode ierr; 52036e030ebSStefano Zampini 52136e030ebSStefano Zampini PetscFunctionBegin; 52236e030ebSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 523674ae819SStefano Zampini PetscValidIntPointer(xadj,3); 524674ae819SStefano Zampini PetscValidIntPointer(xadj,4); 525674ae819SStefano Zampini if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 526674ae819SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 527da1bb401SStefano Zampini } 52836e030ebSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 529575ad6abSStefano Zampini /* free arrays if PCBDDC is not the PC type */ 530575ad6abSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 531575ad6abSStefano Zampini if (!f && copymode == PETSC_OWN_POINTER) { 532575ad6abSStefano Zampini ierr = PetscFree(xadj);CHKERRQ(ierr); 533575ad6abSStefano Zampini ierr = PetscFree(adjncy);CHKERRQ(ierr); 53436e030ebSStefano Zampini } 53536e030ebSStefano Zampini PetscFunctionReturn(0); 53636e030ebSStefano Zampini } 5379c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */ 5381e6b0712SBarry Smith 5399c0446d6SStefano Zampini #undef __FUNCT__ 5409c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 5419c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 5429c0446d6SStefano Zampini { 5439c0446d6SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 5449c0446d6SStefano Zampini PetscInt i; 5459c0446d6SStefano Zampini PetscErrorCode ierr; 5469c0446d6SStefano Zampini 5479c0446d6SStefano Zampini PetscFunctionBegin; 548da1bb401SStefano Zampini /* Destroy ISes if they were already set */ 5499c0446d6SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 5509c0446d6SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 5519c0446d6SStefano Zampini } 552d11ae9bbSstefano_zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 553da1bb401SStefano Zampini /* allocate space then set */ 5549c0446d6SStefano Zampini ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr); 5559c0446d6SStefano Zampini for (i=0;i<n_is;i++) { 556da1bb401SStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 557da1bb401SStefano Zampini pcbddc->ISForDofs[i]=ISForDofs[i]; 5589c0446d6SStefano Zampini } 5599c0446d6SStefano Zampini pcbddc->n_ISForDofs=n_is; 56060d44989SStefano Zampini pcbddc->user_provided_isfordofs = PETSC_TRUE; 5619c0446d6SStefano Zampini PetscFunctionReturn(0); 5629c0446d6SStefano Zampini } 5631e6b0712SBarry Smith 5649c0446d6SStefano Zampini #undef __FUNCT__ 5659c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting" 5669c0446d6SStefano Zampini /*@ 56728509bceSStefano Zampini PCBDDCSetDofsSplitting - Set index sets defining fields of the local Neumann matrix 5689c0446d6SStefano Zampini 5699c0446d6SStefano Zampini Not collective 5709c0446d6SStefano Zampini 5719c0446d6SStefano Zampini Input Parameters: 5729c0446d6SStefano Zampini + pc - the preconditioning context 57328509bceSStefano Zampini - n_is - number of index sets defining the fields 57428509bceSStefano Zampini . ISForDofs - array of IS describing the fields 5759c0446d6SStefano Zampini 5769c0446d6SStefano Zampini Level: intermediate 5779c0446d6SStefano Zampini 5789c0446d6SStefano Zampini Notes: 5799c0446d6SStefano Zampini 5809c0446d6SStefano Zampini .seealso: PCBDDC 5819c0446d6SStefano Zampini @*/ 5829c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 5839c0446d6SStefano Zampini { 5842b510759SStefano Zampini PetscInt i; 5859c0446d6SStefano Zampini PetscErrorCode ierr; 5869c0446d6SStefano Zampini 5879c0446d6SStefano Zampini PetscFunctionBegin; 5889c0446d6SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5892b510759SStefano Zampini for (i=0;i<n_is;i++) { 5902b510759SStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2); 5912b510759SStefano Zampini } 5929c0446d6SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 5939c0446d6SStefano Zampini PetscFunctionReturn(0); 5949c0446d6SStefano Zampini } 595da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 596534831adSStefano Zampini #undef __FUNCT__ 597534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC" 598534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 599534831adSStefano Zampini /* 600534831adSStefano Zampini PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 601534831adSStefano Zampini guess if a transformation of basis approach has been selected. 6029c0446d6SStefano Zampini 603534831adSStefano Zampini Input Parameter: 604534831adSStefano Zampini + pc - the preconditioner contex 605534831adSStefano Zampini 606534831adSStefano Zampini Application Interface Routine: PCPreSolve() 607534831adSStefano Zampini 608534831adSStefano Zampini Notes: 609534831adSStefano Zampini The interface routine PCPreSolve() is not usually called directly by 610534831adSStefano Zampini the user, but instead is called by KSPSolve(). 611534831adSStefano Zampini */ 612534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 613534831adSStefano Zampini { 614534831adSStefano Zampini PetscErrorCode ierr; 615534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 616534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 617534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 618534831adSStefano Zampini Mat temp_mat; 6193972b0daSStefano Zampini IS dirIS; 6203972b0daSStefano Zampini PetscInt dirsize,i,*is_indices; 6213972b0daSStefano Zampini PetscScalar *array_x,*array_diagonal; 6223972b0daSStefano Zampini Vec used_vec; 62392e3dcfbSStefano Zampini PetscBool guess_nonzero,flg,bddc_has_dirichlet_boundaries; 624534831adSStefano Zampini 625534831adSStefano Zampini PetscFunctionBegin; 62685c4d303SStefano Zampini /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 62785c4d303SStefano Zampini if (ksp) { 62885c4d303SStefano Zampini PetscBool iscg; 62985c4d303SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 63085c4d303SStefano Zampini if (!iscg) { 63185c4d303SStefano Zampini ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 63285c4d303SStefano Zampini } 63385c4d303SStefano Zampini } 63485c4d303SStefano Zampini /* Creates parallel work vectors used in presolve */ 63562a6ff1dSStefano Zampini if (!pcbddc->original_rhs) { 63662a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 63762a6ff1dSStefano Zampini } 63862a6ff1dSStefano Zampini if (!pcbddc->temp_solution) { 63962a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 64062a6ff1dSStefano Zampini } 6413972b0daSStefano Zampini if (x) { 6423972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 6433972b0daSStefano Zampini used_vec = x; 6443972b0daSStefano Zampini } else { 6453972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 6463972b0daSStefano Zampini used_vec = pcbddc->temp_solution; 6473972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 6483972b0daSStefano Zampini } 6493972b0daSStefano Zampini /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 6503972b0daSStefano Zampini if (ksp) { 6513972b0daSStefano Zampini ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 6523972b0daSStefano Zampini if (!guess_nonzero) { 6533972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 6543972b0daSStefano Zampini } 6553972b0daSStefano Zampini } 6563308cffdSStefano Zampini 65792e3dcfbSStefano Zampini /* TODO: remove when Dirichlet boundaries will be shared */ 65892e3dcfbSStefano Zampini ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr); 65992e3dcfbSStefano Zampini flg = PETSC_FALSE; 66092e3dcfbSStefano Zampini if (dirIS) flg = PETSC_TRUE; 66192e3dcfbSStefano Zampini ierr = MPI_Allreduce(&flg,&bddc_has_dirichlet_boundaries,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 66292e3dcfbSStefano Zampini 6633972b0daSStefano Zampini /* store the original rhs */ 6643972b0daSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 6653972b0daSStefano Zampini 6663972b0daSStefano Zampini /* Take into account zeroed rows -> change rhs and store solution removed */ 66792e3dcfbSStefano Zampini if (rhs && bddc_has_dirichlet_boundaries) { 6683972b0daSStefano Zampini ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 6693972b0daSStefano Zampini ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 6703972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6713972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6723972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6733972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6743972b0daSStefano Zampini if (dirIS) { 6753972b0daSStefano Zampini ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr); 6763972b0daSStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 6773972b0daSStefano Zampini ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 6783972b0daSStefano Zampini ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 6792fa5cd67SKarl Rupp for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 6803972b0daSStefano Zampini ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 6813972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 6823972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 6833972b0daSStefano Zampini } 6843972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6853972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 686b76ba322SStefano Zampini 6873972b0daSStefano Zampini /* remove the computed solution from the rhs */ 6883972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 6893972b0daSStefano Zampini ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 6903972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 6913308cffdSStefano Zampini } 692b76ba322SStefano Zampini 693b76ba322SStefano Zampini /* store partially computed solution and set initial guess */ 6943972b0daSStefano Zampini if (x) { 6953972b0daSStefano Zampini ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 6963972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 69785c4d303SStefano Zampini if (pcbddc->use_exact_dirichlet_trick) { 698b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 699b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 700b76ba322SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 701b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 702b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 703b76ba322SStefano Zampini if (ksp) { 704b76ba322SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 705b76ba322SStefano Zampini } 706b76ba322SStefano Zampini } 7073972b0daSStefano Zampini } 708b76ba322SStefano Zampini 70992e3dcfbSStefano Zampini /* prepare MatMult and rhs for solver */ 710674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 711b76ba322SStefano Zampini /* swap pointers for local matrices */ 712b76ba322SStefano Zampini temp_mat = matis->A; 713b76ba322SStefano Zampini matis->A = pcbddc->local_mat; 714b76ba322SStefano Zampini pcbddc->local_mat = temp_mat; 71592e3dcfbSStefano Zampini if (rhs) { 716b76ba322SStefano Zampini /* Get local rhs and apply transformation of basis */ 717b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 718b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 719b76ba322SStefano Zampini /* from original basis to modified basis */ 720b76ba322SStefano Zampini ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 721b76ba322SStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 722b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 723b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 724674ae819SStefano Zampini } 72592e3dcfbSStefano Zampini } 72692e3dcfbSStefano Zampini 72792e3dcfbSStefano Zampini /* remove nullspace if present */ 7280bdf917eSStefano Zampini if (ksp && pcbddc->NullSpace) { 729d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 730d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 731b76ba322SStefano Zampini } 7320bdf917eSStefano Zampini ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 733534831adSStefano Zampini PetscFunctionReturn(0); 734534831adSStefano Zampini } 735534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 736534831adSStefano Zampini #undef __FUNCT__ 737534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC" 738534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 739534831adSStefano Zampini /* 740534831adSStefano Zampini PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 741534831adSStefano Zampini approach has been selected. Also, restores rhs to its original state. 742534831adSStefano Zampini 743534831adSStefano Zampini Input Parameter: 744534831adSStefano Zampini + pc - the preconditioner contex 745534831adSStefano Zampini 746534831adSStefano Zampini Application Interface Routine: PCPostSolve() 747534831adSStefano Zampini 748534831adSStefano Zampini Notes: 749534831adSStefano Zampini The interface routine PCPostSolve() is not usually called directly by 750534831adSStefano Zampini the user, but instead is called by KSPSolve(). 751534831adSStefano Zampini */ 752534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 753534831adSStefano Zampini { 754534831adSStefano Zampini PetscErrorCode ierr; 755534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 756534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 757534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 758534831adSStefano Zampini Mat temp_mat; 759534831adSStefano Zampini 760534831adSStefano Zampini PetscFunctionBegin; 761674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 762534831adSStefano Zampini /* swap pointers for local matrices */ 763534831adSStefano Zampini temp_mat = matis->A; 764534831adSStefano Zampini matis->A = pcbddc->local_mat; 765534831adSStefano Zampini pcbddc->local_mat = temp_mat; 7663425bc38SStefano Zampini } 7673308cffdSStefano Zampini if (pcbddc->use_change_of_basis && x) { 768534831adSStefano Zampini /* Get Local boundary and apply transformation of basis to solution vector */ 769534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 770534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 771534831adSStefano Zampini /* from modified basis to original basis */ 772534831adSStefano Zampini ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 773534831adSStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 774534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 775534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 776534831adSStefano Zampini } 7773972b0daSStefano Zampini /* add solution removed in presolve */ 7783425bc38SStefano Zampini if (x) { 7793425bc38SStefano Zampini ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 7803425bc38SStefano Zampini } 781fb223d50SStefano Zampini /* restore rhs to its original state */ 782fb223d50SStefano Zampini if (rhs) { 783fb223d50SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 784fb223d50SStefano Zampini } 785534831adSStefano Zampini PetscFunctionReturn(0); 786534831adSStefano Zampini } 787534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 78853cdbc3dSStefano Zampini #undef __FUNCT__ 78953cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC" 7900c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 7910c7d97c5SJed Brown /* 7920c7d97c5SJed Brown PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 7930c7d97c5SJed Brown by setting data structures and options. 7940c7d97c5SJed Brown 7950c7d97c5SJed Brown Input Parameter: 79653cdbc3dSStefano Zampini + pc - the preconditioner context 7970c7d97c5SJed Brown 7980c7d97c5SJed Brown Application Interface Routine: PCSetUp() 7990c7d97c5SJed Brown 8000c7d97c5SJed Brown Notes: 8010c7d97c5SJed Brown The interface routine PCSetUp() is not usually called directly by 8020c7d97c5SJed Brown the user, but instead is called by PCApply() if necessary. 8030c7d97c5SJed Brown */ 80453cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc) 8050c7d97c5SJed Brown { 8060c7d97c5SJed Brown PetscErrorCode ierr; 8070c7d97c5SJed Brown PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 808*f4ddd8eeSStefano Zampini MatNullSpace nearnullspace; 809*f4ddd8eeSStefano Zampini const Vec *nearnullvecs,*onearnullvecs; 810674ae819SStefano Zampini MatStructure flag; 811*f4ddd8eeSStefano Zampini PetscObjectState state; 812*f4ddd8eeSStefano Zampini PetscInt nnsp_size,onnsp_size; 813674ae819SStefano Zampini PetscBool computeis,computetopography,computesolvers; 814*f4ddd8eeSStefano Zampini PetscBool new_nearnullspace_provided,nnsp_has_cnst,onnsp_has_cnst; 8150c7d97c5SJed Brown 8160c7d97c5SJed Brown PetscFunctionBegin; 817*f4ddd8eeSStefano Zampini /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 818*f4ddd8eeSStefano Zampini /* PCIS does not support MatStructure flags different from SAME_PRECONDITIONER */ 8193b03a366Sstefano_zampini /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 820fcd409f5SStefano Zampini Also, BDDC directly build the Dirichlet problem */ 821*f4ddd8eeSStefano Zampini 822*f4ddd8eeSStefano Zampini /* split work */ 823674ae819SStefano Zampini if (pc->setupcalled) { 824674ae819SStefano Zampini computeis = PETSC_FALSE; 825674ae819SStefano Zampini ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr); 826674ae819SStefano Zampini if (flag == SAME_PRECONDITIONER) { 827*f4ddd8eeSStefano Zampini PetscFunctionReturn(0); 828674ae819SStefano Zampini } else if (flag == SAME_NONZERO_PATTERN) { 829674ae819SStefano Zampini computetopography = PETSC_FALSE; 830674ae819SStefano Zampini computesolvers = PETSC_TRUE; 831674ae819SStefano Zampini } else { /* DIFFERENT_NONZERO_PATTERN */ 832674ae819SStefano Zampini computetopography = PETSC_TRUE; 833674ae819SStefano Zampini computesolvers = PETSC_TRUE; 834674ae819SStefano Zampini } 835674ae819SStefano Zampini } else { 836674ae819SStefano Zampini computeis = PETSC_TRUE; 837674ae819SStefano Zampini computetopography = PETSC_TRUE; 838674ae819SStefano Zampini computesolvers = PETSC_TRUE; 839674ae819SStefano Zampini } 840*f4ddd8eeSStefano Zampini 841*f4ddd8eeSStefano Zampini /* Get stdout for dbg */ 842*f4ddd8eeSStefano Zampini if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) { 843*f4ddd8eeSStefano Zampini ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr); 844*f4ddd8eeSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 845*f4ddd8eeSStefano Zampini if (pcbddc->current_level) { 846*f4ddd8eeSStefano Zampini ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 847*f4ddd8eeSStefano Zampini } 848*f4ddd8eeSStefano Zampini } 849*f4ddd8eeSStefano Zampini 850fcd409f5SStefano Zampini /* Set up all the "iterative substructuring" common block without computing solvers */ 851674ae819SStefano Zampini if (computeis) { 852fcd409f5SStefano Zampini /* HACK INTO PCIS */ 853fcd409f5SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 854fcd409f5SStefano Zampini pcis->computesolvers = PETSC_FALSE; 855674ae819SStefano Zampini ierr = PCISSetUp(pc);CHKERRQ(ierr); 856674ae819SStefano Zampini } 857*f4ddd8eeSStefano Zampini 858*f4ddd8eeSStefano Zampini /* Analyze interface */ 859674ae819SStefano Zampini if (computetopography) { 860674ae819SStefano Zampini ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 861fb8d54d4SStefano Zampini } 862fb8d54d4SStefano Zampini 863*f4ddd8eeSStefano Zampini /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 864fb8d54d4SStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 865*f4ddd8eeSStefano Zampini ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr); 866*f4ddd8eeSStefano Zampini if (pcbddc->onearnullspace) { /* already used nearnullspace */ 867*f4ddd8eeSStefano Zampini if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 868*f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 869*f4ddd8eeSStefano Zampini } else { 870*f4ddd8eeSStefano Zampini /* determine if the two nullspaces are different (should be lightweight) */ 871*f4ddd8eeSStefano Zampini if (nearnullspace != pcbddc->onearnullspace) { 872*f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 873*f4ddd8eeSStefano Zampini } else { /* maybe the user has changed the content of the nearnullspace */ 874*f4ddd8eeSStefano Zampini ierr = MatNullSpaceGetVecs(nearnullspace,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 875*f4ddd8eeSStefano Zampini ierr = MatNullSpaceGetVecs(pcbddc->onearnullspace,&onnsp_has_cnst,&onnsp_size,&onearnullvecs);CHKERRQ(ierr); 876*f4ddd8eeSStefano Zampini if ( (nnsp_has_cnst != onnsp_has_cnst) || (nnsp_size != onnsp_size) ) { 877*f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 878*f4ddd8eeSStefano Zampini } else { /* nullspaces have the same size, so check vectors or their ObjectStateId */ 879*f4ddd8eeSStefano Zampini PetscInt i; 880*f4ddd8eeSStefano Zampini for (i=0;i<nnsp_size;i++) { 881*f4ddd8eeSStefano Zampini ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr); 882*f4ddd8eeSStefano Zampini if (nearnullvecs[i] != onearnullvecs[i] || pcbddc->onearnullvecs_state[i] != state) { 883*f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 884*f4ddd8eeSStefano Zampini break; 885*f4ddd8eeSStefano Zampini } 886*f4ddd8eeSStefano Zampini } 887*f4ddd8eeSStefano Zampini } 888*f4ddd8eeSStefano Zampini } 889*f4ddd8eeSStefano Zampini } 890*f4ddd8eeSStefano Zampini } else { 891*f4ddd8eeSStefano Zampini if (!nearnullspace) { /* both nearnullspaces are null */ 892*f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 893*f4ddd8eeSStefano Zampini } else { /* nearnullspace attached later */ 894*f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 895*f4ddd8eeSStefano Zampini } 896*f4ddd8eeSStefano Zampini } 897*f4ddd8eeSStefano Zampini 898*f4ddd8eeSStefano Zampini /* Setup constraints and related work vectors */ 899*f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 900fb8d54d4SStefano Zampini if (computetopography || new_nearnullspace_provided) { 901674ae819SStefano Zampini ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 902e7b262bdSStefano Zampini /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 903*f4ddd8eeSStefano Zampini ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr); 904*f4ddd8eeSStefano Zampini /* set flag for primal space */ 905*f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_TRUE; 906674ae819SStefano Zampini } 907fb8d54d4SStefano Zampini 908*f4ddd8eeSStefano Zampini if (computesolvers || pcbddc->new_primal_space) { 909674ae819SStefano Zampini /* reset data */ 910674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 91199cc7994SStefano Zampini /* Create coarse and local stuffs */ 91299cc7994SStefano Zampini ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 913674ae819SStefano Zampini ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 9140c7d97c5SJed Brown } 9152b510759SStefano Zampini if (pcbddc->dbg_flag && pcbddc->current_level) { 9162b510759SStefano Zampini ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 9172b510759SStefano Zampini } 9180c7d97c5SJed Brown PetscFunctionReturn(0); 9190c7d97c5SJed Brown } 9200c7d97c5SJed Brown 9210c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 9220c7d97c5SJed Brown /* 9230c7d97c5SJed Brown PCApply_BDDC - Applies the BDDC preconditioner to a vector. 9240c7d97c5SJed Brown 9250c7d97c5SJed Brown Input Parameters: 9260c7d97c5SJed Brown . pc - the preconditioner context 9270c7d97c5SJed Brown . r - input vector (global) 9280c7d97c5SJed Brown 9290c7d97c5SJed Brown Output Parameter: 9300c7d97c5SJed Brown . z - output vector (global) 9310c7d97c5SJed Brown 9320c7d97c5SJed Brown Application Interface Routine: PCApply() 9330c7d97c5SJed Brown */ 9340c7d97c5SJed Brown #undef __FUNCT__ 9350c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC" 93653cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 9370c7d97c5SJed Brown { 9380c7d97c5SJed Brown PC_IS *pcis = (PC_IS*)(pc->data); 9390c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 9400c7d97c5SJed Brown PetscErrorCode ierr; 9413b03a366Sstefano_zampini const PetscScalar one = 1.0; 9423b03a366Sstefano_zampini const PetscScalar m_one = -1.0; 9432617d88aSStefano Zampini const PetscScalar zero = 0.0; 9440c7d97c5SJed Brown 9450c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN 9460c7d97c5SJed Brown NN interface preconditioner changed to BDDC 9478eeda7d8SStefano Zampini Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 9480c7d97c5SJed Brown 9490c7d97c5SJed Brown PetscFunctionBegin; 95085c4d303SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 9510c7d97c5SJed Brown /* First Dirichlet solve */ 9520c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9530c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 95453cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 9550c7d97c5SJed Brown /* 9560c7d97c5SJed Brown Assembling right hand side for BDDC operator 957674ae819SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 958674ae819SStefano Zampini - pcis->vec1_B the interface part of the global vector z 9590c7d97c5SJed Brown */ 9600c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 9610c7d97c5SJed Brown ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 9628eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 9630c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 9640c7d97c5SJed Brown ierr = VecCopy(r,z);CHKERRQ(ierr); 9650c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9660c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 967674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 968b76ba322SStefano Zampini } else { 9690bdf917eSStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 970b76ba322SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 971674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 972b76ba322SStefano Zampini } 973b76ba322SStefano Zampini 9742617d88aSStefano Zampini /* Apply interface preconditioner 9752617d88aSStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 9762617d88aSStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr); 9772617d88aSStefano Zampini 978674ae819SStefano Zampini /* Apply transpose of partition of unity operator */ 979674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 9800c7d97c5SJed Brown 9813b03a366Sstefano_zampini /* Second Dirichlet solve and assembling of output */ 9820c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9830c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9840c7d97c5SJed Brown ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 9858eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 986df187020SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 987df187020SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 988df187020SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY (pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 989df187020SStefano Zampini ierr = VecAXPY (pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 9900c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9910c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9920c7d97c5SJed Brown PetscFunctionReturn(0); 9930c7d97c5SJed Brown } 994da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 995674ae819SStefano Zampini 996da1bb401SStefano Zampini #undef __FUNCT__ 997da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC" 998da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc) 999da1bb401SStefano Zampini { 1000da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1001da1bb401SStefano Zampini PetscErrorCode ierr; 1002da1bb401SStefano Zampini 1003da1bb401SStefano Zampini PetscFunctionBegin; 1004da1bb401SStefano Zampini /* free data created by PCIS */ 1005da1bb401SStefano Zampini ierr = PCISDestroy(pc);CHKERRQ(ierr); 1006674ae819SStefano Zampini /* free BDDC custom data */ 1007674ae819SStefano Zampini ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 1008674ae819SStefano Zampini /* destroy objects related to topography */ 1009674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 1010674ae819SStefano Zampini /* free allocated graph structure */ 1011da1bb401SStefano Zampini ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 1012674ae819SStefano Zampini /* free data for scaling operator */ 1013674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 1014674ae819SStefano Zampini /* free solvers stuff */ 1015674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 101662a6ff1dSStefano Zampini /* free global vectors needed in presolve */ 101762a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 101862a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 10193425bc38SStefano Zampini /* remove functions */ 1020674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 1021bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 10222b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 1023b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 10242b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 1025bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 1026bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 1027bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1028bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 1029bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1030bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 1031bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 1032bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 1033bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 1034bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 1035674ae819SStefano Zampini /* Free the private data structure */ 1036674ae819SStefano Zampini ierr = PetscFree(pc->data);CHKERRQ(ierr); 1037da1bb401SStefano Zampini PetscFunctionReturn(0); 1038da1bb401SStefano Zampini } 10393425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 10401e6b0712SBarry Smith 10413425bc38SStefano Zampini #undef __FUNCT__ 10423425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 10433425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 10443425bc38SStefano Zampini { 1045674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10463425bc38SStefano Zampini PC_IS* pcis; 10473425bc38SStefano Zampini PC_BDDC* pcbddc; 10483425bc38SStefano Zampini PetscErrorCode ierr; 10490c7d97c5SJed Brown 10503425bc38SStefano Zampini PetscFunctionBegin; 10513425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10523425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 10533425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 10543425bc38SStefano Zampini 10553425bc38SStefano Zampini /* change of basis for physical rhs if needed 10563425bc38SStefano Zampini It also changes the rhs in case of dirichlet boundaries */ 10573308cffdSStefano Zampini ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 10583425bc38SStefano Zampini /* store vectors for computation of fetidp final solution */ 10593425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10603425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1061fb223d50SStefano Zampini /* scale rhs since it should be unassembled */ 1062fb223d50SStefano Zampini /* TODO use counter scaling? (also below) */ 10633425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10643425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1065674ae819SStefano Zampini /* Apply partition of unity */ 10663425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1067674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 10688eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 10693425bc38SStefano Zampini /* compute partially subassembled Schur complement right-hand side */ 10703425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10713425bc38SStefano Zampini ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 10723425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 10733425bc38SStefano Zampini ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 10743425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10753425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1076674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 10773425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10783425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10793425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 10803425bc38SStefano Zampini } 10813425bc38SStefano Zampini /* BDDC rhs */ 10823425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 10838eeda7d8SStefano Zampini if (pcbddc->switch_static) { 10843425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10853425bc38SStefano Zampini } 10863425bc38SStefano Zampini /* apply BDDC */ 10873425bc38SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 10883425bc38SStefano Zampini /* Application of B_delta and assembling of rhs for fetidp fluxes */ 10893425bc38SStefano Zampini ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 10903425bc38SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 10913425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10923425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10933425bc38SStefano Zampini /* restore original rhs */ 10943425bc38SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 10953425bc38SStefano Zampini PetscFunctionReturn(0); 10963425bc38SStefano Zampini } 10971e6b0712SBarry Smith 10983425bc38SStefano Zampini #undef __FUNCT__ 10993425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 11003425bc38SStefano Zampini /*@ 110128509bceSStefano Zampini PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system 11023425bc38SStefano Zampini 11033425bc38SStefano Zampini Collective 11043425bc38SStefano Zampini 11053425bc38SStefano Zampini Input Parameters: 110628509bceSStefano Zampini + fetidp_mat - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators 110728509bceSStefano Zampini . standard_rhs - the right-hand side for your linear system 11083425bc38SStefano Zampini 11093425bc38SStefano Zampini Output Parameters: 111028509bceSStefano Zampini - fetidp_flux_rhs - the right-hand side for the FETIDP linear system 11113425bc38SStefano Zampini 11123425bc38SStefano Zampini Level: developer 11133425bc38SStefano Zampini 11143425bc38SStefano Zampini Notes: 11153425bc38SStefano Zampini 111628509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 11173425bc38SStefano Zampini @*/ 11183425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 11193425bc38SStefano Zampini { 1120674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 11213425bc38SStefano Zampini PetscErrorCode ierr; 11223425bc38SStefano Zampini 11233425bc38SStefano Zampini PetscFunctionBegin; 11243425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 11253425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 11263425bc38SStefano Zampini PetscFunctionReturn(0); 11273425bc38SStefano Zampini } 11283425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 11291e6b0712SBarry Smith 11303425bc38SStefano Zampini #undef __FUNCT__ 11313425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 11323425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 11333425bc38SStefano Zampini { 1134674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 11353425bc38SStefano Zampini PC_IS* pcis; 11363425bc38SStefano Zampini PC_BDDC* pcbddc; 11373425bc38SStefano Zampini PetscErrorCode ierr; 11383425bc38SStefano Zampini 11393425bc38SStefano Zampini PetscFunctionBegin; 11403425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 11413425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 11423425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 11433425bc38SStefano Zampini 11443425bc38SStefano Zampini /* apply B_delta^T */ 11453425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11463425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11473425bc38SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 11483425bc38SStefano Zampini /* compute rhs for BDDC application */ 11493425bc38SStefano Zampini ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 11508eeda7d8SStefano Zampini if (pcbddc->switch_static) { 11513425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 11523425bc38SStefano Zampini } 11533425bc38SStefano Zampini /* apply BDDC */ 11543425bc38SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 11553425bc38SStefano Zampini /* put values into standard global vector */ 11563425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11573425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11588eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 11593425bc38SStefano Zampini /* compute values into the interior if solved for the partially subassembled Schur complement */ 11603425bc38SStefano Zampini ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 11613425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 11623425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 11633425bc38SStefano Zampini } 11643425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11653425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11663425bc38SStefano Zampini /* final change of basis if needed 11673425bc38SStefano Zampini Is also sums the dirichlet part removed during RHS assembling */ 11683308cffdSStefano Zampini ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 11693425bc38SStefano Zampini PetscFunctionReturn(0); 11703425bc38SStefano Zampini } 11711e6b0712SBarry Smith 11723425bc38SStefano Zampini #undef __FUNCT__ 11733425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 11743425bc38SStefano Zampini /*@ 117528509bceSStefano Zampini PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system 11763425bc38SStefano Zampini 11773425bc38SStefano Zampini Collective 11783425bc38SStefano Zampini 11793425bc38SStefano Zampini Input Parameters: 118028509bceSStefano Zampini + fetidp_mat - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators 118128509bceSStefano Zampini . fetidp_flux_sol - the solution of the FETIDP linear system 11823425bc38SStefano Zampini 11833425bc38SStefano Zampini Output Parameters: 118428509bceSStefano Zampini - standard_sol - the solution defined on the physical domain 11853425bc38SStefano Zampini 11863425bc38SStefano Zampini Level: developer 11873425bc38SStefano Zampini 11883425bc38SStefano Zampini Notes: 11893425bc38SStefano Zampini 119028509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 11913425bc38SStefano Zampini @*/ 11923425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 11933425bc38SStefano Zampini { 1194674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 11953425bc38SStefano Zampini PetscErrorCode ierr; 11963425bc38SStefano Zampini 11973425bc38SStefano Zampini PetscFunctionBegin; 11983425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 11993425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 12003425bc38SStefano Zampini PetscFunctionReturn(0); 12013425bc38SStefano Zampini } 12023425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 12031e6b0712SBarry Smith 1204f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1205f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1206f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1207f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1208674ae819SStefano Zampini 12093425bc38SStefano Zampini #undef __FUNCT__ 12103425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 12113425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 12123425bc38SStefano Zampini { 1213674ae819SStefano Zampini 1214674ae819SStefano Zampini FETIDPMat_ctx fetidpmat_ctx; 12153425bc38SStefano Zampini Mat newmat; 1216674ae819SStefano Zampini FETIDPPC_ctx fetidppc_ctx; 12173425bc38SStefano Zampini PC newpc; 1218ce94432eSBarry Smith MPI_Comm comm; 12193425bc38SStefano Zampini PetscErrorCode ierr; 12203425bc38SStefano Zampini 12213425bc38SStefano Zampini PetscFunctionBegin; 1222ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 12233425bc38SStefano Zampini /* FETIDP linear matrix */ 12243425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 12253425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 12263425bc38SStefano Zampini ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 12273425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 12283425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 12293425bc38SStefano Zampini ierr = MatSetUp(newmat);CHKERRQ(ierr); 12303425bc38SStefano Zampini /* FETIDP preconditioner */ 12313425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 12323425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 12333425bc38SStefano Zampini ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 12343425bc38SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 12353425bc38SStefano Zampini ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 12363425bc38SStefano Zampini ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 12373425bc38SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 12383425bc38SStefano Zampini ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr); 12393425bc38SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 12403425bc38SStefano Zampini /* return pointers for objects created */ 12413425bc38SStefano Zampini *fetidp_mat=newmat; 12423425bc38SStefano Zampini *fetidp_pc=newpc; 12433425bc38SStefano Zampini PetscFunctionReturn(0); 12443425bc38SStefano Zampini } 12451e6b0712SBarry Smith 12463425bc38SStefano Zampini #undef __FUNCT__ 12473425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 12483425bc38SStefano Zampini /*@ 124928509bceSStefano Zampini PCBDDCCreateFETIDPOperators - Create operators for FETIDP 12503425bc38SStefano Zampini 12513425bc38SStefano Zampini Collective 12523425bc38SStefano Zampini 12533425bc38SStefano Zampini Input Parameters: 125428509bceSStefano Zampini + pc - the BDDC preconditioning context already setup 125528509bceSStefano Zampini 125628509bceSStefano Zampini Output Parameters: 125728509bceSStefano Zampini . fetidp_mat - shell FETIDP matrix object 125828509bceSStefano Zampini . fetidp_pc - shell Dirichlet preconditioner for FETIDP matrix 125928509bceSStefano Zampini 126028509bceSStefano Zampini Options Database Keys: 126128509bceSStefano Zampini - -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers 12623425bc38SStefano Zampini 12633425bc38SStefano Zampini Level: developer 12643425bc38SStefano Zampini 12653425bc38SStefano Zampini Notes: 126628509bceSStefano Zampini Currently the only operation provided for FETIDP matrix is MatMult 12673425bc38SStefano Zampini 12683425bc38SStefano Zampini .seealso: PCBDDC 12693425bc38SStefano Zampini @*/ 12703425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 12713425bc38SStefano Zampini { 12723425bc38SStefano Zampini PetscErrorCode ierr; 12733425bc38SStefano Zampini 12743425bc38SStefano Zampini PetscFunctionBegin; 12753425bc38SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12763425bc38SStefano Zampini if (pc->setupcalled) { 1277516d51a7SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1278f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 12793425bc38SStefano Zampini PetscFunctionReturn(0); 12803425bc38SStefano Zampini } 12810c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 1282da1bb401SStefano Zampini /*MC 1283da1bb401SStefano Zampini PCBDDC - Balancing Domain Decomposition by Constraints. 12840c7d97c5SJed Brown 128528509bceSStefano Zampini An implementation of the BDDC preconditioner based on 128628509bceSStefano Zampini 128728509bceSStefano Zampini .vb 128828509bceSStefano Zampini [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 128928509bceSStefano 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 129028509bceSStefano Zampini [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 129128509bceSStefano Zampini .ve 129228509bceSStefano Zampini 129328509bceSStefano Zampini The matrix to be preconditioned (Pmat) must be of type MATIS. 129428509bceSStefano Zampini 1295b6fdb6dfSStefano Zampini Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 129628509bceSStefano Zampini 129728509bceSStefano Zampini It also works with unsymmetric and indefinite problems. 129828509bceSStefano Zampini 1299b6fdb6dfSStefano 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. 1300b6fdb6dfSStefano Zampini 130128509bceSStefano Zampini Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace 130228509bceSStefano Zampini 130328509bceSStefano 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 130428509bceSStefano Zampini 130528509bceSStefano Zampini Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace. 130628509bceSStefano Zampini 1307b6fdb6dfSStefano 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. 130828509bceSStefano Zampini 130928509bceSStefano Zampini The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object. 131028509bceSStefano Zampini 1311da1bb401SStefano Zampini Options Database Keys: 131228509bceSStefano Zampini 131328509bceSStefano Zampini . -pc_bddc_use_vertices <1> - use or not vertices in primal space 131428509bceSStefano Zampini . -pc_bddc_use_edges <1> - use or not edges in primal space 1315b6fdb6dfSStefano Zampini . -pc_bddc_use_faces <0> - use or not faces in primal space 131628509bceSStefano Zampini . -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only) 131728509bceSStefano Zampini . -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested 131828509bceSStefano Zampini . -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1]) 131928509bceSStefano Zampini . -pc_bddc_levels <0> - maximum number of levels for multilevel 132028509bceSStefano Zampini . -pc_bddc_coarsening_ratio - H/h ratio at the coarser level 132128509bceSStefano Zampini - -pc_bddc_check_level <0> - set verbosity level of debugging output 132228509bceSStefano Zampini 132328509bceSStefano Zampini Options for Dirichlet, Neumann or coarse solver can be set with 132428509bceSStefano Zampini .vb 132528509bceSStefano Zampini -pc_bddc_dirichlet_ 132628509bceSStefano Zampini -pc_bddc_neumann_ 132728509bceSStefano Zampini -pc_bddc_coarse_ 132828509bceSStefano Zampini .ve 132928509bceSStefano Zampini e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg 133028509bceSStefano Zampini 133128509bceSStefano Zampini When using a multilevel approach, solvers' options at the N-th level can be specified as 133228509bceSStefano Zampini .vb 133328509bceSStefano Zampini -pc_bddc_dirichlet_N_ 133428509bceSStefano Zampini -pc_bddc_neumann_N_ 133528509bceSStefano Zampini -pc_bddc_coarse_N_ 133628509bceSStefano Zampini .ve 133728509bceSStefano Zampini Note that level number ranges from the finest 0 to the coarsest N 1338da1bb401SStefano Zampini 1339da1bb401SStefano Zampini Level: intermediate 1340da1bb401SStefano Zampini 1341b6fdb6dfSStefano Zampini Developer notes: 134228509bceSStefano Zampini Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose 1343da1bb401SStefano Zampini 134428509bceSStefano Zampini New deluxe scaling operator will be available soon. 1345da1bb401SStefano Zampini 1346da1bb401SStefano Zampini Contributed by Stefano Zampini 1347da1bb401SStefano Zampini 1348da1bb401SStefano Zampini .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1349da1bb401SStefano Zampini M*/ 1350b2573a8aSBarry Smith 1351da1bb401SStefano Zampini #undef __FUNCT__ 1352da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC" 13538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1354da1bb401SStefano Zampini { 1355da1bb401SStefano Zampini PetscErrorCode ierr; 1356da1bb401SStefano Zampini PC_BDDC *pcbddc; 1357da1bb401SStefano Zampini 1358da1bb401SStefano Zampini PetscFunctionBegin; 1359da1bb401SStefano Zampini /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1360da1bb401SStefano Zampini ierr = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr); 1361da1bb401SStefano Zampini pc->data = (void*)pcbddc; 1362da1bb401SStefano Zampini 1363da1bb401SStefano Zampini /* create PCIS data structure */ 1364da1bb401SStefano Zampini ierr = PCISCreate(pc);CHKERRQ(ierr); 1365da1bb401SStefano Zampini 136647d04d0dSStefano Zampini /* BDDC customization */ 136747d04d0dSStefano Zampini pcbddc->use_vertices = PETSC_TRUE; 136847d04d0dSStefano Zampini pcbddc->use_edges = PETSC_TRUE; 136947d04d0dSStefano Zampini pcbddc->use_faces = PETSC_FALSE; 137047d04d0dSStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 137147d04d0dSStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 137247d04d0dSStefano Zampini pcbddc->switch_static = PETSC_FALSE; 137347d04d0dSStefano Zampini pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 137447d04d0dSStefano Zampini pcbddc->dbg_flag = 0; 137547d04d0dSStefano Zampini 1376*f4ddd8eeSStefano Zampini pcbddc->coarse_size = 0; 1377*f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1378*f4ddd8eeSStefano Zampini pcbddc->global_primal_indices = 0; 1379*f4ddd8eeSStefano Zampini pcbddc->onearnullspace = 0; 1380*f4ddd8eeSStefano Zampini pcbddc->onearnullvecs_state = 0; 1381674ae819SStefano Zampini pcbddc->user_primal_vertices = 0; 13820bdf917eSStefano Zampini pcbddc->NullSpace = 0; 13833972b0daSStefano Zampini pcbddc->temp_solution = 0; 1384534831adSStefano Zampini pcbddc->original_rhs = 0; 1385534831adSStefano Zampini pcbddc->local_mat = 0; 1386534831adSStefano Zampini pcbddc->ChangeOfBasisMatrix = 0; 1387da1bb401SStefano Zampini pcbddc->coarse_vec = 0; 1388da1bb401SStefano Zampini pcbddc->coarse_rhs = 0; 1389da1bb401SStefano Zampini pcbddc->coarse_ksp = 0; 1390da1bb401SStefano Zampini pcbddc->coarse_phi_B = 0; 1391da1bb401SStefano Zampini pcbddc->coarse_phi_D = 0; 139215aaf578SStefano Zampini pcbddc->coarse_psi_B = 0; 139315aaf578SStefano Zampini pcbddc->coarse_psi_D = 0; 1394da1bb401SStefano Zampini pcbddc->vec1_P = 0; 1395da1bb401SStefano Zampini pcbddc->vec1_R = 0; 1396da1bb401SStefano Zampini pcbddc->vec2_R = 0; 1397da1bb401SStefano Zampini pcbddc->local_auxmat1 = 0; 1398da1bb401SStefano Zampini pcbddc->local_auxmat2 = 0; 1399da1bb401SStefano Zampini pcbddc->R_to_B = 0; 1400da1bb401SStefano Zampini pcbddc->R_to_D = 0; 1401da1bb401SStefano Zampini pcbddc->ksp_D = 0; 1402da1bb401SStefano Zampini pcbddc->ksp_R = 0; 1403da1bb401SStefano Zampini pcbddc->NeumannBoundaries = 0; 140460d44989SStefano Zampini pcbddc->user_provided_isfordofs = PETSC_FALSE; 140560d44989SStefano Zampini pcbddc->n_ISForDofs = 0; 1406da1bb401SStefano Zampini pcbddc->ISForDofs = 0; 1407da1bb401SStefano Zampini pcbddc->ConstraintMatrix = 0; 140885c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 140947d04d0dSStefano Zampini pcbddc->coarse_loc_to_glob = 0; 141047d04d0dSStefano Zampini pcbddc->coarsening_ratio = 8; 14114fad6a16SStefano Zampini pcbddc->current_level = 0; 14122b510759SStefano Zampini pcbddc->max_levels = 0; 1413da1bb401SStefano Zampini 1414674ae819SStefano Zampini /* create local graph structure */ 1415674ae819SStefano Zampini ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1416674ae819SStefano Zampini 1417674ae819SStefano Zampini /* scaling */ 1418674ae819SStefano Zampini pcbddc->use_deluxe_scaling = PETSC_FALSE; 1419674ae819SStefano Zampini pcbddc->work_scaling = 0; 1420da1bb401SStefano Zampini 1421da1bb401SStefano Zampini /* function pointers */ 1422da1bb401SStefano Zampini pc->ops->apply = PCApply_BDDC; 1423da1bb401SStefano Zampini pc->ops->applytranspose = 0; 1424da1bb401SStefano Zampini pc->ops->setup = PCSetUp_BDDC; 1425da1bb401SStefano Zampini pc->ops->destroy = PCDestroy_BDDC; 1426da1bb401SStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1427da1bb401SStefano Zampini pc->ops->view = 0; 1428da1bb401SStefano Zampini pc->ops->applyrichardson = 0; 1429da1bb401SStefano Zampini pc->ops->applysymmetricleft = 0; 1430da1bb401SStefano Zampini pc->ops->applysymmetricright = 0; 1431534831adSStefano Zampini pc->ops->presolve = PCPreSolve_BDDC; 1432534831adSStefano Zampini pc->ops->postsolve = PCPostSolve_BDDC; 1433da1bb401SStefano Zampini 1434da1bb401SStefano Zampini /* composing function */ 1435674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1436bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 14372b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1438b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 14392b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1440bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1441bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1442bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1443bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1444bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1445bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 1446bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1447bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1448bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1449bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1450da1bb401SStefano Zampini PetscFunctionReturn(0); 1451da1bb401SStefano Zampini } 14523425bc38SStefano Zampini 1453