153cdbc3dSStefano Zampini /* TODOLIST 2eb97c9d2SStefano Zampini 3eb97c9d2SStefano Zampini ConstraintsSetup 4eb97c9d2SStefano Zampini - tolerances for constraints as an option (take care of single precision!) 5eb97c9d2SStefano Zampini 6eb97c9d2SStefano Zampini Solvers 7a0d3c3abSStefano Zampini - Add support for cholesky for coarse solver (similar to local solvers) 8eb97c9d2SStefano Zampini - Propagate ksp prefixes for solvers to mat objects? 9eb97c9d2SStefano Zampini - Propagate nearnullspace info among levels 10eb97c9d2SStefano Zampini 11eb97c9d2SStefano Zampini User interface 12b9b85e73SStefano Zampini - ** DofSplitting and DM attached to pc? 13eb97c9d2SStefano Zampini 14eb97c9d2SStefano Zampini Debugging output 15b9b85e73SStefano Zampini - * Better management of verbosity levels of debugging output 16eb97c9d2SStefano Zampini 17eb97c9d2SStefano Zampini Build 18eb97c9d2SStefano Zampini - make runexe59 19eb97c9d2SStefano Zampini 20eb97c9d2SStefano Zampini Extra 21b9b85e73SStefano Zampini - ** GetRid of PCBDDCApplySchur, use MatSchur instead 22b9b85e73SStefano Zampini - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)? 23eb97c9d2SStefano Zampini - add support for computing h,H and related using coordinates? 24c0b83709SStefano Zampini - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap) 25eb97c9d2SStefano Zampini - Better management in PCIS code 26eb97c9d2SStefano Zampini - BDDC with MG framework? 27eb97c9d2SStefano Zampini 28eb97c9d2SStefano Zampini FETIDP 29eb97c9d2SStefano Zampini - Move FETIDP code to its own classes 30eb97c9d2SStefano Zampini 31eb97c9d2SStefano Zampini MATIS related operations contained in BDDC code 32eb97c9d2SStefano Zampini - Provide general case for subassembling 33b9b85e73SStefano Zampini - *** Preallocation routines in MatISGetMPIAXAIJ 34eb97c9d2SStefano Zampini 3553cdbc3dSStefano Zampini */ 360c7d97c5SJed Brown 3753cdbc3dSStefano Zampini /* ---------------------------------------------------------------------------------------------------------------------------------------------- 380c7d97c5SJed Brown Implementation of BDDC preconditioner based on: 390c7d97c5SJed Brown C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 4053cdbc3dSStefano Zampini ---------------------------------------------------------------------------------------------------------------------------------------------- */ 4153cdbc3dSStefano Zampini 42ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 43ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 443b03a366Sstefano_zampini #include <petscblaslapack.h> 45674ae819SStefano Zampini 460c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 470c7d97c5SJed Brown #undef __FUNCT__ 480c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC" 490c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc) 500c7d97c5SJed Brown { 510c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 520c7d97c5SJed Brown PetscErrorCode ierr; 530c7d97c5SJed Brown 540c7d97c5SJed Brown PetscFunctionBegin; 550c7d97c5SJed Brown ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 568eeda7d8SStefano Zampini /* Verbose debugging */ 578eeda7d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 588eeda7d8SStefano Zampini /* Primal space cumstomization */ 5908a5cf49SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);CHKERRQ(ierr); 608eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr); 618eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr); 628eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr); 638eeda7d8SStefano Zampini /* Change of basis */ 64b9b85e73SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not internal change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr); 65b9b85e73SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not internal change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr); 66674ae819SStefano Zampini if (!pcbddc->use_change_of_basis) { 67674ae819SStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 68674ae819SStefano Zampini } 698eeda7d8SStefano Zampini /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 708eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);CHKERRQ(ierr); 710298fd71SBarry Smith ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 722b510759SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 73323d395dSStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);CHKERRQ(ierr); 74674ae819SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr); 75*1cef56d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_deluxe_threshold","Deluxe subproblems (Schur principal minors) smaller than this value are explicilty computed (-1 computes all)","none",pcbddc->deluxe_threshold,&pcbddc->deluxe_threshold,NULL);CHKERRQ(ierr); 7634a97f8cSStefano Zampini ierr = PetscOptionsBool("-pc_bddc_deluxe_rebuild","Whether or not the interface graph for deluxe has to be rebuilt (i.e. use a standard definition of the interface)","none",pcbddc->deluxe_rebuild,&pcbddc->deluxe_rebuild,NULL);CHKERRQ(ierr); 7734a97f8cSStefano Zampini ierr = PetscOptionsInt("-pc_bddc_deluxe_layers","Number of dofs' layers for the application of deluxe cheap version (i.e. -1 uses all dofs)","none",pcbddc->deluxe_layers,&pcbddc->deluxe_layers,NULL);CHKERRQ(ierr); 7834a97f8cSStefano Zampini ierr = PetscOptionsBool("-pc_bddc_deluxe_use_useradj","Whether or not the CSR graph specified by the user should be used for computing layers (default is to use adj of local mat)","none",pcbddc->deluxe_use_useradj,&pcbddc->deluxe_use_useradj,NULL);CHKERRQ(ierr); 790c7d97c5SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 800c7d97c5SJed Brown PetscFunctionReturn(0); 810c7d97c5SJed Brown } 820c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 83674ae819SStefano Zampini #undef __FUNCT__ 84b9b85e73SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisLocalMat_BDDC" 85b9b85e73SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisLocalMat_BDDC(PC pc, Mat change) 86b9b85e73SStefano Zampini { 87b9b85e73SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 88b9b85e73SStefano Zampini PetscErrorCode ierr; 89b9b85e73SStefano Zampini 90b9b85e73SStefano Zampini PetscFunctionBegin; 91b9b85e73SStefano Zampini ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr); 92b9b85e73SStefano Zampini ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr); 93b9b85e73SStefano Zampini pcbddc->user_ChangeOfBasisMatrix = change; 94b9b85e73SStefano Zampini PetscFunctionReturn(0); 95b9b85e73SStefano Zampini } 96b9b85e73SStefano Zampini #undef __FUNCT__ 97b9b85e73SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisLocalMat" 98b9b85e73SStefano Zampini /*@ 99b9b85e73SStefano Zampini PCBDDCSetChangeOfBasisLocalMat - Set user defined change of basis for local boundary dofs 100b9b85e73SStefano Zampini 101b9b85e73SStefano Zampini Collective on PC 102b9b85e73SStefano Zampini 103b9b85e73SStefano Zampini Input Parameters: 104b9b85e73SStefano Zampini + pc - the preconditioning context 105b9b85e73SStefano Zampini - change - the local change of basis matrix, either in local (internal + boundary) or in local boundary numbering 106b9b85e73SStefano Zampini 107b9b85e73SStefano Zampini Level: intermediate 108b9b85e73SStefano Zampini 109b9b85e73SStefano Zampini Notes: 110b9b85e73SStefano Zampini 111b9b85e73SStefano Zampini .seealso: PCBDDC 112b9b85e73SStefano Zampini @*/ 113b9b85e73SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisLocalMat(PC pc, Mat change) 114b9b85e73SStefano Zampini { 115b9b85e73SStefano Zampini PetscErrorCode ierr; 116b9b85e73SStefano Zampini 117b9b85e73SStefano Zampini PetscFunctionBegin; 118b9b85e73SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 119b9b85e73SStefano Zampini PetscValidHeaderSpecific(change,MAT_CLASSID,2); 120b9b85e73SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisLocalMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr); 121b9b85e73SStefano Zampini PetscFunctionReturn(0); 122b9b85e73SStefano Zampini } 123b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */ 124b9b85e73SStefano Zampini #undef __FUNCT__ 125674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 126674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 127674ae819SStefano Zampini { 128674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 129674ae819SStefano Zampini PetscErrorCode ierr; 1301e6b0712SBarry Smith 131674ae819SStefano Zampini PetscFunctionBegin; 132674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 133674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 134674ae819SStefano Zampini pcbddc->user_primal_vertices = PrimalVertices; 135674ae819SStefano Zampini PetscFunctionReturn(0); 136674ae819SStefano Zampini } 137674ae819SStefano Zampini #undef __FUNCT__ 138674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 139674ae819SStefano Zampini /*@ 14028509bceSStefano Zampini PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 141674ae819SStefano Zampini 142674ae819SStefano Zampini Not collective 143674ae819SStefano Zampini 144674ae819SStefano Zampini Input Parameters: 145674ae819SStefano Zampini + pc - the preconditioning context 14628509bceSStefano Zampini - PrimalVertices - index set of primal vertices in local numbering 147674ae819SStefano Zampini 148674ae819SStefano Zampini Level: intermediate 149674ae819SStefano Zampini 150674ae819SStefano Zampini Notes: 151674ae819SStefano Zampini 152674ae819SStefano Zampini .seealso: PCBDDC 153674ae819SStefano Zampini @*/ 154674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 155674ae819SStefano Zampini { 156674ae819SStefano Zampini PetscErrorCode ierr; 157674ae819SStefano Zampini 158674ae819SStefano Zampini PetscFunctionBegin; 159674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 160674ae819SStefano Zampini PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 161674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 162674ae819SStefano Zampini PetscFunctionReturn(0); 163674ae819SStefano Zampini } 164674ae819SStefano Zampini /* -------------------------------------------------------------------------- */ 1650c7d97c5SJed Brown #undef __FUNCT__ 1664fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 1674fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 1684fad6a16SStefano Zampini { 1694fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1704fad6a16SStefano Zampini 1714fad6a16SStefano Zampini PetscFunctionBegin; 1724fad6a16SStefano Zampini pcbddc->coarsening_ratio = k; 1734fad6a16SStefano Zampini PetscFunctionReturn(0); 1744fad6a16SStefano Zampini } 1751e6b0712SBarry Smith 1764fad6a16SStefano Zampini #undef __FUNCT__ 1774fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio" 1784fad6a16SStefano Zampini /*@ 17928509bceSStefano Zampini PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 1804fad6a16SStefano Zampini 1814fad6a16SStefano Zampini Logically collective on PC 1824fad6a16SStefano Zampini 1834fad6a16SStefano Zampini Input Parameters: 1844fad6a16SStefano Zampini + pc - the preconditioning context 18528509bceSStefano Zampini - k - coarsening ratio (H/h at the coarser level) 1864fad6a16SStefano Zampini 18728509bceSStefano Zampini Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 1884fad6a16SStefano Zampini 1894fad6a16SStefano Zampini Level: intermediate 1904fad6a16SStefano Zampini 1914fad6a16SStefano Zampini Notes: 1924fad6a16SStefano Zampini 1934fad6a16SStefano Zampini .seealso: PCBDDC 1944fad6a16SStefano Zampini @*/ 1954fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 1964fad6a16SStefano Zampini { 1974fad6a16SStefano Zampini PetscErrorCode ierr; 1984fad6a16SStefano Zampini 1994fad6a16SStefano Zampini PetscFunctionBegin; 2004fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2012b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,k,2); 2024fad6a16SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 2034fad6a16SStefano Zampini PetscFunctionReturn(0); 2044fad6a16SStefano Zampini } 2052b510759SStefano Zampini 206b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 2072b510759SStefano Zampini #undef __FUNCT__ 208b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 209b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 210b8ffe317SStefano Zampini { 211b8ffe317SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 212b8ffe317SStefano Zampini 213b8ffe317SStefano Zampini PetscFunctionBegin; 21485c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = flg; 215b8ffe317SStefano Zampini PetscFunctionReturn(0); 216b8ffe317SStefano Zampini } 217b8ffe317SStefano Zampini 218b8ffe317SStefano Zampini #undef __FUNCT__ 219b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 220b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 2212b510759SStefano Zampini { 2222b510759SStefano Zampini PetscErrorCode ierr; 2232b510759SStefano Zampini 2242b510759SStefano Zampini PetscFunctionBegin; 2252b510759SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 226b8ffe317SStefano Zampini PetscValidLogicalCollectiveBool(pc,flg,2); 227b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 2282b510759SStefano Zampini PetscFunctionReturn(0); 2292b510759SStefano Zampini } 2301e6b0712SBarry Smith 2314fad6a16SStefano Zampini #undef __FUNCT__ 2322b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC" 2332b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 2344fad6a16SStefano Zampini { 2354fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2364fad6a16SStefano Zampini 2374fad6a16SStefano Zampini PetscFunctionBegin; 2382b510759SStefano Zampini pcbddc->current_level = level; 2394fad6a16SStefano Zampini PetscFunctionReturn(0); 2404fad6a16SStefano Zampini } 2411e6b0712SBarry Smith 2424fad6a16SStefano Zampini #undef __FUNCT__ 243b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel" 244b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 245b8ffe317SStefano Zampini { 246b8ffe317SStefano Zampini PetscErrorCode ierr; 247b8ffe317SStefano Zampini 248b8ffe317SStefano Zampini PetscFunctionBegin; 249b8ffe317SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 250b8ffe317SStefano Zampini PetscValidLogicalCollectiveInt(pc,level,2); 251b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 252b8ffe317SStefano Zampini PetscFunctionReturn(0); 253b8ffe317SStefano Zampini } 254b8ffe317SStefano Zampini 255b8ffe317SStefano Zampini #undef __FUNCT__ 2562b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC" 2572b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 2582b510759SStefano Zampini { 2592b510759SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2602b510759SStefano Zampini 2612b510759SStefano Zampini PetscFunctionBegin; 2622b510759SStefano Zampini pcbddc->max_levels = levels; 2632b510759SStefano Zampini PetscFunctionReturn(0); 2642b510759SStefano Zampini } 2652b510759SStefano Zampini 2662b510759SStefano Zampini #undef __FUNCT__ 2672b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels" 2684fad6a16SStefano Zampini /*@ 26928509bceSStefano Zampini PCBDDCSetLevels - Sets the maximum number of levels for multilevel 2704fad6a16SStefano Zampini 2714fad6a16SStefano Zampini Logically collective on PC 2724fad6a16SStefano Zampini 2734fad6a16SStefano Zampini Input Parameters: 2744fad6a16SStefano Zampini + pc - the preconditioning context 27528509bceSStefano Zampini - levels - the maximum number of levels (max 9) 2764fad6a16SStefano Zampini 27728509bceSStefano Zampini Default value is 0, i.e. traditional one-level BDDC 2784fad6a16SStefano Zampini 2794fad6a16SStefano Zampini Level: intermediate 2804fad6a16SStefano Zampini 2814fad6a16SStefano Zampini Notes: 2824fad6a16SStefano Zampini 2834fad6a16SStefano Zampini .seealso: PCBDDC 2844fad6a16SStefano Zampini @*/ 2852b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 2864fad6a16SStefano Zampini { 2874fad6a16SStefano Zampini PetscErrorCode ierr; 2884fad6a16SStefano Zampini 2894fad6a16SStefano Zampini PetscFunctionBegin; 2904fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2912b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,levels,2); 292312be037SStefano Zampini if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n"); 2932b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 2944fad6a16SStefano Zampini PetscFunctionReturn(0); 2954fad6a16SStefano Zampini } 2964fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */ 2971e6b0712SBarry Smith 2984fad6a16SStefano Zampini #undef __FUNCT__ 2990bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 3000bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 3010bdf917eSStefano Zampini { 3020bdf917eSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3030bdf917eSStefano Zampini PetscErrorCode ierr; 3040bdf917eSStefano Zampini 3050bdf917eSStefano Zampini PetscFunctionBegin; 3060bdf917eSStefano Zampini ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 3070bdf917eSStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 3080bdf917eSStefano Zampini pcbddc->NullSpace = NullSpace; 3090bdf917eSStefano Zampini PetscFunctionReturn(0); 3100bdf917eSStefano Zampini } 3111e6b0712SBarry Smith 3120bdf917eSStefano Zampini #undef __FUNCT__ 3130bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace" 3140bdf917eSStefano Zampini /*@ 31528509bceSStefano Zampini PCBDDCSetNullSpace - Set nullspace for BDDC operator 3160bdf917eSStefano Zampini 3170bdf917eSStefano Zampini Logically collective on PC and MatNullSpace 3180bdf917eSStefano Zampini 3190bdf917eSStefano Zampini Input Parameters: 3200bdf917eSStefano Zampini + pc - the preconditioning context 32128509bceSStefano Zampini - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 3220bdf917eSStefano Zampini 3230bdf917eSStefano Zampini Level: intermediate 3240bdf917eSStefano Zampini 3250bdf917eSStefano Zampini Notes: 3260bdf917eSStefano Zampini 3270bdf917eSStefano Zampini .seealso: PCBDDC 3280bdf917eSStefano Zampini @*/ 3290bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 3300bdf917eSStefano Zampini { 3310bdf917eSStefano Zampini PetscErrorCode ierr; 3320bdf917eSStefano Zampini 3330bdf917eSStefano Zampini PetscFunctionBegin; 3340bdf917eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 335674ae819SStefano Zampini PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 3362b510759SStefano Zampini PetscCheckSameComm(pc,1,NullSpace,2); 3370bdf917eSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 3380bdf917eSStefano Zampini PetscFunctionReturn(0); 3390bdf917eSStefano Zampini } 3400bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */ 3411e6b0712SBarry Smith 3420bdf917eSStefano Zampini #undef __FUNCT__ 3433b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 3443b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 3453b03a366Sstefano_zampini { 3463b03a366Sstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3473b03a366Sstefano_zampini PetscErrorCode ierr; 3483b03a366Sstefano_zampini 3493b03a366Sstefano_zampini PetscFunctionBegin; 350785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 351785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 3523b03a366Sstefano_zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 35336e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 35436e030ebSStefano Zampini pcbddc->DirichletBoundaries = DirichletBoundaries; 355fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 3563b03a366Sstefano_zampini PetscFunctionReturn(0); 3573b03a366Sstefano_zampini } 3581e6b0712SBarry Smith 3593b03a366Sstefano_zampini #undef __FUNCT__ 3603b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 3613b03a366Sstefano_zampini /*@ 36228509bceSStefano Zampini PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 3633b03a366Sstefano_zampini 364785d1243SStefano Zampini Collective 3653b03a366Sstefano_zampini 3663b03a366Sstefano_zampini Input Parameters: 3673b03a366Sstefano_zampini + pc - the preconditioning context 368785d1243SStefano Zampini - DirichletBoundaries - parallel IS defining the Dirichlet boundaries 3693b03a366Sstefano_zampini 3703b03a366Sstefano_zampini Level: intermediate 3713b03a366Sstefano_zampini 372785d1243SStefano Zampini Notes: Any process can list any global node 3733b03a366Sstefano_zampini 3743b03a366Sstefano_zampini .seealso: PCBDDC 3753b03a366Sstefano_zampini @*/ 3763b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 3773b03a366Sstefano_zampini { 3783b03a366Sstefano_zampini PetscErrorCode ierr; 3793b03a366Sstefano_zampini 3803b03a366Sstefano_zampini PetscFunctionBegin; 3813b03a366Sstefano_zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 382674ae819SStefano Zampini PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 383785d1243SStefano Zampini PetscCheckSameComm(pc,1,DirichletBoundaries,2); 3843b03a366Sstefano_zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 3853b03a366Sstefano_zampini PetscFunctionReturn(0); 3863b03a366Sstefano_zampini } 3873b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */ 3881e6b0712SBarry Smith 3893b03a366Sstefano_zampini #undef __FUNCT__ 39082ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC" 39182ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries) 3920c7d97c5SJed Brown { 3930c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3940c7d97c5SJed Brown PetscErrorCode ierr; 3950c7d97c5SJed Brown 3960c7d97c5SJed Brown PetscFunctionBegin; 397785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 398785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 3990c7d97c5SJed Brown ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 4000c7d97c5SJed Brown ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 401785d1243SStefano Zampini pcbddc->DirichletBoundariesLocal = DirichletBoundaries; 4027d24d155SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4030c7d97c5SJed Brown PetscFunctionReturn(0); 4040c7d97c5SJed Brown } 4050c7d97c5SJed Brown 4060c7d97c5SJed Brown #undef __FUNCT__ 40782ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal" 4089c0446d6SStefano Zampini /*@ 40982ba6b80SStefano Zampini PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering. 4109c0446d6SStefano Zampini 411785d1243SStefano Zampini Collective 41253cdbc3dSStefano Zampini 41353cdbc3dSStefano Zampini Input Parameters: 41453cdbc3dSStefano Zampini + pc - the preconditioning context 41582ba6b80SStefano Zampini - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering) 41653cdbc3dSStefano Zampini 41753cdbc3dSStefano Zampini Level: intermediate 41853cdbc3dSStefano Zampini 4199c0446d6SStefano Zampini Notes: 42053cdbc3dSStefano Zampini 42153cdbc3dSStefano Zampini .seealso: PCBDDC 42253cdbc3dSStefano Zampini @*/ 42382ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries) 4240c7d97c5SJed Brown { 4250c7d97c5SJed Brown PetscErrorCode ierr; 4260c7d97c5SJed Brown 4270c7d97c5SJed Brown PetscFunctionBegin; 4280c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4290c7d97c5SJed Brown PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 43082ba6b80SStefano Zampini PetscCheckSameComm(pc,1,DirichletBoundaries,2); 43182ba6b80SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 4320c7d97c5SJed Brown PetscFunctionReturn(0); 4330c7d97c5SJed Brown } 4340c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 4350c7d97c5SJed Brown 4360c7d97c5SJed Brown #undef __FUNCT__ 4370c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 43853cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 4390c7d97c5SJed Brown { 4400c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 44153cdbc3dSStefano Zampini PetscErrorCode ierr; 4420c7d97c5SJed Brown 4430c7d97c5SJed Brown PetscFunctionBegin; 444785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 445785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 44653cdbc3dSStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 44736e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 44836e030ebSStefano Zampini pcbddc->NeumannBoundaries = NeumannBoundaries; 449fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4500c7d97c5SJed Brown PetscFunctionReturn(0); 4510c7d97c5SJed Brown } 4521e6b0712SBarry Smith 4530c7d97c5SJed Brown #undef __FUNCT__ 4540c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 45557527edcSJed Brown /*@ 45628509bceSStefano Zampini PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 45757527edcSJed Brown 458785d1243SStefano Zampini Collective 45957527edcSJed Brown 46057527edcSJed Brown Input Parameters: 46157527edcSJed Brown + pc - the preconditioning context 462785d1243SStefano Zampini - NeumannBoundaries - parallel IS defining the Neumann boundaries 46357527edcSJed Brown 46457527edcSJed Brown Level: intermediate 46557527edcSJed Brown 466785d1243SStefano Zampini Notes: Any process can list any global node 46757527edcSJed Brown 46857527edcSJed Brown .seealso: PCBDDC 46957527edcSJed Brown @*/ 47053cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 4710c7d97c5SJed Brown { 4720c7d97c5SJed Brown PetscErrorCode ierr; 4730c7d97c5SJed Brown 4740c7d97c5SJed Brown PetscFunctionBegin; 4750c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 476674ae819SStefano Zampini PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 477785d1243SStefano Zampini PetscCheckSameComm(pc,1,NeumannBoundaries,2); 47853cdbc3dSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 47953cdbc3dSStefano Zampini PetscFunctionReturn(0); 48053cdbc3dSStefano Zampini } 48153cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 4821e6b0712SBarry Smith 48353cdbc3dSStefano Zampini #undef __FUNCT__ 48482ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC" 48582ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries) 4860c7d97c5SJed Brown { 4870c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 4880c7d97c5SJed Brown PetscErrorCode ierr; 4890c7d97c5SJed Brown 4900c7d97c5SJed Brown PetscFunctionBegin; 491785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 492785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 4930c7d97c5SJed Brown ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 4940c7d97c5SJed Brown ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 495785d1243SStefano Zampini pcbddc->NeumannBoundariesLocal = NeumannBoundaries; 4967d24d155SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4970c7d97c5SJed Brown PetscFunctionReturn(0); 4980c7d97c5SJed Brown } 4990c7d97c5SJed Brown 5000c7d97c5SJed Brown #undef __FUNCT__ 50182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal" 5020c7d97c5SJed Brown /*@ 50382ba6b80SStefano Zampini PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering. 5040c7d97c5SJed Brown 505785d1243SStefano Zampini Collective 5060c7d97c5SJed Brown 5070c7d97c5SJed Brown Input Parameters: 5080c7d97c5SJed Brown + pc - the preconditioning context 50982ba6b80SStefano Zampini - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering) 5100c7d97c5SJed Brown 5110c7d97c5SJed Brown Level: intermediate 5120c7d97c5SJed Brown 5130c7d97c5SJed Brown Notes: 5140c7d97c5SJed Brown 5150c7d97c5SJed Brown .seealso: PCBDDC 5160c7d97c5SJed Brown @*/ 51782ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries) 5180c7d97c5SJed Brown { 5190c7d97c5SJed Brown PetscErrorCode ierr; 5200c7d97c5SJed Brown 5210c7d97c5SJed Brown PetscFunctionBegin; 5220c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5230c7d97c5SJed Brown PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 52482ba6b80SStefano Zampini PetscCheckSameComm(pc,1,NeumannBoundaries,2); 52582ba6b80SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 52653cdbc3dSStefano Zampini PetscFunctionReturn(0); 52753cdbc3dSStefano Zampini } 52853cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 52953cdbc3dSStefano Zampini 53053cdbc3dSStefano Zampini #undef __FUNCT__ 531da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 532da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 533da1bb401SStefano Zampini { 534da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 535da1bb401SStefano Zampini 536da1bb401SStefano Zampini PetscFunctionBegin; 537da1bb401SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundaries; 538da1bb401SStefano Zampini PetscFunctionReturn(0); 539da1bb401SStefano Zampini } 5401e6b0712SBarry Smith 541da1bb401SStefano Zampini #undef __FUNCT__ 542da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 543da1bb401SStefano Zampini /*@ 544785d1243SStefano Zampini PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries 545da1bb401SStefano Zampini 546785d1243SStefano Zampini Collective 547785d1243SStefano Zampini 548785d1243SStefano Zampini Input Parameters: 549785d1243SStefano Zampini . pc - the preconditioning context 550785d1243SStefano Zampini 551785d1243SStefano Zampini Output Parameters: 552785d1243SStefano Zampini . DirichletBoundaries - index set defining the Dirichlet boundaries 553785d1243SStefano Zampini 554785d1243SStefano Zampini Level: intermediate 555785d1243SStefano Zampini 556785d1243SStefano Zampini Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries 557785d1243SStefano Zampini 558785d1243SStefano Zampini .seealso: PCBDDC 559785d1243SStefano Zampini @*/ 560785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 561785d1243SStefano Zampini { 562785d1243SStefano Zampini PetscErrorCode ierr; 563785d1243SStefano Zampini 564785d1243SStefano Zampini PetscFunctionBegin; 565785d1243SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 566785d1243SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 567785d1243SStefano Zampini PetscFunctionReturn(0); 568785d1243SStefano Zampini } 569785d1243SStefano Zampini /* -------------------------------------------------------------------------- */ 570785d1243SStefano Zampini 571785d1243SStefano Zampini #undef __FUNCT__ 572785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC" 573785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries) 574785d1243SStefano Zampini { 575785d1243SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 576785d1243SStefano Zampini 577785d1243SStefano Zampini PetscFunctionBegin; 578785d1243SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundariesLocal; 579785d1243SStefano Zampini PetscFunctionReturn(0); 580785d1243SStefano Zampini } 581785d1243SStefano Zampini 582785d1243SStefano Zampini #undef __FUNCT__ 58382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal" 584da1bb401SStefano Zampini /*@ 58582ba6b80SStefano Zampini PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering) 586da1bb401SStefano Zampini 587785d1243SStefano Zampini Collective 588da1bb401SStefano Zampini 589da1bb401SStefano Zampini Input Parameters: 59028509bceSStefano Zampini . pc - the preconditioning context 591da1bb401SStefano Zampini 592da1bb401SStefano Zampini Output Parameters: 59328509bceSStefano Zampini . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 594da1bb401SStefano Zampini 595da1bb401SStefano Zampini Level: intermediate 596da1bb401SStefano Zampini 597da1bb401SStefano Zampini Notes: 598da1bb401SStefano Zampini 599da1bb401SStefano Zampini .seealso: PCBDDC 600da1bb401SStefano Zampini @*/ 60182ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries) 602da1bb401SStefano Zampini { 603da1bb401SStefano Zampini PetscErrorCode ierr; 604da1bb401SStefano Zampini 605da1bb401SStefano Zampini PetscFunctionBegin; 606da1bb401SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 60782ba6b80SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 608da1bb401SStefano Zampini PetscFunctionReturn(0); 609da1bb401SStefano Zampini } 610da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 6111e6b0712SBarry Smith 612da1bb401SStefano Zampini #undef __FUNCT__ 61353cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 61453cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 61553cdbc3dSStefano Zampini { 61653cdbc3dSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 61753cdbc3dSStefano Zampini 61853cdbc3dSStefano Zampini PetscFunctionBegin; 61953cdbc3dSStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundaries; 62053cdbc3dSStefano Zampini PetscFunctionReturn(0); 62153cdbc3dSStefano Zampini } 6221e6b0712SBarry Smith 62353cdbc3dSStefano Zampini #undef __FUNCT__ 62453cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 62553cdbc3dSStefano Zampini /*@ 626785d1243SStefano Zampini PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries 62753cdbc3dSStefano Zampini 628785d1243SStefano Zampini Collective 629785d1243SStefano Zampini 630785d1243SStefano Zampini Input Parameters: 631785d1243SStefano Zampini . pc - the preconditioning context 632785d1243SStefano Zampini 633785d1243SStefano Zampini Output Parameters: 634785d1243SStefano Zampini . NeumannBoundaries - index set defining the Neumann boundaries 635785d1243SStefano Zampini 636785d1243SStefano Zampini Level: intermediate 637785d1243SStefano Zampini 638785d1243SStefano Zampini Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries 639785d1243SStefano Zampini 640785d1243SStefano Zampini .seealso: PCBDDC 641785d1243SStefano Zampini @*/ 642785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 643785d1243SStefano Zampini { 644785d1243SStefano Zampini PetscErrorCode ierr; 645785d1243SStefano Zampini 646785d1243SStefano Zampini PetscFunctionBegin; 647785d1243SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 648785d1243SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 649785d1243SStefano Zampini PetscFunctionReturn(0); 650785d1243SStefano Zampini } 651785d1243SStefano Zampini /* -------------------------------------------------------------------------- */ 652785d1243SStefano Zampini 653785d1243SStefano Zampini #undef __FUNCT__ 654785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC" 655785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries) 656785d1243SStefano Zampini { 657785d1243SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 658785d1243SStefano Zampini 659785d1243SStefano Zampini PetscFunctionBegin; 660785d1243SStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundariesLocal; 661785d1243SStefano Zampini PetscFunctionReturn(0); 662785d1243SStefano Zampini } 663785d1243SStefano Zampini 664785d1243SStefano Zampini #undef __FUNCT__ 66582ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal" 66653cdbc3dSStefano Zampini /*@ 66782ba6b80SStefano Zampini PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering) 66853cdbc3dSStefano Zampini 669785d1243SStefano Zampini Collective 67053cdbc3dSStefano Zampini 67153cdbc3dSStefano Zampini Input Parameters: 67228509bceSStefano Zampini . pc - the preconditioning context 67353cdbc3dSStefano Zampini 67453cdbc3dSStefano Zampini Output Parameters: 67528509bceSStefano Zampini . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 67653cdbc3dSStefano Zampini 67753cdbc3dSStefano Zampini Level: intermediate 67853cdbc3dSStefano Zampini 67953cdbc3dSStefano Zampini Notes: 68053cdbc3dSStefano Zampini 68153cdbc3dSStefano Zampini .seealso: PCBDDC 68253cdbc3dSStefano Zampini @*/ 68382ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries) 68453cdbc3dSStefano Zampini { 68553cdbc3dSStefano Zampini PetscErrorCode ierr; 68653cdbc3dSStefano Zampini 68753cdbc3dSStefano Zampini PetscFunctionBegin; 68853cdbc3dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 68982ba6b80SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 6900c7d97c5SJed Brown PetscFunctionReturn(0); 6910c7d97c5SJed Brown } 69236e030ebSStefano Zampini /* -------------------------------------------------------------------------- */ 6931e6b0712SBarry Smith 69436e030ebSStefano Zampini #undef __FUNCT__ 695da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 6961a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 69736e030ebSStefano Zampini { 69836e030ebSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 699da1bb401SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 700da1bb401SStefano Zampini PetscErrorCode ierr; 70136e030ebSStefano Zampini 70236e030ebSStefano Zampini PetscFunctionBegin; 703674ae819SStefano Zampini /* free old CSR */ 704674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 705fb223d50SStefano Zampini /* TODO: PCBDDCGraphSetAdjacency */ 706674ae819SStefano Zampini /* get CSR into graph structure */ 707da1bb401SStefano Zampini if (copymode == PETSC_COPY_VALUES) { 708785e854fSJed Brown ierr = PetscMalloc1((nvtxs+1),&mat_graph->xadj);CHKERRQ(ierr); 709785e854fSJed Brown ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr); 710674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 711674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 712da1bb401SStefano Zampini } else if (copymode == PETSC_OWN_POINTER) { 7131a83f524SJed Brown mat_graph->xadj = (PetscInt*)xadj; 7141a83f524SJed Brown mat_graph->adjncy = (PetscInt*)adjncy; 715674ae819SStefano Zampini } 716575ad6abSStefano Zampini mat_graph->nvtxs_csr = nvtxs; 71736e030ebSStefano Zampini PetscFunctionReturn(0); 71836e030ebSStefano Zampini } 7191e6b0712SBarry Smith 72036e030ebSStefano Zampini #undef __FUNCT__ 721da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 72236e030ebSStefano Zampini /*@ 72328509bceSStefano Zampini PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix 72436e030ebSStefano Zampini 72536e030ebSStefano Zampini Not collective 72636e030ebSStefano Zampini 72736e030ebSStefano Zampini Input Parameters: 72836e030ebSStefano Zampini + pc - the preconditioning context 72928509bceSStefano Zampini . nvtxs - number of local vertices of the graph (i.e., the local size of your problem) 73028509bceSStefano Zampini . xadj, adjncy - the CSR graph 73128509bceSStefano Zampini - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 73236e030ebSStefano Zampini 73336e030ebSStefano Zampini Level: intermediate 73436e030ebSStefano Zampini 73536e030ebSStefano Zampini Notes: 73636e030ebSStefano Zampini 73728509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode 73836e030ebSStefano Zampini @*/ 7391a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 74036e030ebSStefano Zampini { 741575ad6abSStefano Zampini void (*f)(void) = 0; 74236e030ebSStefano Zampini PetscErrorCode ierr; 74336e030ebSStefano Zampini 74436e030ebSStefano Zampini PetscFunctionBegin; 74536e030ebSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 746674ae819SStefano Zampini PetscValidIntPointer(xadj,3); 747674ae819SStefano Zampini PetscValidIntPointer(xadj,4); 748674ae819SStefano Zampini if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 749674ae819SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 750da1bb401SStefano Zampini } 75136e030ebSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 752575ad6abSStefano Zampini /* free arrays if PCBDDC is not the PC type */ 753575ad6abSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 754575ad6abSStefano Zampini if (!f && copymode == PETSC_OWN_POINTER) { 755575ad6abSStefano Zampini ierr = PetscFree(xadj);CHKERRQ(ierr); 756575ad6abSStefano Zampini ierr = PetscFree(adjncy);CHKERRQ(ierr); 75736e030ebSStefano Zampini } 75836e030ebSStefano Zampini PetscFunctionReturn(0); 75936e030ebSStefano Zampini } 7609c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */ 7611e6b0712SBarry Smith 7629c0446d6SStefano Zampini #undef __FUNCT__ 76363602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC" 76463602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 76563602bcaSStefano Zampini { 76663602bcaSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 76763602bcaSStefano Zampini PetscInt i; 76863602bcaSStefano Zampini PetscErrorCode ierr; 76963602bcaSStefano Zampini 77063602bcaSStefano Zampini PetscFunctionBegin; 77163602bcaSStefano Zampini /* Destroy ISes if they were already set */ 77263602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 77363602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 77463602bcaSStefano Zampini } 77563602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 77663602bcaSStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 77763602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 77863602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 77963602bcaSStefano Zampini } 78063602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 78163602bcaSStefano Zampini pcbddc->n_ISForDofs = 0; 78263602bcaSStefano Zampini /* allocate space then set */ 78363602bcaSStefano Zampini ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr); 78463602bcaSStefano Zampini for (i=0;i<n_is;i++) { 78563602bcaSStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 78663602bcaSStefano Zampini pcbddc->ISForDofsLocal[i]=ISForDofs[i]; 78763602bcaSStefano Zampini } 78863602bcaSStefano Zampini pcbddc->n_ISForDofsLocal=n_is; 78963602bcaSStefano Zampini if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 7901a8a16d1SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 79163602bcaSStefano Zampini PetscFunctionReturn(0); 79263602bcaSStefano Zampini } 79363602bcaSStefano Zampini 79463602bcaSStefano Zampini #undef __FUNCT__ 79563602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal" 79663602bcaSStefano Zampini /*@ 79763602bcaSStefano Zampini PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix 79863602bcaSStefano Zampini 79963602bcaSStefano Zampini Collective 80063602bcaSStefano Zampini 80163602bcaSStefano Zampini Input Parameters: 80263602bcaSStefano Zampini + pc - the preconditioning context 80363602bcaSStefano Zampini - n_is - number of index sets defining the fields 80463602bcaSStefano Zampini . ISForDofs - array of IS describing the fields in local ordering 80563602bcaSStefano Zampini 80663602bcaSStefano Zampini Level: intermediate 80763602bcaSStefano Zampini 80863602bcaSStefano Zampini Notes: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to a different field. 80963602bcaSStefano Zampini 81063602bcaSStefano Zampini .seealso: PCBDDC 81163602bcaSStefano Zampini @*/ 81263602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[]) 81363602bcaSStefano Zampini { 81463602bcaSStefano Zampini PetscInt i; 81563602bcaSStefano Zampini PetscErrorCode ierr; 81663602bcaSStefano Zampini 81763602bcaSStefano Zampini PetscFunctionBegin; 81863602bcaSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 81963602bcaSStefano Zampini PetscValidLogicalCollectiveInt(pc,n_is,2); 82063602bcaSStefano Zampini for (i=0;i<n_is;i++) { 82163602bcaSStefano Zampini PetscCheckSameComm(pc,1,ISForDofs[i],3); 82263602bcaSStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 82363602bcaSStefano Zampini } 82463602bcaSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 82563602bcaSStefano Zampini PetscFunctionReturn(0); 82663602bcaSStefano Zampini } 82763602bcaSStefano Zampini /* -------------------------------------------------------------------------- */ 82863602bcaSStefano Zampini 82963602bcaSStefano Zampini #undef __FUNCT__ 8309c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 8319c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 8329c0446d6SStefano Zampini { 8339c0446d6SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 8349c0446d6SStefano Zampini PetscInt i; 8359c0446d6SStefano Zampini PetscErrorCode ierr; 8369c0446d6SStefano Zampini 8379c0446d6SStefano Zampini PetscFunctionBegin; 838da1bb401SStefano Zampini /* Destroy ISes if they were already set */ 8399c0446d6SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 8409c0446d6SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 8419c0446d6SStefano Zampini } 842d11ae9bbSstefano_zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 84363602bcaSStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 84463602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 84563602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 84663602bcaSStefano Zampini } 84763602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 84863602bcaSStefano Zampini pcbddc->n_ISForDofsLocal = 0; 849da1bb401SStefano Zampini /* allocate space then set */ 850785e854fSJed Brown ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr); 8519c0446d6SStefano Zampini for (i=0;i<n_is;i++) { 852da1bb401SStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 853da1bb401SStefano Zampini pcbddc->ISForDofs[i]=ISForDofs[i]; 8549c0446d6SStefano Zampini } 8559c0446d6SStefano Zampini pcbddc->n_ISForDofs=n_is; 85663602bcaSStefano Zampini if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 8571a8a16d1SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 8589c0446d6SStefano Zampini PetscFunctionReturn(0); 8599c0446d6SStefano Zampini } 8601e6b0712SBarry Smith 8619c0446d6SStefano Zampini #undef __FUNCT__ 8629c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting" 8639c0446d6SStefano Zampini /*@ 86463602bcaSStefano Zampini PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix 8659c0446d6SStefano Zampini 86663602bcaSStefano Zampini Collective 8679c0446d6SStefano Zampini 8689c0446d6SStefano Zampini Input Parameters: 8699c0446d6SStefano Zampini + pc - the preconditioning context 87028509bceSStefano Zampini - n_is - number of index sets defining the fields 87163602bcaSStefano Zampini . ISForDofs - array of IS describing the fields in global ordering 8729c0446d6SStefano Zampini 8739c0446d6SStefano Zampini Level: intermediate 8749c0446d6SStefano Zampini 87563602bcaSStefano Zampini Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field. 8769c0446d6SStefano Zampini 8779c0446d6SStefano Zampini .seealso: PCBDDC 8789c0446d6SStefano Zampini @*/ 8799c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 8809c0446d6SStefano Zampini { 8812b510759SStefano Zampini PetscInt i; 8829c0446d6SStefano Zampini PetscErrorCode ierr; 8839c0446d6SStefano Zampini 8849c0446d6SStefano Zampini PetscFunctionBegin; 8859c0446d6SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 88663602bcaSStefano Zampini PetscValidLogicalCollectiveInt(pc,n_is,2); 8872b510759SStefano Zampini for (i=0;i<n_is;i++) { 88863602bcaSStefano Zampini PetscCheckSameComm(pc,1,ISForDofs[i],3); 88963602bcaSStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 8902b510759SStefano Zampini } 8919c0446d6SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 8929c0446d6SStefano Zampini PetscFunctionReturn(0); 8939c0446d6SStefano Zampini } 894da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 895534831adSStefano Zampini #undef __FUNCT__ 896534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC" 897534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 898534831adSStefano Zampini /* 899534831adSStefano Zampini PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 900534831adSStefano Zampini guess if a transformation of basis approach has been selected. 9019c0446d6SStefano Zampini 902534831adSStefano Zampini Input Parameter: 903534831adSStefano Zampini + pc - the preconditioner contex 904534831adSStefano Zampini 905534831adSStefano Zampini Application Interface Routine: PCPreSolve() 906534831adSStefano Zampini 907534831adSStefano Zampini Notes: 908534831adSStefano Zampini The interface routine PCPreSolve() is not usually called directly by 909534831adSStefano Zampini the user, but instead is called by KSPSolve(). 910534831adSStefano Zampini */ 911534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 912534831adSStefano Zampini { 913534831adSStefano Zampini PetscErrorCode ierr; 914534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 915534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 916534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 917534831adSStefano Zampini Mat temp_mat; 9183972b0daSStefano Zampini IS dirIS; 9193972b0daSStefano Zampini Vec used_vec; 92082ba6b80SStefano Zampini PetscBool guess_nonzero; 921534831adSStefano Zampini 922534831adSStefano Zampini PetscFunctionBegin; 92385c4d303SStefano Zampini /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 92485c4d303SStefano Zampini if (ksp) { 92585c4d303SStefano Zampini PetscBool iscg; 92685c4d303SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 92785c4d303SStefano Zampini if (!iscg) { 92885c4d303SStefano Zampini ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 92985c4d303SStefano Zampini } 93085c4d303SStefano Zampini } 93185c4d303SStefano Zampini /* Creates parallel work vectors used in presolve */ 93262a6ff1dSStefano Zampini if (!pcbddc->original_rhs) { 93362a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 93462a6ff1dSStefano Zampini } 93562a6ff1dSStefano Zampini if (!pcbddc->temp_solution) { 93662a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 93762a6ff1dSStefano Zampini } 9383972b0daSStefano Zampini if (x) { 9393972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 9403972b0daSStefano Zampini used_vec = x; 9413972b0daSStefano Zampini } else { 9423972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 9433972b0daSStefano Zampini used_vec = pcbddc->temp_solution; 9443972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 9453972b0daSStefano Zampini } 9463972b0daSStefano Zampini /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 9473972b0daSStefano Zampini if (ksp) { 9483972b0daSStefano Zampini ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 9493972b0daSStefano Zampini if (!guess_nonzero) { 9503972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 9513972b0daSStefano Zampini } 9523972b0daSStefano Zampini } 9533308cffdSStefano Zampini 9543972b0daSStefano Zampini /* store the original rhs */ 9553972b0daSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 9563972b0daSStefano Zampini 9573972b0daSStefano Zampini /* Take into account zeroed rows -> change rhs and store solution removed */ 958785d1243SStefano Zampini /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */ 95982ba6b80SStefano Zampini ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr); 96082ba6b80SStefano Zampini if (rhs && dirIS) { 961785d1243SStefano Zampini PetscInt dirsize,i,*is_indices; 962785d1243SStefano Zampini PetscScalar *array_x,*array_diagonal; 963785d1243SStefano Zampini 9643972b0daSStefano Zampini ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 9653972b0daSStefano Zampini ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 9663972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9673972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9683972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9693972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 97082ba6b80SStefano Zampini ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr); 9713972b0daSStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 9723972b0daSStefano Zampini ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 9733972b0daSStefano Zampini ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 9742fa5cd67SKarl Rupp for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 9753972b0daSStefano Zampini ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 9763972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 9773972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 9783972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9793972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 980b76ba322SStefano Zampini 9813972b0daSStefano Zampini /* remove the computed solution from the rhs */ 9823972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 9833972b0daSStefano Zampini ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 9843972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 9853308cffdSStefano Zampini } 986b76ba322SStefano Zampini 987b76ba322SStefano Zampini /* store partially computed solution and set initial guess */ 9883972b0daSStefano Zampini if (x) { 9893972b0daSStefano Zampini ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 9903972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 99185c4d303SStefano Zampini if (pcbddc->use_exact_dirichlet_trick) { 992b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 993b76ba322SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 994b76ba322SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 995b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 996b76ba322SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 997b76ba322SStefano Zampini if (ksp) { 998b76ba322SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 999b76ba322SStefano Zampini } 1000b76ba322SStefano Zampini } 10013972b0daSStefano Zampini } 1002b76ba322SStefano Zampini 100392e3dcfbSStefano Zampini /* prepare MatMult and rhs for solver */ 1004b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix) { 1005b76ba322SStefano Zampini /* swap pointers for local matrices */ 1006b76ba322SStefano Zampini temp_mat = matis->A; 1007b76ba322SStefano Zampini matis->A = pcbddc->local_mat; 1008b76ba322SStefano Zampini pcbddc->local_mat = temp_mat; 100934a97f8cSStefano Zampini if (rhs) { /* TODO the change of basis is fragile when the change involves dofs which are not on the same subset of the interface (i.e. Toselli's for curl-curl) */ 1010b76ba322SStefano Zampini /* Get local rhs and apply transformation of basis */ 1011b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1012b76ba322SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 101334a97f8cSStefano Zampini ierr = VecPointwiseMult(pcis->vec2_B,pcis->vec2_B,pcis->D);CHKERRQ(ierr); 1014b76ba322SStefano Zampini /* from original basis to modified basis */ 1015b76ba322SStefano Zampini ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 101634a97f8cSStefano Zampini /* put back modified values into the global vec using ADD_VALUES copy mode */ 101734a97f8cSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 101834a97f8cSStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec2_B,rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1019674ae819SStefano Zampini } 102092e3dcfbSStefano Zampini } 102192e3dcfbSStefano Zampini 102292e3dcfbSStefano Zampini /* remove nullspace if present */ 10230bdf917eSStefano Zampini if (ksp && pcbddc->NullSpace) { 1024d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 1025d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 1026b76ba322SStefano Zampini } 10270bdf917eSStefano Zampini ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 1028534831adSStefano Zampini PetscFunctionReturn(0); 1029534831adSStefano Zampini } 1030534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 1031534831adSStefano Zampini #undef __FUNCT__ 1032534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC" 1033534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 1034534831adSStefano Zampini /* 1035534831adSStefano Zampini PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 1036534831adSStefano Zampini approach has been selected. Also, restores rhs to its original state. 1037534831adSStefano Zampini 1038534831adSStefano Zampini Input Parameter: 1039534831adSStefano Zampini + pc - the preconditioner contex 1040534831adSStefano Zampini 1041534831adSStefano Zampini Application Interface Routine: PCPostSolve() 1042534831adSStefano Zampini 1043534831adSStefano Zampini Notes: 1044534831adSStefano Zampini The interface routine PCPostSolve() is not usually called directly by 1045534831adSStefano Zampini the user, but instead is called by KSPSolve(). 1046534831adSStefano Zampini */ 1047534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1048534831adSStefano Zampini { 1049534831adSStefano Zampini PetscErrorCode ierr; 1050534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1051534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 1052534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 1053534831adSStefano Zampini Mat temp_mat; 1054534831adSStefano Zampini 1055534831adSStefano Zampini PetscFunctionBegin; 1056b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix) { 1057534831adSStefano Zampini /* swap pointers for local matrices */ 1058534831adSStefano Zampini temp_mat = matis->A; 1059534831adSStefano Zampini matis->A = pcbddc->local_mat; 1060534831adSStefano Zampini pcbddc->local_mat = temp_mat; 10613425bc38SStefano Zampini } 1062b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix && x) { 1063534831adSStefano Zampini /* Get Local boundary and apply transformation of basis to solution vector */ 1064534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1065534831adSStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1066534831adSStefano Zampini /* from modified basis to original basis */ 1067534831adSStefano Zampini ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 1068534831adSStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 1069534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1070534831adSStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1071534831adSStefano Zampini } 10723972b0daSStefano Zampini /* add solution removed in presolve */ 10733425bc38SStefano Zampini if (x) { 10743425bc38SStefano Zampini ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 10753425bc38SStefano Zampini } 1076fb223d50SStefano Zampini /* restore rhs to its original state */ 1077fb223d50SStefano Zampini if (rhs) { 1078fb223d50SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 1079fb223d50SStefano Zampini } 1080534831adSStefano Zampini PetscFunctionReturn(0); 1081534831adSStefano Zampini } 1082534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 108353cdbc3dSStefano Zampini #undef __FUNCT__ 108453cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC" 10850c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 10860c7d97c5SJed Brown /* 10870c7d97c5SJed Brown PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 10880c7d97c5SJed Brown by setting data structures and options. 10890c7d97c5SJed Brown 10900c7d97c5SJed Brown Input Parameter: 109153cdbc3dSStefano Zampini + pc - the preconditioner context 10920c7d97c5SJed Brown 10930c7d97c5SJed Brown Application Interface Routine: PCSetUp() 10940c7d97c5SJed Brown 10950c7d97c5SJed Brown Notes: 10960c7d97c5SJed Brown The interface routine PCSetUp() is not usually called directly by 10970c7d97c5SJed Brown the user, but instead is called by PCApply() if necessary. 10980c7d97c5SJed Brown */ 109953cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc) 11000c7d97c5SJed Brown { 11010c7d97c5SJed Brown PetscErrorCode ierr; 11020c7d97c5SJed Brown PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 1103f4ddd8eeSStefano Zampini MatNullSpace nearnullspace; 1104674ae819SStefano Zampini PetscBool computeis,computetopography,computesolvers; 1105165b64e2SStefano Zampini PetscBool new_nearnullspace_provided; 11060c7d97c5SJed Brown 11070c7d97c5SJed Brown PetscFunctionBegin; 1108f4ddd8eeSStefano Zampini /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 11093b03a366Sstefano_zampini /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1110fcd409f5SStefano Zampini Also, BDDC directly build the Dirichlet problem */ 1111f4ddd8eeSStefano Zampini /* split work */ 1112674ae819SStefano Zampini if (pc->setupcalled) { 1113674ae819SStefano Zampini computeis = PETSC_FALSE; 1114d1e9a80fSBarry Smith if (pc->flag == SAME_NONZERO_PATTERN) { 1115674ae819SStefano Zampini computetopography = PETSC_FALSE; 1116674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1117674ae819SStefano Zampini } else { /* DIFFERENT_NONZERO_PATTERN */ 1118674ae819SStefano Zampini computetopography = PETSC_TRUE; 1119674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1120674ae819SStefano Zampini } 1121674ae819SStefano Zampini } else { 1122674ae819SStefano Zampini computeis = PETSC_TRUE; 1123674ae819SStefano Zampini computetopography = PETSC_TRUE; 1124674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1125674ae819SStefano Zampini } 1126fb180af4SStefano Zampini if (pcbddc->recompute_topography) { 1127fb180af4SStefano Zampini computetopography = PETSC_TRUE; 1128fb180af4SStefano Zampini } 1129f4ddd8eeSStefano Zampini 1130f4ddd8eeSStefano Zampini /* Get stdout for dbg */ 113170cf5478SStefano Zampini if (pcbddc->dbg_flag) { 113270cf5478SStefano Zampini if (!pcbddc->dbg_viewer) { 113358a03d70SStefano Zampini pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)); 1134f4ddd8eeSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 1135f4ddd8eeSStefano Zampini } 113658a03d70SStefano Zampini ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1137f4ddd8eeSStefano Zampini } 1138f4ddd8eeSStefano Zampini 1139fcd409f5SStefano Zampini /* Set up all the "iterative substructuring" common block without computing solvers */ 1140674ae819SStefano Zampini if (computeis) { 1141fcd409f5SStefano Zampini /* HACK INTO PCIS */ 1142fcd409f5SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 1143fcd409f5SStefano Zampini pcis->computesolvers = PETSC_FALSE; 1144674ae819SStefano Zampini ierr = PCISSetUp(pc);CHKERRQ(ierr); 114539e2fb2aSStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr); 1146674ae819SStefano Zampini } 1147f4ddd8eeSStefano Zampini 1148b9b85e73SStefano Zampini /* check user defined change of basis (if any) */ 1149b9b85e73SStefano Zampini if (pcbddc->user_ChangeOfBasisMatrix) { 1150b9b85e73SStefano Zampini PC_IS* pcis= (PC_IS*)pc->data; 1151b9b85e73SStefano Zampini PetscInt n,m; 1152b9b85e73SStefano Zampini ierr = MatGetSize(pcbddc->user_ChangeOfBasisMatrix,&n,&m);CHKERRQ(ierr); 1153b9b85e73SStefano Zampini if (n != m) { 1154b9b85e73SStefano Zampini SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Change of basis matrix should be square %d != %d\n",n,m); 1155b9b85e73SStefano Zampini } else if (n != pcis->n_B && n != pcis->n) { 1156b9b85e73SStefano Zampini SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Size of change of basis matrix %d differs from boundary size %d and local size %d\n",n,pcis->n_B,pcis->n); 1157b9b85e73SStefano Zampini } 1158b9b85e73SStefano Zampini /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */ 1159b9b85e73SStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 1160b9b85e73SStefano Zampini } 1161f4ddd8eeSStefano Zampini /* Analyze interface */ 1162674ae819SStefano Zampini if (computetopography) { 1163674ae819SStefano Zampini ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 116434a97f8cSStefano Zampini /* Schurs on subsets should be reset */ 116534a97f8cSStefano Zampini if (pcbddc->deluxe_ctx) { 116634a97f8cSStefano Zampini ierr = PCBDDCSubSchursReset(pcbddc->deluxe_ctx->sub_schurs);CHKERRQ(ierr); 116734a97f8cSStefano Zampini } 1168674ae819SStefano Zampini } 1169fb8d54d4SStefano Zampini 1170f4ddd8eeSStefano Zampini /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1171fb8d54d4SStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 1172f4ddd8eeSStefano Zampini ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr); 1173f4ddd8eeSStefano Zampini if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1174f4ddd8eeSStefano Zampini if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1175f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1176f4ddd8eeSStefano Zampini } else { 1177f4ddd8eeSStefano Zampini /* determine if the two nullspaces are different (should be lightweight) */ 1178f4ddd8eeSStefano Zampini if (nearnullspace != pcbddc->onearnullspace) { 1179f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1180165b64e2SStefano Zampini } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1181f4ddd8eeSStefano Zampini PetscInt i; 1182165b64e2SStefano Zampini const Vec *nearnullvecs; 1183165b64e2SStefano Zampini PetscObjectState state; 1184165b64e2SStefano Zampini PetscInt nnsp_size; 1185165b64e2SStefano Zampini ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 1186f4ddd8eeSStefano Zampini for (i=0;i<nnsp_size;i++) { 1187f4ddd8eeSStefano Zampini ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr); 1188165b64e2SStefano Zampini if (pcbddc->onearnullvecs_state[i] != state) { 1189f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1190f4ddd8eeSStefano Zampini break; 1191f4ddd8eeSStefano Zampini } 1192f4ddd8eeSStefano Zampini } 1193f4ddd8eeSStefano Zampini } 1194f4ddd8eeSStefano Zampini } 1195f4ddd8eeSStefano Zampini } else { 1196f4ddd8eeSStefano Zampini if (!nearnullspace) { /* both nearnullspaces are null */ 1197f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 1198f4ddd8eeSStefano Zampini } else { /* nearnullspace attached later */ 1199f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1200f4ddd8eeSStefano Zampini } 1201f4ddd8eeSStefano Zampini } 1202f4ddd8eeSStefano Zampini 1203f4ddd8eeSStefano Zampini /* Setup constraints and related work vectors */ 1204727cdba6SStefano Zampini /* reset primal space flags */ 1205f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1206727cdba6SStefano Zampini pcbddc->new_primal_space_local = PETSC_FALSE; 1207fb8d54d4SStefano Zampini if (computetopography || new_nearnullspace_provided) { 1208727cdba6SStefano Zampini /* It also sets the primal space flags */ 1209674ae819SStefano Zampini ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 1210e7b262bdSStefano Zampini /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1211f4ddd8eeSStefano Zampini ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr); 121234a97f8cSStefano Zampini /* Schurs on subsets should be reset */ 121334a97f8cSStefano Zampini if (pcbddc->deluxe_ctx) { 121434a97f8cSStefano Zampini ierr = PCBDDCSubSchursReset(pcbddc->deluxe_ctx->sub_schurs);CHKERRQ(ierr); 121534a97f8cSStefano Zampini } 1216674ae819SStefano Zampini } 1217fb8d54d4SStefano Zampini 1218f4ddd8eeSStefano Zampini if (computesolvers || pcbddc->new_primal_space) { 121999cc7994SStefano Zampini /* Create coarse and local stuffs */ 122099cc7994SStefano Zampini ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 122134a97f8cSStefano Zampini /* Create scaling operators */ 1222674ae819SStefano Zampini ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 12230c7d97c5SJed Brown } 122470cf5478SStefano Zampini 122558a03d70SStefano Zampini if (pcbddc->dbg_flag) { 122658a03d70SStefano Zampini ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 12272b510759SStefano Zampini } 12280c7d97c5SJed Brown PetscFunctionReturn(0); 12290c7d97c5SJed Brown } 12300c7d97c5SJed Brown 12310c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 12320c7d97c5SJed Brown /* 123350efa1b5SStefano Zampini PCApply_BDDC - Applies the BDDC operator to a vector. 12340c7d97c5SJed Brown 12350c7d97c5SJed Brown Input Parameters: 12360c7d97c5SJed Brown . pc - the preconditioner context 12370c7d97c5SJed Brown . r - input vector (global) 12380c7d97c5SJed Brown 12390c7d97c5SJed Brown Output Parameter: 12400c7d97c5SJed Brown . z - output vector (global) 12410c7d97c5SJed Brown 12420c7d97c5SJed Brown Application Interface Routine: PCApply() 12430c7d97c5SJed Brown */ 12440c7d97c5SJed Brown #undef __FUNCT__ 12450c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC" 124653cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 12470c7d97c5SJed Brown { 12480c7d97c5SJed Brown PC_IS *pcis = (PC_IS*)(pc->data); 12490c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 12500c7d97c5SJed Brown PetscErrorCode ierr; 12513b03a366Sstefano_zampini const PetscScalar one = 1.0; 12523b03a366Sstefano_zampini const PetscScalar m_one = -1.0; 12532617d88aSStefano Zampini const PetscScalar zero = 0.0; 12540c7d97c5SJed Brown 12550c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN 12560c7d97c5SJed Brown NN interface preconditioner changed to BDDC 12578eeda7d8SStefano Zampini Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 12580c7d97c5SJed Brown 12590c7d97c5SJed Brown PetscFunctionBegin; 126085c4d303SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 12610c7d97c5SJed Brown /* First Dirichlet solve */ 12620c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12630c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 126453cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 12650c7d97c5SJed Brown /* 12660c7d97c5SJed Brown Assembling right hand side for BDDC operator 1267674ae819SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 1268674ae819SStefano Zampini - pcis->vec1_B the interface part of the global vector z 12690c7d97c5SJed Brown */ 12700c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 12710c7d97c5SJed Brown ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 12728eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 12730c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 12740c7d97c5SJed Brown ierr = VecCopy(r,z);CHKERRQ(ierr); 12750c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12760c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1277674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1278b76ba322SStefano Zampini } else { 12790bdf917eSStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 1280b76ba322SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 1281674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1282b76ba322SStefano Zampini } 1283b76ba322SStefano Zampini 12842617d88aSStefano Zampini /* Apply interface preconditioner 12852617d88aSStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1286dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr); 12872617d88aSStefano Zampini 1288674ae819SStefano Zampini /* Apply transpose of partition of unity operator */ 1289674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 12900c7d97c5SJed Brown 12913b03a366Sstefano_zampini /* Second Dirichlet solve and assembling of output */ 12920c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12930c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12940c7d97c5SJed Brown ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 12958eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1296df187020SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1297df187020SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1298df187020SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1299df187020SStefano Zampini ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 13000c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 13010c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 13020c7d97c5SJed Brown PetscFunctionReturn(0); 13030c7d97c5SJed Brown } 130450efa1b5SStefano Zampini 130550efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */ 130650efa1b5SStefano Zampini /* 130750efa1b5SStefano Zampini PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector. 130850efa1b5SStefano Zampini 130950efa1b5SStefano Zampini Input Parameters: 131050efa1b5SStefano Zampini . pc - the preconditioner context 131150efa1b5SStefano Zampini . r - input vector (global) 131250efa1b5SStefano Zampini 131350efa1b5SStefano Zampini Output Parameter: 131450efa1b5SStefano Zampini . z - output vector (global) 131550efa1b5SStefano Zampini 131650efa1b5SStefano Zampini Application Interface Routine: PCApplyTranspose() 131750efa1b5SStefano Zampini */ 131850efa1b5SStefano Zampini #undef __FUNCT__ 131950efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC" 132050efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z) 132150efa1b5SStefano Zampini { 132250efa1b5SStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 132350efa1b5SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 132450efa1b5SStefano Zampini PetscErrorCode ierr; 132550efa1b5SStefano Zampini const PetscScalar one = 1.0; 132650efa1b5SStefano Zampini const PetscScalar m_one = -1.0; 132750efa1b5SStefano Zampini const PetscScalar zero = 0.0; 132850efa1b5SStefano Zampini 132950efa1b5SStefano Zampini PetscFunctionBegin; 133050efa1b5SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 133150efa1b5SStefano Zampini /* First Dirichlet solve */ 133250efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 133350efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 133450efa1b5SStefano Zampini ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 133550efa1b5SStefano Zampini /* 133650efa1b5SStefano Zampini Assembling right hand side for BDDC operator 133750efa1b5SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 133850efa1b5SStefano Zampini - pcis->vec1_B the interface part of the global vector z 133950efa1b5SStefano Zampini */ 134050efa1b5SStefano Zampini ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 134150efa1b5SStefano Zampini ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 134250efa1b5SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 134350efa1b5SStefano Zampini ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 134450efa1b5SStefano Zampini ierr = VecCopy(r,z);CHKERRQ(ierr); 134550efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 134650efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 134750efa1b5SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 134850efa1b5SStefano Zampini } else { 134950efa1b5SStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 135050efa1b5SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 135150efa1b5SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 135250efa1b5SStefano Zampini } 135350efa1b5SStefano Zampini 135450efa1b5SStefano Zampini /* Apply interface preconditioner 135550efa1b5SStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1356dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr); 135750efa1b5SStefano Zampini 135850efa1b5SStefano Zampini /* Apply transpose of partition of unity operator */ 135950efa1b5SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 136050efa1b5SStefano Zampini 136150efa1b5SStefano Zampini /* Second Dirichlet solve and assembling of output */ 136250efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 136350efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 136450efa1b5SStefano Zampini ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 136550efa1b5SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1366b0147a47SStefano Zampini ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1367b0147a47SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1368b0147a47SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1369b0147a47SStefano Zampini ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 137050efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 137150efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 137250efa1b5SStefano Zampini PetscFunctionReturn(0); 137350efa1b5SStefano Zampini } 1374da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 1375674ae819SStefano Zampini 1376da1bb401SStefano Zampini #undef __FUNCT__ 1377da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC" 1378da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc) 1379da1bb401SStefano Zampini { 1380da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1381da1bb401SStefano Zampini PetscErrorCode ierr; 1382da1bb401SStefano Zampini 1383da1bb401SStefano Zampini PetscFunctionBegin; 1384da1bb401SStefano Zampini /* free data created by PCIS */ 1385da1bb401SStefano Zampini ierr = PCISDestroy(pc);CHKERRQ(ierr); 1386674ae819SStefano Zampini /* free BDDC custom data */ 1387674ae819SStefano Zampini ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 1388674ae819SStefano Zampini /* destroy objects related to topography */ 1389674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 1390674ae819SStefano Zampini /* free allocated graph structure */ 1391da1bb401SStefano Zampini ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 139234a97f8cSStefano Zampini /* destroy objects for scaling operator */ 1393674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 139434a97f8cSStefano Zampini ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr); 1395674ae819SStefano Zampini /* free solvers stuff */ 1396674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 139762a6ff1dSStefano Zampini /* free global vectors needed in presolve */ 139862a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 139962a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 140039e2fb2aSStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr); 14013425bc38SStefano Zampini /* remove functions */ 1402b9b85e73SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisLocalMat_C",NULL);CHKERRQ(ierr); 1403674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 1404bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 14052b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 1406b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 14072b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 1408bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 1409bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 141082ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1411bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 141282ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1413bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 141482ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1415bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1416785d1243SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1417bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 141863602bcaSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr); 1419bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 1420bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 1421bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 1422bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 1423674ae819SStefano Zampini /* Free the private data structure */ 1424674ae819SStefano Zampini ierr = PetscFree(pc->data);CHKERRQ(ierr); 1425da1bb401SStefano Zampini PetscFunctionReturn(0); 1426da1bb401SStefano Zampini } 14273425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 14281e6b0712SBarry Smith 14293425bc38SStefano Zampini #undef __FUNCT__ 14303425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 14313425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 14323425bc38SStefano Zampini { 1433674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 14343425bc38SStefano Zampini PC_IS* pcis; 14353425bc38SStefano Zampini PC_BDDC* pcbddc; 14363425bc38SStefano Zampini PetscErrorCode ierr; 14370c7d97c5SJed Brown 14383425bc38SStefano Zampini PetscFunctionBegin; 14393425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 14403425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 14413425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 14423425bc38SStefano Zampini 14433425bc38SStefano Zampini /* change of basis for physical rhs if needed 14443425bc38SStefano Zampini It also changes the rhs in case of dirichlet boundaries */ 14453308cffdSStefano Zampini ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 14463425bc38SStefano Zampini /* store vectors for computation of fetidp final solution */ 14473425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14483425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1449fb223d50SStefano Zampini /* scale rhs since it should be unassembled */ 1450fb223d50SStefano Zampini /* TODO use counter scaling? (also below) */ 14513425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14523425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1453674ae819SStefano Zampini /* Apply partition of unity */ 14543425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1455674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 14568eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 14573425bc38SStefano Zampini /* compute partially subassembled Schur complement right-hand side */ 14583425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 14593425bc38SStefano Zampini ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 14603425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 14613425bc38SStefano Zampini ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 14623425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14633425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1464674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 14653425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14663425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14673425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 14683425bc38SStefano Zampini } 14693425bc38SStefano Zampini /* BDDC rhs */ 14703425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 14718eeda7d8SStefano Zampini if (pcbddc->switch_static) { 14723425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 14733425bc38SStefano Zampini } 14743425bc38SStefano Zampini /* apply BDDC */ 1475dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 14763425bc38SStefano Zampini /* Application of B_delta and assembling of rhs for fetidp fluxes */ 14773425bc38SStefano Zampini ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 14783425bc38SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 14793425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14803425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14813425bc38SStefano Zampini /* restore original rhs */ 14823425bc38SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 14833425bc38SStefano Zampini PetscFunctionReturn(0); 14843425bc38SStefano Zampini } 14851e6b0712SBarry Smith 14863425bc38SStefano Zampini #undef __FUNCT__ 14873425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 14883425bc38SStefano Zampini /*@ 148928509bceSStefano Zampini PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system 14903425bc38SStefano Zampini 14913425bc38SStefano Zampini Collective 14923425bc38SStefano Zampini 14933425bc38SStefano Zampini Input Parameters: 149428509bceSStefano Zampini + fetidp_mat - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators 149528509bceSStefano Zampini . standard_rhs - the right-hand side for your linear system 14963425bc38SStefano Zampini 14973425bc38SStefano Zampini Output Parameters: 149828509bceSStefano Zampini - fetidp_flux_rhs - the right-hand side for the FETIDP linear system 14993425bc38SStefano Zampini 15003425bc38SStefano Zampini Level: developer 15013425bc38SStefano Zampini 15023425bc38SStefano Zampini Notes: 15033425bc38SStefano Zampini 150428509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 15053425bc38SStefano Zampini @*/ 15063425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 15073425bc38SStefano Zampini { 1508674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 15093425bc38SStefano Zampini PetscErrorCode ierr; 15103425bc38SStefano Zampini 15113425bc38SStefano Zampini PetscFunctionBegin; 15123425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15133425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 15143425bc38SStefano Zampini PetscFunctionReturn(0); 15153425bc38SStefano Zampini } 15163425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 15171e6b0712SBarry Smith 15183425bc38SStefano Zampini #undef __FUNCT__ 15193425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 15203425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 15213425bc38SStefano Zampini { 1522674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 15233425bc38SStefano Zampini PC_IS* pcis; 15243425bc38SStefano Zampini PC_BDDC* pcbddc; 15253425bc38SStefano Zampini PetscErrorCode ierr; 15263425bc38SStefano Zampini 15273425bc38SStefano Zampini PetscFunctionBegin; 15283425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15293425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 15303425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 15313425bc38SStefano Zampini 15323425bc38SStefano Zampini /* apply B_delta^T */ 15333425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15343425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15353425bc38SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 15363425bc38SStefano Zampini /* compute rhs for BDDC application */ 15373425bc38SStefano Zampini ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 15388eeda7d8SStefano Zampini if (pcbddc->switch_static) { 15393425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 15403425bc38SStefano Zampini } 15413425bc38SStefano Zampini /* apply BDDC */ 1542dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 15433425bc38SStefano Zampini /* put values into standard global vector */ 15443425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15453425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15468eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 15473425bc38SStefano Zampini /* compute values into the interior if solved for the partially subassembled Schur complement */ 15483425bc38SStefano Zampini ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 15493425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 15503425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 15513425bc38SStefano Zampini } 15523425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15533425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15543425bc38SStefano Zampini /* final change of basis if needed 15553425bc38SStefano Zampini Is also sums the dirichlet part removed during RHS assembling */ 15563308cffdSStefano Zampini ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 15573425bc38SStefano Zampini PetscFunctionReturn(0); 15583425bc38SStefano Zampini } 15591e6b0712SBarry Smith 15603425bc38SStefano Zampini #undef __FUNCT__ 15613425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 15623425bc38SStefano Zampini /*@ 156328509bceSStefano Zampini PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system 15643425bc38SStefano Zampini 15653425bc38SStefano Zampini Collective 15663425bc38SStefano Zampini 15673425bc38SStefano Zampini Input Parameters: 156828509bceSStefano Zampini + fetidp_mat - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators 156928509bceSStefano Zampini . fetidp_flux_sol - the solution of the FETIDP linear system 15703425bc38SStefano Zampini 15713425bc38SStefano Zampini Output Parameters: 157228509bceSStefano Zampini - standard_sol - the solution defined on the physical domain 15733425bc38SStefano Zampini 15743425bc38SStefano Zampini Level: developer 15753425bc38SStefano Zampini 15763425bc38SStefano Zampini Notes: 15773425bc38SStefano Zampini 157828509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 15793425bc38SStefano Zampini @*/ 15803425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 15813425bc38SStefano Zampini { 1582674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 15833425bc38SStefano Zampini PetscErrorCode ierr; 15843425bc38SStefano Zampini 15853425bc38SStefano Zampini PetscFunctionBegin; 15863425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15873425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 15883425bc38SStefano Zampini PetscFunctionReturn(0); 15893425bc38SStefano Zampini } 15903425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 15911e6b0712SBarry Smith 1592f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1593edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec); 1594f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1595f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1596edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec); 1597f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1598674ae819SStefano Zampini 15993425bc38SStefano Zampini #undef __FUNCT__ 16003425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 16013425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 16023425bc38SStefano Zampini { 1603674ae819SStefano Zampini 1604674ae819SStefano Zampini FETIDPMat_ctx fetidpmat_ctx; 16053425bc38SStefano Zampini Mat newmat; 1606674ae819SStefano Zampini FETIDPPC_ctx fetidppc_ctx; 16073425bc38SStefano Zampini PC newpc; 1608ce94432eSBarry Smith MPI_Comm comm; 16093425bc38SStefano Zampini PetscErrorCode ierr; 16103425bc38SStefano Zampini 16113425bc38SStefano Zampini PetscFunctionBegin; 1612ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 16133425bc38SStefano Zampini /* FETIDP linear matrix */ 16143425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 16153425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 16163425bc38SStefano Zampini ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 16173425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 1618edf7251bSStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr); 16193425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 16203425bc38SStefano Zampini ierr = MatSetUp(newmat);CHKERRQ(ierr); 16213425bc38SStefano Zampini /* FETIDP preconditioner */ 16223425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 16233425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 16243425bc38SStefano Zampini ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 16253425bc38SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 16263425bc38SStefano Zampini ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 16273425bc38SStefano Zampini ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 1628edf7251bSStefano Zampini ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr); 16293425bc38SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 163023ee1639SBarry Smith ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr); 16313425bc38SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 16323425bc38SStefano Zampini /* return pointers for objects created */ 16333425bc38SStefano Zampini *fetidp_mat=newmat; 16343425bc38SStefano Zampini *fetidp_pc=newpc; 16353425bc38SStefano Zampini PetscFunctionReturn(0); 16363425bc38SStefano Zampini } 16371e6b0712SBarry Smith 16383425bc38SStefano Zampini #undef __FUNCT__ 16393425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 16403425bc38SStefano Zampini /*@ 164128509bceSStefano Zampini PCBDDCCreateFETIDPOperators - Create operators for FETIDP 16423425bc38SStefano Zampini 16433425bc38SStefano Zampini Collective 16443425bc38SStefano Zampini 16453425bc38SStefano Zampini Input Parameters: 164628509bceSStefano Zampini + pc - the BDDC preconditioning context already setup 164728509bceSStefano Zampini 164828509bceSStefano Zampini Output Parameters: 164928509bceSStefano Zampini . fetidp_mat - shell FETIDP matrix object 165028509bceSStefano Zampini . fetidp_pc - shell Dirichlet preconditioner for FETIDP matrix 165128509bceSStefano Zampini 165228509bceSStefano Zampini Options Database Keys: 165328509bceSStefano Zampini - -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers 16543425bc38SStefano Zampini 16553425bc38SStefano Zampini Level: developer 16563425bc38SStefano Zampini 16573425bc38SStefano Zampini Notes: 165828509bceSStefano Zampini Currently the only operation provided for FETIDP matrix is MatMult 16593425bc38SStefano Zampini 16603425bc38SStefano Zampini .seealso: PCBDDC 16613425bc38SStefano Zampini @*/ 16623425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 16633425bc38SStefano Zampini { 16643425bc38SStefano Zampini PetscErrorCode ierr; 16653425bc38SStefano Zampini 16663425bc38SStefano Zampini PetscFunctionBegin; 16673425bc38SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16683425bc38SStefano Zampini if (pc->setupcalled) { 1669516d51a7SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1670f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 16713425bc38SStefano Zampini PetscFunctionReturn(0); 16723425bc38SStefano Zampini } 16730c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 1674da1bb401SStefano Zampini /*MC 1675da1bb401SStefano Zampini PCBDDC - Balancing Domain Decomposition by Constraints. 16760c7d97c5SJed Brown 167728509bceSStefano Zampini An implementation of the BDDC preconditioner based on 167828509bceSStefano Zampini 167928509bceSStefano Zampini .vb 168028509bceSStefano Zampini [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 168128509bceSStefano 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 168228509bceSStefano Zampini [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 168328509bceSStefano Zampini .ve 168428509bceSStefano Zampini 168528509bceSStefano Zampini The matrix to be preconditioned (Pmat) must be of type MATIS. 168628509bceSStefano Zampini 1687b6fdb6dfSStefano Zampini Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 168828509bceSStefano Zampini 168928509bceSStefano Zampini It also works with unsymmetric and indefinite problems. 169028509bceSStefano Zampini 1691b6fdb6dfSStefano 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. 1692b6fdb6dfSStefano Zampini 169328509bceSStefano Zampini Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace 169428509bceSStefano Zampini 16956a818285SBarry Smith 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() 169628509bceSStefano Zampini 16976a818285SBarry Smith Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). 169828509bceSStefano Zampini 1699b6fdb6dfSStefano 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. 170028509bceSStefano Zampini 170128509bceSStefano Zampini The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object. 170228509bceSStefano Zampini 1703da1bb401SStefano Zampini Options Database Keys: 170428509bceSStefano Zampini 170528509bceSStefano Zampini . -pc_bddc_use_vertices <1> - use or not vertices in primal space 170628509bceSStefano Zampini . -pc_bddc_use_edges <1> - use or not edges in primal space 1707b6fdb6dfSStefano Zampini . -pc_bddc_use_faces <0> - use or not faces in primal space 170828509bceSStefano Zampini . -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only) 170928509bceSStefano Zampini . -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested 171028509bceSStefano Zampini . -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1]) 171128509bceSStefano Zampini . -pc_bddc_levels <0> - maximum number of levels for multilevel 171228509bceSStefano Zampini . -pc_bddc_coarsening_ratio - H/h ratio at the coarser level 171328509bceSStefano Zampini - -pc_bddc_check_level <0> - set verbosity level of debugging output 171428509bceSStefano Zampini 171528509bceSStefano Zampini Options for Dirichlet, Neumann or coarse solver can be set with 171628509bceSStefano Zampini .vb 171728509bceSStefano Zampini -pc_bddc_dirichlet_ 171828509bceSStefano Zampini -pc_bddc_neumann_ 171928509bceSStefano Zampini -pc_bddc_coarse_ 172028509bceSStefano Zampini .ve 172128509bceSStefano Zampini e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg 172228509bceSStefano Zampini 172328509bceSStefano Zampini When using a multilevel approach, solvers' options at the N-th level can be specified as 172428509bceSStefano Zampini .vb 1725312be037SStefano Zampini -pc_bddc_dirichlet_lN_ 1726312be037SStefano Zampini -pc_bddc_neumann_lN_ 1727312be037SStefano Zampini -pc_bddc_coarse_lN_ 172828509bceSStefano Zampini .ve 1729312be037SStefano Zampini Note that level number ranges from the finest 0 to the coarsest N. 1730da1bb401SStefano Zampini 1731da1bb401SStefano Zampini Level: intermediate 1732da1bb401SStefano Zampini 1733b6fdb6dfSStefano Zampini Developer notes: 1734da1bb401SStefano Zampini 173528509bceSStefano Zampini New deluxe scaling operator will be available soon. 1736da1bb401SStefano Zampini 1737da1bb401SStefano Zampini Contributed by Stefano Zampini 1738da1bb401SStefano Zampini 1739da1bb401SStefano Zampini .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1740da1bb401SStefano Zampini M*/ 1741b2573a8aSBarry Smith 1742da1bb401SStefano Zampini #undef __FUNCT__ 1743da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC" 17448cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1745da1bb401SStefano Zampini { 1746da1bb401SStefano Zampini PetscErrorCode ierr; 1747da1bb401SStefano Zampini PC_BDDC *pcbddc; 1748da1bb401SStefano Zampini 1749da1bb401SStefano Zampini PetscFunctionBegin; 1750da1bb401SStefano Zampini /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1751b00a9115SJed Brown ierr = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr); 1752da1bb401SStefano Zampini pc->data = (void*)pcbddc; 1753da1bb401SStefano Zampini 1754da1bb401SStefano Zampini /* create PCIS data structure */ 1755da1bb401SStefano Zampini ierr = PCISCreate(pc);CHKERRQ(ierr); 1756da1bb401SStefano Zampini 175747d04d0dSStefano Zampini /* BDDC customization */ 175808a5cf49SStefano Zampini pcbddc->use_local_adj = PETSC_TRUE; 175947d04d0dSStefano Zampini pcbddc->use_vertices = PETSC_TRUE; 176047d04d0dSStefano Zampini pcbddc->use_edges = PETSC_TRUE; 176147d04d0dSStefano Zampini pcbddc->use_faces = PETSC_FALSE; 176247d04d0dSStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 176347d04d0dSStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 176447d04d0dSStefano Zampini pcbddc->switch_static = PETSC_FALSE; 176547d04d0dSStefano Zampini pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 176647d04d0dSStefano Zampini pcbddc->dbg_flag = 0; 1767b9d89cd5SStefano Zampini /* private */ 1768b9d89cd5SStefano Zampini pcbddc->issym = PETSC_FALSE; 176939e2fb2aSStefano Zampini pcbddc->BtoNmap = 0; 1770727cdba6SStefano Zampini pcbddc->local_primal_size = 0; 1771e9189074SStefano Zampini pcbddc->n_vertices = 0; 1772e9189074SStefano Zampini pcbddc->n_actual_vertices = 0; 1773e9189074SStefano Zampini pcbddc->n_constraints = 0; 1774727cdba6SStefano Zampini pcbddc->primal_indices_local_idxs = 0; 1775fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_FALSE; 177668457ee5SStefano Zampini pcbddc->coarse_size = -1; 1777f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1778727cdba6SStefano Zampini pcbddc->new_primal_space_local = PETSC_FALSE; 1779f4ddd8eeSStefano Zampini pcbddc->global_primal_indices = 0; 1780f4ddd8eeSStefano Zampini pcbddc->onearnullspace = 0; 1781f4ddd8eeSStefano Zampini pcbddc->onearnullvecs_state = 0; 1782674ae819SStefano Zampini pcbddc->user_primal_vertices = 0; 17830bdf917eSStefano Zampini pcbddc->NullSpace = 0; 17843972b0daSStefano Zampini pcbddc->temp_solution = 0; 1785534831adSStefano Zampini pcbddc->original_rhs = 0; 1786534831adSStefano Zampini pcbddc->local_mat = 0; 1787534831adSStefano Zampini pcbddc->ChangeOfBasisMatrix = 0; 1788b9b85e73SStefano Zampini pcbddc->user_ChangeOfBasisMatrix = 0; 1789da1bb401SStefano Zampini pcbddc->coarse_vec = 0; 1790da1bb401SStefano Zampini pcbddc->coarse_rhs = 0; 1791da1bb401SStefano Zampini pcbddc->coarse_ksp = 0; 1792da1bb401SStefano Zampini pcbddc->coarse_phi_B = 0; 1793da1bb401SStefano Zampini pcbddc->coarse_phi_D = 0; 179415aaf578SStefano Zampini pcbddc->coarse_psi_B = 0; 179515aaf578SStefano Zampini pcbddc->coarse_psi_D = 0; 1796da1bb401SStefano Zampini pcbddc->vec1_P = 0; 1797da1bb401SStefano Zampini pcbddc->vec1_R = 0; 1798da1bb401SStefano Zampini pcbddc->vec2_R = 0; 1799da1bb401SStefano Zampini pcbddc->local_auxmat1 = 0; 1800da1bb401SStefano Zampini pcbddc->local_auxmat2 = 0; 1801da1bb401SStefano Zampini pcbddc->R_to_B = 0; 1802da1bb401SStefano Zampini pcbddc->R_to_D = 0; 1803da1bb401SStefano Zampini pcbddc->ksp_D = 0; 1804da1bb401SStefano Zampini pcbddc->ksp_R = 0; 1805da1bb401SStefano Zampini pcbddc->NeumannBoundaries = 0; 1806785d1243SStefano Zampini pcbddc->NeumannBoundariesLocal = 0; 1807785d1243SStefano Zampini pcbddc->DirichletBoundaries = 0; 1808785d1243SStefano Zampini pcbddc->DirichletBoundariesLocal = 0; 180960d44989SStefano Zampini pcbddc->user_provided_isfordofs = PETSC_FALSE; 181060d44989SStefano Zampini pcbddc->n_ISForDofs = 0; 181163602bcaSStefano Zampini pcbddc->n_ISForDofsLocal = 0; 1812da1bb401SStefano Zampini pcbddc->ISForDofs = 0; 181363602bcaSStefano Zampini pcbddc->ISForDofsLocal = 0; 1814da1bb401SStefano Zampini pcbddc->ConstraintMatrix = 0; 181585c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 181647d04d0dSStefano Zampini pcbddc->coarse_loc_to_glob = 0; 181747d04d0dSStefano Zampini pcbddc->coarsening_ratio = 8; 18184fad6a16SStefano Zampini pcbddc->current_level = 0; 18192b510759SStefano Zampini pcbddc->max_levels = 0; 1820323d395dSStefano Zampini pcbddc->use_coarse_estimates = PETSC_FALSE; 1821f3bde8b3SStefano Zampini pcbddc->coarse_subassembling = 0; 1822323d395dSStefano Zampini pcbddc->coarse_subassembling_init = 0; 1823da1bb401SStefano Zampini 1824674ae819SStefano Zampini /* create local graph structure */ 1825674ae819SStefano Zampini ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1826674ae819SStefano Zampini 1827674ae819SStefano Zampini /* scaling */ 1828674ae819SStefano Zampini pcbddc->work_scaling = 0; 182934a97f8cSStefano Zampini pcbddc->use_deluxe_scaling = PETSC_FALSE; 183034a97f8cSStefano Zampini pcbddc->deluxe_threshold = 1; 183134a97f8cSStefano Zampini pcbddc->deluxe_rebuild = PETSC_FALSE; 183234a97f8cSStefano Zampini pcbddc->deluxe_layers = -1; 183334a97f8cSStefano Zampini pcbddc->deluxe_use_useradj = PETSC_FALSE; 183434a97f8cSStefano Zampini pcbddc->deluxe_compute_rowadj = PETSC_TRUE; 1835da1bb401SStefano Zampini 1836da1bb401SStefano Zampini /* function pointers */ 1837da1bb401SStefano Zampini pc->ops->apply = PCApply_BDDC; 183893bd9ae7SStefano Zampini pc->ops->applytranspose = PCApplyTranspose_BDDC; 1839da1bb401SStefano Zampini pc->ops->setup = PCSetUp_BDDC; 1840da1bb401SStefano Zampini pc->ops->destroy = PCDestroy_BDDC; 1841da1bb401SStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1842da1bb401SStefano Zampini pc->ops->view = 0; 1843da1bb401SStefano Zampini pc->ops->applyrichardson = 0; 1844da1bb401SStefano Zampini pc->ops->applysymmetricleft = 0; 1845da1bb401SStefano Zampini pc->ops->applysymmetricright = 0; 1846534831adSStefano Zampini pc->ops->presolve = PCPreSolve_BDDC; 1847534831adSStefano Zampini pc->ops->postsolve = PCPostSolve_BDDC; 1848da1bb401SStefano Zampini 1849da1bb401SStefano Zampini /* composing function */ 1850b9b85e73SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisLocalMat_C",PCBDDCSetChangeOfBasisLocalMat_BDDC);CHKERRQ(ierr); 1851674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1852bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 18532b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1854b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 18552b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1856bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1857bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 185882ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1859bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 186082ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1861bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 186282ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1863bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 186482ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1865bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 186663602bcaSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr); 1867bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1868bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1869bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1870bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1871da1bb401SStefano Zampini PetscFunctionReturn(0); 1872da1bb401SStefano Zampini } 18733425bc38SStefano Zampini 1874