153cdbc3dSStefano Zampini /* TODOLIST 2eb97c9d2SStefano Zampini 3eb97c9d2SStefano Zampini ConstraintsSetup 4eb97c9d2SStefano Zampini - tolerances for constraints as an option (take care of single precision!) 5c1c8e736SStefano Zampini - Can MAT_IGNORE_ZERO_ENTRIES be used for Constraints Matrix? 6eb97c9d2SStefano Zampini 7eb97c9d2SStefano Zampini Solvers 8eb97c9d2SStefano Zampini - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers) 9eb97c9d2SStefano Zampini - Propagate ksp prefixes for solvers to mat objects? 10eb97c9d2SStefano Zampini - Propagate nearnullspace info among levels 11eb97c9d2SStefano Zampini 12eb97c9d2SStefano Zampini User interface 13eb97c9d2SStefano Zampini - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access) 14eb97c9d2SStefano Zampini - DofSplitting and DM attached to pc? 15eb97c9d2SStefano Zampini 16eb97c9d2SStefano Zampini Debugging output 17eb97c9d2SStefano Zampini - Better management of verbosity levels of debugging output 18eb97c9d2SStefano Zampini 19eb97c9d2SStefano Zampini Build 20eb97c9d2SStefano Zampini - make runexe59 21eb97c9d2SStefano Zampini 22eb97c9d2SStefano Zampini Extra 23eb97c9d2SStefano Zampini - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)? 24eb97c9d2SStefano Zampini - add support for computing h,H and related using coordinates? 25c0b83709SStefano Zampini - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap) 26eb97c9d2SStefano Zampini - Better management in PCIS code 27eb97c9d2SStefano Zampini - BDDC with MG framework? 28eb97c9d2SStefano Zampini 29eb97c9d2SStefano Zampini FETIDP 30eb97c9d2SStefano Zampini - Move FETIDP code to its own classes 31eb97c9d2SStefano Zampini 32eb97c9d2SStefano Zampini MATIS related operations contained in BDDC code 33eb97c9d2SStefano Zampini - Provide general case for subassembling 3439e2fb2aSStefano Zampini - Preallocation routines in MatISGetMPIAXAIJ 35eb97c9d2SStefano Zampini 3653cdbc3dSStefano Zampini */ 370c7d97c5SJed Brown 3853cdbc3dSStefano Zampini /* ---------------------------------------------------------------------------------------------------------------------------------------------- 390c7d97c5SJed Brown Implementation of BDDC preconditioner based on: 400c7d97c5SJed Brown C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 4153cdbc3dSStefano Zampini ---------------------------------------------------------------------------------------------------------------------------------------------- */ 4253cdbc3dSStefano Zampini 43ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 44ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 453b03a366Sstefano_zampini #include <petscblaslapack.h> 46674ae819SStefano Zampini 470c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 480c7d97c5SJed Brown #undef __FUNCT__ 490c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC" 500c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc) 510c7d97c5SJed Brown { 520c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 530c7d97c5SJed Brown PetscErrorCode ierr; 540c7d97c5SJed Brown 550c7d97c5SJed Brown PetscFunctionBegin; 560c7d97c5SJed Brown ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 578eeda7d8SStefano Zampini /* Verbose debugging */ 588eeda7d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 598eeda7d8SStefano Zampini /* Primal space cumstomization */ 60*08a5cf49SStefano 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); 618eeda7d8SStefano 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); 628eeda7d8SStefano 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); 638eeda7d8SStefano 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); 648eeda7d8SStefano Zampini /* Change of basis */ 658eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr); 668eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr); 67674ae819SStefano Zampini if (!pcbddc->use_change_of_basis) { 68674ae819SStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 69674ae819SStefano Zampini } 708eeda7d8SStefano Zampini /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 718eeda7d8SStefano 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); 720298fd71SBarry Smith ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 732b510759SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 74323d395dSStefano 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); 75674ae819SStefano 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); 760c7d97c5SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 770c7d97c5SJed Brown PetscFunctionReturn(0); 780c7d97c5SJed Brown } 790c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 80674ae819SStefano Zampini #undef __FUNCT__ 81674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 82674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 83674ae819SStefano Zampini { 84674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 85674ae819SStefano Zampini PetscErrorCode ierr; 861e6b0712SBarry Smith 87674ae819SStefano Zampini PetscFunctionBegin; 88674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 89674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 90674ae819SStefano Zampini pcbddc->user_primal_vertices = PrimalVertices; 91674ae819SStefano Zampini PetscFunctionReturn(0); 92674ae819SStefano Zampini } 93674ae819SStefano Zampini #undef __FUNCT__ 94674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 95674ae819SStefano Zampini /*@ 9628509bceSStefano Zampini PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 97674ae819SStefano Zampini 98674ae819SStefano Zampini Not collective 99674ae819SStefano Zampini 100674ae819SStefano Zampini Input Parameters: 101674ae819SStefano Zampini + pc - the preconditioning context 10228509bceSStefano Zampini - PrimalVertices - index set of primal vertices in local numbering 103674ae819SStefano Zampini 104674ae819SStefano Zampini Level: intermediate 105674ae819SStefano Zampini 106674ae819SStefano Zampini Notes: 107674ae819SStefano Zampini 108674ae819SStefano Zampini .seealso: PCBDDC 109674ae819SStefano Zampini @*/ 110674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 111674ae819SStefano Zampini { 112674ae819SStefano Zampini PetscErrorCode ierr; 113674ae819SStefano Zampini 114674ae819SStefano Zampini PetscFunctionBegin; 115674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 116674ae819SStefano Zampini PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 117674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 118674ae819SStefano Zampini PetscFunctionReturn(0); 119674ae819SStefano Zampini } 120674ae819SStefano Zampini /* -------------------------------------------------------------------------- */ 1210c7d97c5SJed Brown #undef __FUNCT__ 1224fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 1234fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 1244fad6a16SStefano Zampini { 1254fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1264fad6a16SStefano Zampini 1274fad6a16SStefano Zampini PetscFunctionBegin; 1284fad6a16SStefano Zampini pcbddc->coarsening_ratio = k; 1294fad6a16SStefano Zampini PetscFunctionReturn(0); 1304fad6a16SStefano Zampini } 1311e6b0712SBarry Smith 1324fad6a16SStefano Zampini #undef __FUNCT__ 1334fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio" 1344fad6a16SStefano Zampini /*@ 13528509bceSStefano Zampini PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 1364fad6a16SStefano Zampini 1374fad6a16SStefano Zampini Logically collective on PC 1384fad6a16SStefano Zampini 1394fad6a16SStefano Zampini Input Parameters: 1404fad6a16SStefano Zampini + pc - the preconditioning context 14128509bceSStefano Zampini - k - coarsening ratio (H/h at the coarser level) 1424fad6a16SStefano Zampini 14328509bceSStefano Zampini Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 1444fad6a16SStefano Zampini 1454fad6a16SStefano Zampini Level: intermediate 1464fad6a16SStefano Zampini 1474fad6a16SStefano Zampini Notes: 1484fad6a16SStefano Zampini 1494fad6a16SStefano Zampini .seealso: PCBDDC 1504fad6a16SStefano Zampini @*/ 1514fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 1524fad6a16SStefano Zampini { 1534fad6a16SStefano Zampini PetscErrorCode ierr; 1544fad6a16SStefano Zampini 1554fad6a16SStefano Zampini PetscFunctionBegin; 1564fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1572b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,k,2); 1584fad6a16SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 1594fad6a16SStefano Zampini PetscFunctionReturn(0); 1604fad6a16SStefano Zampini } 1612b510759SStefano Zampini 162b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 1632b510759SStefano Zampini #undef __FUNCT__ 164b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 165b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 166b8ffe317SStefano Zampini { 167b8ffe317SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 168b8ffe317SStefano Zampini 169b8ffe317SStefano Zampini PetscFunctionBegin; 17085c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = flg; 171b8ffe317SStefano Zampini PetscFunctionReturn(0); 172b8ffe317SStefano Zampini } 173b8ffe317SStefano Zampini 174b8ffe317SStefano Zampini #undef __FUNCT__ 175b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 176b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 1772b510759SStefano Zampini { 1782b510759SStefano Zampini PetscErrorCode ierr; 1792b510759SStefano Zampini 1802b510759SStefano Zampini PetscFunctionBegin; 1812b510759SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 182b8ffe317SStefano Zampini PetscValidLogicalCollectiveBool(pc,flg,2); 183b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1842b510759SStefano Zampini PetscFunctionReturn(0); 1852b510759SStefano Zampini } 1861e6b0712SBarry Smith 1874fad6a16SStefano Zampini #undef __FUNCT__ 1882b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC" 1892b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 1904fad6a16SStefano Zampini { 1914fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1924fad6a16SStefano Zampini 1934fad6a16SStefano Zampini PetscFunctionBegin; 1942b510759SStefano Zampini pcbddc->current_level = level; 1954fad6a16SStefano Zampini PetscFunctionReturn(0); 1964fad6a16SStefano Zampini } 1971e6b0712SBarry Smith 1984fad6a16SStefano Zampini #undef __FUNCT__ 199b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel" 200b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 201b8ffe317SStefano Zampini { 202b8ffe317SStefano Zampini PetscErrorCode ierr; 203b8ffe317SStefano Zampini 204b8ffe317SStefano Zampini PetscFunctionBegin; 205b8ffe317SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 206b8ffe317SStefano Zampini PetscValidLogicalCollectiveInt(pc,level,2); 207b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 208b8ffe317SStefano Zampini PetscFunctionReturn(0); 209b8ffe317SStefano Zampini } 210b8ffe317SStefano Zampini 211b8ffe317SStefano Zampini #undef __FUNCT__ 2122b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC" 2132b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 2142b510759SStefano Zampini { 2152b510759SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2162b510759SStefano Zampini 2172b510759SStefano Zampini PetscFunctionBegin; 2182b510759SStefano Zampini pcbddc->max_levels = levels; 2192b510759SStefano Zampini PetscFunctionReturn(0); 2202b510759SStefano Zampini } 2212b510759SStefano Zampini 2222b510759SStefano Zampini #undef __FUNCT__ 2232b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels" 2244fad6a16SStefano Zampini /*@ 22528509bceSStefano Zampini PCBDDCSetLevels - Sets the maximum number of levels for multilevel 2264fad6a16SStefano Zampini 2274fad6a16SStefano Zampini Logically collective on PC 2284fad6a16SStefano Zampini 2294fad6a16SStefano Zampini Input Parameters: 2304fad6a16SStefano Zampini + pc - the preconditioning context 23128509bceSStefano Zampini - levels - the maximum number of levels (max 9) 2324fad6a16SStefano Zampini 23328509bceSStefano Zampini Default value is 0, i.e. traditional one-level BDDC 2344fad6a16SStefano Zampini 2354fad6a16SStefano Zampini Level: intermediate 2364fad6a16SStefano Zampini 2374fad6a16SStefano Zampini Notes: 2384fad6a16SStefano Zampini 2394fad6a16SStefano Zampini .seealso: PCBDDC 2404fad6a16SStefano Zampini @*/ 2412b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 2424fad6a16SStefano Zampini { 2434fad6a16SStefano Zampini PetscErrorCode ierr; 2444fad6a16SStefano Zampini 2454fad6a16SStefano Zampini PetscFunctionBegin; 2464fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2472b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,levels,2); 248312be037SStefano Zampini if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n"); 2492b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 2504fad6a16SStefano Zampini PetscFunctionReturn(0); 2514fad6a16SStefano Zampini } 2524fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */ 2531e6b0712SBarry Smith 2544fad6a16SStefano Zampini #undef __FUNCT__ 2550bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 2560bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 2570bdf917eSStefano Zampini { 2580bdf917eSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2590bdf917eSStefano Zampini PetscErrorCode ierr; 2600bdf917eSStefano Zampini 2610bdf917eSStefano Zampini PetscFunctionBegin; 2620bdf917eSStefano Zampini ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 2630bdf917eSStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 2640bdf917eSStefano Zampini pcbddc->NullSpace = NullSpace; 2650bdf917eSStefano Zampini PetscFunctionReturn(0); 2660bdf917eSStefano Zampini } 2671e6b0712SBarry Smith 2680bdf917eSStefano Zampini #undef __FUNCT__ 2690bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace" 2700bdf917eSStefano Zampini /*@ 27128509bceSStefano Zampini PCBDDCSetNullSpace - Set nullspace for BDDC operator 2720bdf917eSStefano Zampini 2730bdf917eSStefano Zampini Logically collective on PC and MatNullSpace 2740bdf917eSStefano Zampini 2750bdf917eSStefano Zampini Input Parameters: 2760bdf917eSStefano Zampini + pc - the preconditioning context 27728509bceSStefano Zampini - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 2780bdf917eSStefano Zampini 2790bdf917eSStefano Zampini Level: intermediate 2800bdf917eSStefano Zampini 2810bdf917eSStefano Zampini Notes: 2820bdf917eSStefano Zampini 2830bdf917eSStefano Zampini .seealso: PCBDDC 2840bdf917eSStefano Zampini @*/ 2850bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 2860bdf917eSStefano Zampini { 2870bdf917eSStefano Zampini PetscErrorCode ierr; 2880bdf917eSStefano Zampini 2890bdf917eSStefano Zampini PetscFunctionBegin; 2900bdf917eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 291674ae819SStefano Zampini PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 2922b510759SStefano Zampini PetscCheckSameComm(pc,1,NullSpace,2); 2930bdf917eSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 2940bdf917eSStefano Zampini PetscFunctionReturn(0); 2950bdf917eSStefano Zampini } 2960bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */ 2971e6b0712SBarry Smith 2980bdf917eSStefano Zampini #undef __FUNCT__ 2993b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 3003b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 3013b03a366Sstefano_zampini { 3023b03a366Sstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3033b03a366Sstefano_zampini PetscErrorCode ierr; 3043b03a366Sstefano_zampini 3053b03a366Sstefano_zampini PetscFunctionBegin; 306785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 307785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 3083b03a366Sstefano_zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 30936e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 31036e030ebSStefano Zampini pcbddc->DirichletBoundaries = DirichletBoundaries; 311fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 3123b03a366Sstefano_zampini PetscFunctionReturn(0); 3133b03a366Sstefano_zampini } 3141e6b0712SBarry Smith 3153b03a366Sstefano_zampini #undef __FUNCT__ 3163b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 3173b03a366Sstefano_zampini /*@ 31828509bceSStefano Zampini PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 3193b03a366Sstefano_zampini 320785d1243SStefano Zampini Collective 3213b03a366Sstefano_zampini 3223b03a366Sstefano_zampini Input Parameters: 3233b03a366Sstefano_zampini + pc - the preconditioning context 324785d1243SStefano Zampini - DirichletBoundaries - parallel IS defining the Dirichlet boundaries 3253b03a366Sstefano_zampini 3263b03a366Sstefano_zampini Level: intermediate 3273b03a366Sstefano_zampini 328785d1243SStefano Zampini Notes: Any process can list any global node 3293b03a366Sstefano_zampini 3303b03a366Sstefano_zampini .seealso: PCBDDC 3313b03a366Sstefano_zampini @*/ 3323b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 3333b03a366Sstefano_zampini { 3343b03a366Sstefano_zampini PetscErrorCode ierr; 3353b03a366Sstefano_zampini 3363b03a366Sstefano_zampini PetscFunctionBegin; 3373b03a366Sstefano_zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 338674ae819SStefano Zampini PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 339785d1243SStefano Zampini PetscCheckSameComm(pc,1,DirichletBoundaries,2); 3403b03a366Sstefano_zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 3413b03a366Sstefano_zampini PetscFunctionReturn(0); 3423b03a366Sstefano_zampini } 3433b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */ 3441e6b0712SBarry Smith 3453b03a366Sstefano_zampini #undef __FUNCT__ 34682ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC" 34782ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries) 3480c7d97c5SJed Brown { 3490c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3500c7d97c5SJed Brown PetscErrorCode ierr; 3510c7d97c5SJed Brown 3520c7d97c5SJed Brown PetscFunctionBegin; 353785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 354785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 3550c7d97c5SJed Brown ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 3560c7d97c5SJed Brown ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 357785d1243SStefano Zampini pcbddc->DirichletBoundariesLocal = DirichletBoundaries; 3587d24d155SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 3590c7d97c5SJed Brown PetscFunctionReturn(0); 3600c7d97c5SJed Brown } 3610c7d97c5SJed Brown 3620c7d97c5SJed Brown #undef __FUNCT__ 36382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal" 3649c0446d6SStefano Zampini /*@ 36582ba6b80SStefano Zampini PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering. 3669c0446d6SStefano Zampini 367785d1243SStefano Zampini Collective 36853cdbc3dSStefano Zampini 36953cdbc3dSStefano Zampini Input Parameters: 37053cdbc3dSStefano Zampini + pc - the preconditioning context 37182ba6b80SStefano Zampini - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering) 37253cdbc3dSStefano Zampini 37353cdbc3dSStefano Zampini Level: intermediate 37453cdbc3dSStefano Zampini 3759c0446d6SStefano Zampini Notes: 37653cdbc3dSStefano Zampini 37753cdbc3dSStefano Zampini .seealso: PCBDDC 37853cdbc3dSStefano Zampini @*/ 37982ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries) 3800c7d97c5SJed Brown { 3810c7d97c5SJed Brown PetscErrorCode ierr; 3820c7d97c5SJed Brown 3830c7d97c5SJed Brown PetscFunctionBegin; 3840c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3850c7d97c5SJed Brown PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 38682ba6b80SStefano Zampini PetscCheckSameComm(pc,1,DirichletBoundaries,2); 38782ba6b80SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 3880c7d97c5SJed Brown PetscFunctionReturn(0); 3890c7d97c5SJed Brown } 3900c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 3910c7d97c5SJed Brown 3920c7d97c5SJed Brown #undef __FUNCT__ 3930c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 39453cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 3950c7d97c5SJed Brown { 3960c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 39753cdbc3dSStefano Zampini PetscErrorCode ierr; 3980c7d97c5SJed Brown 3990c7d97c5SJed Brown PetscFunctionBegin; 400785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 401785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 40253cdbc3dSStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 40336e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 40436e030ebSStefano Zampini pcbddc->NeumannBoundaries = NeumannBoundaries; 405fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4060c7d97c5SJed Brown PetscFunctionReturn(0); 4070c7d97c5SJed Brown } 4081e6b0712SBarry Smith 4090c7d97c5SJed Brown #undef __FUNCT__ 4100c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 41157527edcSJed Brown /*@ 41228509bceSStefano Zampini PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 41357527edcSJed Brown 414785d1243SStefano Zampini Collective 41557527edcSJed Brown 41657527edcSJed Brown Input Parameters: 41757527edcSJed Brown + pc - the preconditioning context 418785d1243SStefano Zampini - NeumannBoundaries - parallel IS defining the Neumann boundaries 41957527edcSJed Brown 42057527edcSJed Brown Level: intermediate 42157527edcSJed Brown 422785d1243SStefano Zampini Notes: Any process can list any global node 42357527edcSJed Brown 42457527edcSJed Brown .seealso: PCBDDC 42557527edcSJed Brown @*/ 42653cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 4270c7d97c5SJed Brown { 4280c7d97c5SJed Brown PetscErrorCode ierr; 4290c7d97c5SJed Brown 4300c7d97c5SJed Brown PetscFunctionBegin; 4310c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 432674ae819SStefano Zampini PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 433785d1243SStefano Zampini PetscCheckSameComm(pc,1,NeumannBoundaries,2); 43453cdbc3dSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 43553cdbc3dSStefano Zampini PetscFunctionReturn(0); 43653cdbc3dSStefano Zampini } 43753cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 4381e6b0712SBarry Smith 43953cdbc3dSStefano Zampini #undef __FUNCT__ 44082ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC" 44182ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries) 4420c7d97c5SJed Brown { 4430c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 4440c7d97c5SJed Brown PetscErrorCode ierr; 4450c7d97c5SJed Brown 4460c7d97c5SJed Brown PetscFunctionBegin; 447785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 448785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 4490c7d97c5SJed Brown ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 4500c7d97c5SJed Brown ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 451785d1243SStefano Zampini pcbddc->NeumannBoundariesLocal = NeumannBoundaries; 4527d24d155SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4530c7d97c5SJed Brown PetscFunctionReturn(0); 4540c7d97c5SJed Brown } 4550c7d97c5SJed Brown 4560c7d97c5SJed Brown #undef __FUNCT__ 45782ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal" 4580c7d97c5SJed Brown /*@ 45982ba6b80SStefano Zampini PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering. 4600c7d97c5SJed Brown 461785d1243SStefano Zampini Collective 4620c7d97c5SJed Brown 4630c7d97c5SJed Brown Input Parameters: 4640c7d97c5SJed Brown + pc - the preconditioning context 46582ba6b80SStefano Zampini - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering) 4660c7d97c5SJed Brown 4670c7d97c5SJed Brown Level: intermediate 4680c7d97c5SJed Brown 4690c7d97c5SJed Brown Notes: 4700c7d97c5SJed Brown 4710c7d97c5SJed Brown .seealso: PCBDDC 4720c7d97c5SJed Brown @*/ 47382ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries) 4740c7d97c5SJed Brown { 4750c7d97c5SJed Brown PetscErrorCode ierr; 4760c7d97c5SJed Brown 4770c7d97c5SJed Brown PetscFunctionBegin; 4780c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4790c7d97c5SJed Brown PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 48082ba6b80SStefano Zampini PetscCheckSameComm(pc,1,NeumannBoundaries,2); 48182ba6b80SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 48253cdbc3dSStefano Zampini PetscFunctionReturn(0); 48353cdbc3dSStefano Zampini } 48453cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 48553cdbc3dSStefano Zampini 48653cdbc3dSStefano Zampini #undef __FUNCT__ 487da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 488da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 489da1bb401SStefano Zampini { 490da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 491da1bb401SStefano Zampini 492da1bb401SStefano Zampini PetscFunctionBegin; 493da1bb401SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundaries; 494da1bb401SStefano Zampini PetscFunctionReturn(0); 495da1bb401SStefano Zampini } 4961e6b0712SBarry Smith 497da1bb401SStefano Zampini #undef __FUNCT__ 498da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 499da1bb401SStefano Zampini /*@ 500785d1243SStefano Zampini PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries 501da1bb401SStefano Zampini 502785d1243SStefano Zampini Collective 503785d1243SStefano Zampini 504785d1243SStefano Zampini Input Parameters: 505785d1243SStefano Zampini . pc - the preconditioning context 506785d1243SStefano Zampini 507785d1243SStefano Zampini Output Parameters: 508785d1243SStefano Zampini . DirichletBoundaries - index set defining the Dirichlet boundaries 509785d1243SStefano Zampini 510785d1243SStefano Zampini Level: intermediate 511785d1243SStefano Zampini 512785d1243SStefano Zampini Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries 513785d1243SStefano Zampini 514785d1243SStefano Zampini .seealso: PCBDDC 515785d1243SStefano Zampini @*/ 516785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 517785d1243SStefano Zampini { 518785d1243SStefano Zampini PetscErrorCode ierr; 519785d1243SStefano Zampini 520785d1243SStefano Zampini PetscFunctionBegin; 521785d1243SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 522785d1243SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 523785d1243SStefano Zampini PetscFunctionReturn(0); 524785d1243SStefano Zampini } 525785d1243SStefano Zampini /* -------------------------------------------------------------------------- */ 526785d1243SStefano Zampini 527785d1243SStefano Zampini #undef __FUNCT__ 528785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC" 529785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries) 530785d1243SStefano Zampini { 531785d1243SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 532785d1243SStefano Zampini 533785d1243SStefano Zampini PetscFunctionBegin; 534785d1243SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundariesLocal; 535785d1243SStefano Zampini PetscFunctionReturn(0); 536785d1243SStefano Zampini } 537785d1243SStefano Zampini 538785d1243SStefano Zampini #undef __FUNCT__ 53982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal" 540da1bb401SStefano Zampini /*@ 54182ba6b80SStefano Zampini PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering) 542da1bb401SStefano Zampini 543785d1243SStefano Zampini Collective 544da1bb401SStefano Zampini 545da1bb401SStefano Zampini Input Parameters: 54628509bceSStefano Zampini . pc - the preconditioning context 547da1bb401SStefano Zampini 548da1bb401SStefano Zampini Output Parameters: 54928509bceSStefano Zampini . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 550da1bb401SStefano Zampini 551da1bb401SStefano Zampini Level: intermediate 552da1bb401SStefano Zampini 553da1bb401SStefano Zampini Notes: 554da1bb401SStefano Zampini 555da1bb401SStefano Zampini .seealso: PCBDDC 556da1bb401SStefano Zampini @*/ 55782ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries) 558da1bb401SStefano Zampini { 559da1bb401SStefano Zampini PetscErrorCode ierr; 560da1bb401SStefano Zampini 561da1bb401SStefano Zampini PetscFunctionBegin; 562da1bb401SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 56382ba6b80SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 564da1bb401SStefano Zampini PetscFunctionReturn(0); 565da1bb401SStefano Zampini } 566da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 5671e6b0712SBarry Smith 568da1bb401SStefano Zampini #undef __FUNCT__ 56953cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 57053cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 57153cdbc3dSStefano Zampini { 57253cdbc3dSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 57353cdbc3dSStefano Zampini 57453cdbc3dSStefano Zampini PetscFunctionBegin; 57553cdbc3dSStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundaries; 57653cdbc3dSStefano Zampini PetscFunctionReturn(0); 57753cdbc3dSStefano Zampini } 5781e6b0712SBarry Smith 57953cdbc3dSStefano Zampini #undef __FUNCT__ 58053cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 58153cdbc3dSStefano Zampini /*@ 582785d1243SStefano Zampini PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries 58353cdbc3dSStefano Zampini 584785d1243SStefano Zampini Collective 585785d1243SStefano Zampini 586785d1243SStefano Zampini Input Parameters: 587785d1243SStefano Zampini . pc - the preconditioning context 588785d1243SStefano Zampini 589785d1243SStefano Zampini Output Parameters: 590785d1243SStefano Zampini . NeumannBoundaries - index set defining the Neumann boundaries 591785d1243SStefano Zampini 592785d1243SStefano Zampini Level: intermediate 593785d1243SStefano Zampini 594785d1243SStefano Zampini Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries 595785d1243SStefano Zampini 596785d1243SStefano Zampini .seealso: PCBDDC 597785d1243SStefano Zampini @*/ 598785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 599785d1243SStefano Zampini { 600785d1243SStefano Zampini PetscErrorCode ierr; 601785d1243SStefano Zampini 602785d1243SStefano Zampini PetscFunctionBegin; 603785d1243SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 604785d1243SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 605785d1243SStefano Zampini PetscFunctionReturn(0); 606785d1243SStefano Zampini } 607785d1243SStefano Zampini /* -------------------------------------------------------------------------- */ 608785d1243SStefano Zampini 609785d1243SStefano Zampini #undef __FUNCT__ 610785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC" 611785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries) 612785d1243SStefano Zampini { 613785d1243SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 614785d1243SStefano Zampini 615785d1243SStefano Zampini PetscFunctionBegin; 616785d1243SStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundariesLocal; 617785d1243SStefano Zampini PetscFunctionReturn(0); 618785d1243SStefano Zampini } 619785d1243SStefano Zampini 620785d1243SStefano Zampini #undef __FUNCT__ 62182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal" 62253cdbc3dSStefano Zampini /*@ 62382ba6b80SStefano Zampini PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering) 62453cdbc3dSStefano Zampini 625785d1243SStefano Zampini Collective 62653cdbc3dSStefano Zampini 62753cdbc3dSStefano Zampini Input Parameters: 62828509bceSStefano Zampini . pc - the preconditioning context 62953cdbc3dSStefano Zampini 63053cdbc3dSStefano Zampini Output Parameters: 63128509bceSStefano Zampini . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 63253cdbc3dSStefano Zampini 63353cdbc3dSStefano Zampini Level: intermediate 63453cdbc3dSStefano Zampini 63553cdbc3dSStefano Zampini Notes: 63653cdbc3dSStefano Zampini 63753cdbc3dSStefano Zampini .seealso: PCBDDC 63853cdbc3dSStefano Zampini @*/ 63982ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries) 64053cdbc3dSStefano Zampini { 64153cdbc3dSStefano Zampini PetscErrorCode ierr; 64253cdbc3dSStefano Zampini 64353cdbc3dSStefano Zampini PetscFunctionBegin; 64453cdbc3dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 64582ba6b80SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 6460c7d97c5SJed Brown PetscFunctionReturn(0); 6470c7d97c5SJed Brown } 64836e030ebSStefano Zampini /* -------------------------------------------------------------------------- */ 6491e6b0712SBarry Smith 65036e030ebSStefano Zampini #undef __FUNCT__ 651da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 6521a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 65336e030ebSStefano Zampini { 65436e030ebSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 655da1bb401SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 656da1bb401SStefano Zampini PetscErrorCode ierr; 65736e030ebSStefano Zampini 65836e030ebSStefano Zampini PetscFunctionBegin; 659674ae819SStefano Zampini /* free old CSR */ 660674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 661fb223d50SStefano Zampini /* TODO: PCBDDCGraphSetAdjacency */ 662674ae819SStefano Zampini /* get CSR into graph structure */ 663da1bb401SStefano Zampini if (copymode == PETSC_COPY_VALUES) { 664785e854fSJed Brown ierr = PetscMalloc1((nvtxs+1),&mat_graph->xadj);CHKERRQ(ierr); 665785e854fSJed Brown ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr); 666674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 667674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 668da1bb401SStefano Zampini } else if (copymode == PETSC_OWN_POINTER) { 6691a83f524SJed Brown mat_graph->xadj = (PetscInt*)xadj; 6701a83f524SJed Brown mat_graph->adjncy = (PetscInt*)adjncy; 671674ae819SStefano Zampini } 672575ad6abSStefano Zampini mat_graph->nvtxs_csr = nvtxs; 67336e030ebSStefano Zampini PetscFunctionReturn(0); 67436e030ebSStefano Zampini } 6751e6b0712SBarry Smith 67636e030ebSStefano Zampini #undef __FUNCT__ 677da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 67836e030ebSStefano Zampini /*@ 67928509bceSStefano Zampini PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix 68036e030ebSStefano Zampini 68136e030ebSStefano Zampini Not collective 68236e030ebSStefano Zampini 68336e030ebSStefano Zampini Input Parameters: 68436e030ebSStefano Zampini + pc - the preconditioning context 68528509bceSStefano Zampini . nvtxs - number of local vertices of the graph (i.e., the local size of your problem) 68628509bceSStefano Zampini . xadj, adjncy - the CSR graph 68728509bceSStefano Zampini - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 68836e030ebSStefano Zampini 68936e030ebSStefano Zampini Level: intermediate 69036e030ebSStefano Zampini 69136e030ebSStefano Zampini Notes: 69236e030ebSStefano Zampini 69328509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode 69436e030ebSStefano Zampini @*/ 6951a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 69636e030ebSStefano Zampini { 697575ad6abSStefano Zampini void (*f)(void) = 0; 69836e030ebSStefano Zampini PetscErrorCode ierr; 69936e030ebSStefano Zampini 70036e030ebSStefano Zampini PetscFunctionBegin; 70136e030ebSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 702674ae819SStefano Zampini PetscValidIntPointer(xadj,3); 703674ae819SStefano Zampini PetscValidIntPointer(xadj,4); 704674ae819SStefano Zampini if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 705674ae819SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 706da1bb401SStefano Zampini } 70736e030ebSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 708575ad6abSStefano Zampini /* free arrays if PCBDDC is not the PC type */ 709575ad6abSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 710575ad6abSStefano Zampini if (!f && copymode == PETSC_OWN_POINTER) { 711575ad6abSStefano Zampini ierr = PetscFree(xadj);CHKERRQ(ierr); 712575ad6abSStefano Zampini ierr = PetscFree(adjncy);CHKERRQ(ierr); 71336e030ebSStefano Zampini } 71436e030ebSStefano Zampini PetscFunctionReturn(0); 71536e030ebSStefano Zampini } 7169c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */ 7171e6b0712SBarry Smith 7189c0446d6SStefano Zampini #undef __FUNCT__ 71963602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC" 72063602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 72163602bcaSStefano Zampini { 72263602bcaSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 72363602bcaSStefano Zampini PetscInt i; 72463602bcaSStefano Zampini PetscErrorCode ierr; 72563602bcaSStefano Zampini 72663602bcaSStefano Zampini PetscFunctionBegin; 72763602bcaSStefano Zampini /* Destroy ISes if they were already set */ 72863602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 72963602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 73063602bcaSStefano Zampini } 73163602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 73263602bcaSStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 73363602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 73463602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 73563602bcaSStefano Zampini } 73663602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 73763602bcaSStefano Zampini pcbddc->n_ISForDofs = 0; 73863602bcaSStefano Zampini /* allocate space then set */ 73963602bcaSStefano Zampini ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr); 74063602bcaSStefano Zampini for (i=0;i<n_is;i++) { 74163602bcaSStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 74263602bcaSStefano Zampini pcbddc->ISForDofsLocal[i]=ISForDofs[i]; 74363602bcaSStefano Zampini } 74463602bcaSStefano Zampini pcbddc->n_ISForDofsLocal=n_is; 74563602bcaSStefano Zampini if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 7461a8a16d1SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 74763602bcaSStefano Zampini PetscFunctionReturn(0); 74863602bcaSStefano Zampini } 74963602bcaSStefano Zampini 75063602bcaSStefano Zampini #undef __FUNCT__ 75163602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal" 75263602bcaSStefano Zampini /*@ 75363602bcaSStefano Zampini PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix 75463602bcaSStefano Zampini 75563602bcaSStefano Zampini Collective 75663602bcaSStefano Zampini 75763602bcaSStefano Zampini Input Parameters: 75863602bcaSStefano Zampini + pc - the preconditioning context 75963602bcaSStefano Zampini - n_is - number of index sets defining the fields 76063602bcaSStefano Zampini . ISForDofs - array of IS describing the fields in local ordering 76163602bcaSStefano Zampini 76263602bcaSStefano Zampini Level: intermediate 76363602bcaSStefano Zampini 76463602bcaSStefano 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. 76563602bcaSStefano Zampini 76663602bcaSStefano Zampini .seealso: PCBDDC 76763602bcaSStefano Zampini @*/ 76863602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[]) 76963602bcaSStefano Zampini { 77063602bcaSStefano Zampini PetscInt i; 77163602bcaSStefano Zampini PetscErrorCode ierr; 77263602bcaSStefano Zampini 77363602bcaSStefano Zampini PetscFunctionBegin; 77463602bcaSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 77563602bcaSStefano Zampini PetscValidLogicalCollectiveInt(pc,n_is,2); 77663602bcaSStefano Zampini for (i=0;i<n_is;i++) { 77763602bcaSStefano Zampini PetscCheckSameComm(pc,1,ISForDofs[i],3); 77863602bcaSStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 77963602bcaSStefano Zampini } 78063602bcaSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 78163602bcaSStefano Zampini PetscFunctionReturn(0); 78263602bcaSStefano Zampini } 78363602bcaSStefano Zampini /* -------------------------------------------------------------------------- */ 78463602bcaSStefano Zampini 78563602bcaSStefano Zampini #undef __FUNCT__ 7869c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 7879c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 7889c0446d6SStefano Zampini { 7899c0446d6SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 7909c0446d6SStefano Zampini PetscInt i; 7919c0446d6SStefano Zampini PetscErrorCode ierr; 7929c0446d6SStefano Zampini 7939c0446d6SStefano Zampini PetscFunctionBegin; 794da1bb401SStefano Zampini /* Destroy ISes if they were already set */ 7959c0446d6SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 7969c0446d6SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 7979c0446d6SStefano Zampini } 798d11ae9bbSstefano_zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 79963602bcaSStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 80063602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 80163602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 80263602bcaSStefano Zampini } 80363602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 80463602bcaSStefano Zampini pcbddc->n_ISForDofsLocal = 0; 805da1bb401SStefano Zampini /* allocate space then set */ 806785e854fSJed Brown ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr); 8079c0446d6SStefano Zampini for (i=0;i<n_is;i++) { 808da1bb401SStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 809da1bb401SStefano Zampini pcbddc->ISForDofs[i]=ISForDofs[i]; 8109c0446d6SStefano Zampini } 8119c0446d6SStefano Zampini pcbddc->n_ISForDofs=n_is; 81263602bcaSStefano Zampini if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 8131a8a16d1SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 8149c0446d6SStefano Zampini PetscFunctionReturn(0); 8159c0446d6SStefano Zampini } 8161e6b0712SBarry Smith 8179c0446d6SStefano Zampini #undef __FUNCT__ 8189c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting" 8199c0446d6SStefano Zampini /*@ 82063602bcaSStefano Zampini PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix 8219c0446d6SStefano Zampini 82263602bcaSStefano Zampini Collective 8239c0446d6SStefano Zampini 8249c0446d6SStefano Zampini Input Parameters: 8259c0446d6SStefano Zampini + pc - the preconditioning context 82628509bceSStefano Zampini - n_is - number of index sets defining the fields 82763602bcaSStefano Zampini . ISForDofs - array of IS describing the fields in global ordering 8289c0446d6SStefano Zampini 8299c0446d6SStefano Zampini Level: intermediate 8309c0446d6SStefano Zampini 83163602bcaSStefano Zampini Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field. 8329c0446d6SStefano Zampini 8339c0446d6SStefano Zampini .seealso: PCBDDC 8349c0446d6SStefano Zampini @*/ 8359c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 8369c0446d6SStefano Zampini { 8372b510759SStefano Zampini PetscInt i; 8389c0446d6SStefano Zampini PetscErrorCode ierr; 8399c0446d6SStefano Zampini 8409c0446d6SStefano Zampini PetscFunctionBegin; 8419c0446d6SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 84263602bcaSStefano Zampini PetscValidLogicalCollectiveInt(pc,n_is,2); 8432b510759SStefano Zampini for (i=0;i<n_is;i++) { 84463602bcaSStefano Zampini PetscCheckSameComm(pc,1,ISForDofs[i],3); 84563602bcaSStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 8462b510759SStefano Zampini } 8479c0446d6SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 8489c0446d6SStefano Zampini PetscFunctionReturn(0); 8499c0446d6SStefano Zampini } 850da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 851534831adSStefano Zampini #undef __FUNCT__ 852534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC" 853534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 854534831adSStefano Zampini /* 855534831adSStefano Zampini PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 856534831adSStefano Zampini guess if a transformation of basis approach has been selected. 8579c0446d6SStefano Zampini 858534831adSStefano Zampini Input Parameter: 859534831adSStefano Zampini + pc - the preconditioner contex 860534831adSStefano Zampini 861534831adSStefano Zampini Application Interface Routine: PCPreSolve() 862534831adSStefano Zampini 863534831adSStefano Zampini Notes: 864534831adSStefano Zampini The interface routine PCPreSolve() is not usually called directly by 865534831adSStefano Zampini the user, but instead is called by KSPSolve(). 866534831adSStefano Zampini */ 867534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 868534831adSStefano Zampini { 869534831adSStefano Zampini PetscErrorCode ierr; 870534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 871534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 872534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 873534831adSStefano Zampini Mat temp_mat; 8743972b0daSStefano Zampini IS dirIS; 8753972b0daSStefano Zampini Vec used_vec; 87682ba6b80SStefano Zampini PetscBool guess_nonzero; 877534831adSStefano Zampini 878534831adSStefano Zampini PetscFunctionBegin; 87985c4d303SStefano Zampini /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 88085c4d303SStefano Zampini if (ksp) { 88185c4d303SStefano Zampini PetscBool iscg; 88285c4d303SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 88385c4d303SStefano Zampini if (!iscg) { 88485c4d303SStefano Zampini ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 88585c4d303SStefano Zampini } 88685c4d303SStefano Zampini } 88785c4d303SStefano Zampini /* Creates parallel work vectors used in presolve */ 88862a6ff1dSStefano Zampini if (!pcbddc->original_rhs) { 88962a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 89062a6ff1dSStefano Zampini } 89162a6ff1dSStefano Zampini if (!pcbddc->temp_solution) { 89262a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 89362a6ff1dSStefano Zampini } 8943972b0daSStefano Zampini if (x) { 8953972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 8963972b0daSStefano Zampini used_vec = x; 8973972b0daSStefano Zampini } else { 8983972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 8993972b0daSStefano Zampini used_vec = pcbddc->temp_solution; 9003972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 9013972b0daSStefano Zampini } 9023972b0daSStefano Zampini /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 9033972b0daSStefano Zampini if (ksp) { 9043972b0daSStefano Zampini ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 9053972b0daSStefano Zampini if (!guess_nonzero) { 9063972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 9073972b0daSStefano Zampini } 9083972b0daSStefano Zampini } 9093308cffdSStefano Zampini 9103972b0daSStefano Zampini /* store the original rhs */ 9113972b0daSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 9123972b0daSStefano Zampini 9133972b0daSStefano Zampini /* Take into account zeroed rows -> change rhs and store solution removed */ 914785d1243SStefano Zampini /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */ 91582ba6b80SStefano Zampini ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr); 91682ba6b80SStefano Zampini if (rhs && dirIS) { 917785d1243SStefano Zampini PetscInt dirsize,i,*is_indices; 918785d1243SStefano Zampini PetscScalar *array_x,*array_diagonal; 919785d1243SStefano Zampini 9203972b0daSStefano Zampini ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 9213972b0daSStefano Zampini ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 9223972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9233972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9243972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9253972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 92682ba6b80SStefano Zampini ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr); 9273972b0daSStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 9283972b0daSStefano Zampini ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 9293972b0daSStefano Zampini ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 9302fa5cd67SKarl Rupp for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 9313972b0daSStefano Zampini ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 9323972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 9333972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 9343972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9353972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 936b76ba322SStefano Zampini 9373972b0daSStefano Zampini /* remove the computed solution from the rhs */ 9383972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 9393972b0daSStefano Zampini ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 9403972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 9413308cffdSStefano Zampini } 942b76ba322SStefano Zampini 943b76ba322SStefano Zampini /* store partially computed solution and set initial guess */ 9443972b0daSStefano Zampini if (x) { 9453972b0daSStefano Zampini ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 9463972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 94785c4d303SStefano Zampini if (pcbddc->use_exact_dirichlet_trick) { 948b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 949b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 950b76ba322SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 951b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 952b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 953b76ba322SStefano Zampini if (ksp) { 954b76ba322SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 955b76ba322SStefano Zampini } 956b76ba322SStefano Zampini } 9573972b0daSStefano Zampini } 958b76ba322SStefano Zampini 95992e3dcfbSStefano Zampini /* prepare MatMult and rhs for solver */ 960674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 961b76ba322SStefano Zampini /* swap pointers for local matrices */ 962b76ba322SStefano Zampini temp_mat = matis->A; 963b76ba322SStefano Zampini matis->A = pcbddc->local_mat; 964b76ba322SStefano Zampini pcbddc->local_mat = temp_mat; 96592e3dcfbSStefano Zampini if (rhs) { 966b76ba322SStefano Zampini /* Get local rhs and apply transformation of basis */ 967b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 968b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 969b76ba322SStefano Zampini /* from original basis to modified basis */ 970b76ba322SStefano Zampini ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 971b76ba322SStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 972b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 973b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 974674ae819SStefano Zampini } 97592e3dcfbSStefano Zampini } 97692e3dcfbSStefano Zampini 97792e3dcfbSStefano Zampini /* remove nullspace if present */ 9780bdf917eSStefano Zampini if (ksp && pcbddc->NullSpace) { 979d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 980d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 981b76ba322SStefano Zampini } 9820bdf917eSStefano Zampini ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 983534831adSStefano Zampini PetscFunctionReturn(0); 984534831adSStefano Zampini } 985534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 986534831adSStefano Zampini #undef __FUNCT__ 987534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC" 988534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 989534831adSStefano Zampini /* 990534831adSStefano Zampini PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 991534831adSStefano Zampini approach has been selected. Also, restores rhs to its original state. 992534831adSStefano Zampini 993534831adSStefano Zampini Input Parameter: 994534831adSStefano Zampini + pc - the preconditioner contex 995534831adSStefano Zampini 996534831adSStefano Zampini Application Interface Routine: PCPostSolve() 997534831adSStefano Zampini 998534831adSStefano Zampini Notes: 999534831adSStefano Zampini The interface routine PCPostSolve() is not usually called directly by 1000534831adSStefano Zampini the user, but instead is called by KSPSolve(). 1001534831adSStefano Zampini */ 1002534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1003534831adSStefano Zampini { 1004534831adSStefano Zampini PetscErrorCode ierr; 1005534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1006534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 1007534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 1008534831adSStefano Zampini Mat temp_mat; 1009534831adSStefano Zampini 1010534831adSStefano Zampini PetscFunctionBegin; 1011674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 1012534831adSStefano Zampini /* swap pointers for local matrices */ 1013534831adSStefano Zampini temp_mat = matis->A; 1014534831adSStefano Zampini matis->A = pcbddc->local_mat; 1015534831adSStefano Zampini pcbddc->local_mat = temp_mat; 10163425bc38SStefano Zampini } 10173308cffdSStefano Zampini if (pcbddc->use_change_of_basis && x) { 1018534831adSStefano Zampini /* Get Local boundary and apply transformation of basis to solution vector */ 1019534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1020534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1021534831adSStefano Zampini /* from modified basis to original basis */ 1022534831adSStefano Zampini ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 1023534831adSStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 1024534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1025534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1026534831adSStefano Zampini } 10273972b0daSStefano Zampini /* add solution removed in presolve */ 10283425bc38SStefano Zampini if (x) { 10293425bc38SStefano Zampini ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 10303425bc38SStefano Zampini } 1031fb223d50SStefano Zampini /* restore rhs to its original state */ 1032fb223d50SStefano Zampini if (rhs) { 1033fb223d50SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 1034fb223d50SStefano Zampini } 1035534831adSStefano Zampini PetscFunctionReturn(0); 1036534831adSStefano Zampini } 1037534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 103853cdbc3dSStefano Zampini #undef __FUNCT__ 103953cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC" 10400c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 10410c7d97c5SJed Brown /* 10420c7d97c5SJed Brown PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 10430c7d97c5SJed Brown by setting data structures and options. 10440c7d97c5SJed Brown 10450c7d97c5SJed Brown Input Parameter: 104653cdbc3dSStefano Zampini + pc - the preconditioner context 10470c7d97c5SJed Brown 10480c7d97c5SJed Brown Application Interface Routine: PCSetUp() 10490c7d97c5SJed Brown 10500c7d97c5SJed Brown Notes: 10510c7d97c5SJed Brown The interface routine PCSetUp() is not usually called directly by 10520c7d97c5SJed Brown the user, but instead is called by PCApply() if necessary. 10530c7d97c5SJed Brown */ 105453cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc) 10550c7d97c5SJed Brown { 10560c7d97c5SJed Brown PetscErrorCode ierr; 10570c7d97c5SJed Brown PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 1058f4ddd8eeSStefano Zampini MatNullSpace nearnullspace; 1059674ae819SStefano Zampini PetscBool computeis,computetopography,computesolvers; 1060165b64e2SStefano Zampini PetscBool new_nearnullspace_provided; 10610c7d97c5SJed Brown 10620c7d97c5SJed Brown PetscFunctionBegin; 1063f4ddd8eeSStefano Zampini /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 10643b03a366Sstefano_zampini /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1065fcd409f5SStefano Zampini Also, BDDC directly build the Dirichlet problem */ 1066f4ddd8eeSStefano Zampini 1067f4ddd8eeSStefano Zampini /* split work */ 1068674ae819SStefano Zampini if (pc->setupcalled) { 1069674ae819SStefano Zampini computeis = PETSC_FALSE; 1070d1e9a80fSBarry Smith if (pc->flag == SAME_NONZERO_PATTERN) { 1071674ae819SStefano Zampini computetopography = PETSC_FALSE; 1072674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1073674ae819SStefano Zampini } else { /* DIFFERENT_NONZERO_PATTERN */ 1074674ae819SStefano Zampini computetopography = PETSC_TRUE; 1075674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1076674ae819SStefano Zampini } 1077674ae819SStefano Zampini } else { 1078674ae819SStefano Zampini computeis = PETSC_TRUE; 1079674ae819SStefano Zampini computetopography = PETSC_TRUE; 1080674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1081674ae819SStefano Zampini } 1082fb180af4SStefano Zampini if (pcbddc->recompute_topography) { 1083fb180af4SStefano Zampini computetopography = PETSC_TRUE; 1084fb180af4SStefano Zampini } 1085f4ddd8eeSStefano Zampini 1086f4ddd8eeSStefano Zampini /* Get stdout for dbg */ 108770cf5478SStefano Zampini if (pcbddc->dbg_flag) { 108870cf5478SStefano Zampini if (!pcbddc->dbg_viewer) { 108958a03d70SStefano Zampini pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)); 1090f4ddd8eeSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 109170cf5478SStefano Zampini } 109258a03d70SStefano Zampini ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1093f4ddd8eeSStefano Zampini } 1094f4ddd8eeSStefano Zampini 1095fcd409f5SStefano Zampini /* Set up all the "iterative substructuring" common block without computing solvers */ 1096674ae819SStefano Zampini if (computeis) { 1097fcd409f5SStefano Zampini /* HACK INTO PCIS */ 1098fcd409f5SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 1099fcd409f5SStefano Zampini pcis->computesolvers = PETSC_FALSE; 1100674ae819SStefano Zampini ierr = PCISSetUp(pc);CHKERRQ(ierr); 110139e2fb2aSStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr); 1102674ae819SStefano Zampini } 1103f4ddd8eeSStefano Zampini 1104f4ddd8eeSStefano Zampini /* Analyze interface */ 1105674ae819SStefano Zampini if (computetopography) { 1106674ae819SStefano Zampini ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 1107fb8d54d4SStefano Zampini } 1108fb8d54d4SStefano Zampini 1109f4ddd8eeSStefano Zampini /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1110fb8d54d4SStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 1111f4ddd8eeSStefano Zampini ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr); 1112f4ddd8eeSStefano Zampini if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1113f4ddd8eeSStefano Zampini if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1114f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1115f4ddd8eeSStefano Zampini } else { 1116f4ddd8eeSStefano Zampini /* determine if the two nullspaces are different (should be lightweight) */ 1117f4ddd8eeSStefano Zampini if (nearnullspace != pcbddc->onearnullspace) { 1118f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1119165b64e2SStefano Zampini } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1120f4ddd8eeSStefano Zampini PetscInt i; 1121165b64e2SStefano Zampini const Vec *nearnullvecs; 1122165b64e2SStefano Zampini PetscObjectState state; 1123165b64e2SStefano Zampini PetscInt nnsp_size; 1124165b64e2SStefano Zampini ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 1125f4ddd8eeSStefano Zampini for (i=0;i<nnsp_size;i++) { 1126f4ddd8eeSStefano Zampini ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr); 1127165b64e2SStefano Zampini if (pcbddc->onearnullvecs_state[i] != state) { 1128f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1129f4ddd8eeSStefano Zampini break; 1130f4ddd8eeSStefano Zampini } 1131f4ddd8eeSStefano Zampini } 1132f4ddd8eeSStefano Zampini } 1133f4ddd8eeSStefano Zampini } 1134f4ddd8eeSStefano Zampini } else { 1135f4ddd8eeSStefano Zampini if (!nearnullspace) { /* both nearnullspaces are null */ 1136f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 1137f4ddd8eeSStefano Zampini } else { /* nearnullspace attached later */ 1138f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1139f4ddd8eeSStefano Zampini } 1140f4ddd8eeSStefano Zampini } 1141f4ddd8eeSStefano Zampini 1142f4ddd8eeSStefano Zampini /* Setup constraints and related work vectors */ 1143727cdba6SStefano Zampini /* reset primal space flags */ 1144f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1145727cdba6SStefano Zampini pcbddc->new_primal_space_local = PETSC_FALSE; 1146fb8d54d4SStefano Zampini if (computetopography || new_nearnullspace_provided) { 1147727cdba6SStefano Zampini /* It also sets the primal space flags */ 1148674ae819SStefano Zampini ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 1149e7b262bdSStefano Zampini /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1150f4ddd8eeSStefano Zampini ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr); 1151674ae819SStefano Zampini } 1152fb8d54d4SStefano Zampini 1153f4ddd8eeSStefano Zampini if (computesolvers || pcbddc->new_primal_space) { 1154674ae819SStefano Zampini /* reset data */ 1155674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 115699cc7994SStefano Zampini /* Create coarse and local stuffs */ 115799cc7994SStefano Zampini ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 1158674ae819SStefano Zampini ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 11590c7d97c5SJed Brown } 116070cf5478SStefano Zampini 116158a03d70SStefano Zampini if (pcbddc->dbg_flag) { 116258a03d70SStefano Zampini ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 11632b510759SStefano Zampini } 11640c7d97c5SJed Brown PetscFunctionReturn(0); 11650c7d97c5SJed Brown } 11660c7d97c5SJed Brown 11670c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 11680c7d97c5SJed Brown /* 116950efa1b5SStefano Zampini PCApply_BDDC - Applies the BDDC operator to a vector. 11700c7d97c5SJed Brown 11710c7d97c5SJed Brown Input Parameters: 11720c7d97c5SJed Brown . pc - the preconditioner context 11730c7d97c5SJed Brown . r - input vector (global) 11740c7d97c5SJed Brown 11750c7d97c5SJed Brown Output Parameter: 11760c7d97c5SJed Brown . z - output vector (global) 11770c7d97c5SJed Brown 11780c7d97c5SJed Brown Application Interface Routine: PCApply() 11790c7d97c5SJed Brown */ 11800c7d97c5SJed Brown #undef __FUNCT__ 11810c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC" 118253cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 11830c7d97c5SJed Brown { 11840c7d97c5SJed Brown PC_IS *pcis = (PC_IS*)(pc->data); 11850c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 11860c7d97c5SJed Brown PetscErrorCode ierr; 11873b03a366Sstefano_zampini const PetscScalar one = 1.0; 11883b03a366Sstefano_zampini const PetscScalar m_one = -1.0; 11892617d88aSStefano Zampini const PetscScalar zero = 0.0; 11900c7d97c5SJed Brown 11910c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN 11920c7d97c5SJed Brown NN interface preconditioner changed to BDDC 11938eeda7d8SStefano Zampini Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 11940c7d97c5SJed Brown 11950c7d97c5SJed Brown PetscFunctionBegin; 119685c4d303SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 11970c7d97c5SJed Brown /* First Dirichlet solve */ 11980c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 11990c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 120053cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 12010c7d97c5SJed Brown /* 12020c7d97c5SJed Brown Assembling right hand side for BDDC operator 1203674ae819SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 1204674ae819SStefano Zampini - pcis->vec1_B the interface part of the global vector z 12050c7d97c5SJed Brown */ 12060c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 12070c7d97c5SJed Brown ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 12088eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 12090c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 12100c7d97c5SJed Brown ierr = VecCopy(r,z);CHKERRQ(ierr); 12110c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12120c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1213674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1214b76ba322SStefano Zampini } else { 12150bdf917eSStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 1216b76ba322SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 1217674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1218b76ba322SStefano Zampini } 1219b76ba322SStefano Zampini 12202617d88aSStefano Zampini /* Apply interface preconditioner 12212617d88aSStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1222dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr); 12232617d88aSStefano Zampini 1224674ae819SStefano Zampini /* Apply transpose of partition of unity operator */ 1225674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 12260c7d97c5SJed Brown 12273b03a366Sstefano_zampini /* Second Dirichlet solve and assembling of output */ 12280c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12290c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12300c7d97c5SJed Brown ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 12318eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1232df187020SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1233df187020SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1234df187020SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1235df187020SStefano Zampini ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 12360c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12370c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12380c7d97c5SJed Brown PetscFunctionReturn(0); 12390c7d97c5SJed Brown } 124050efa1b5SStefano Zampini 124150efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */ 124250efa1b5SStefano Zampini /* 124350efa1b5SStefano Zampini PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector. 124450efa1b5SStefano Zampini 124550efa1b5SStefano Zampini Input Parameters: 124650efa1b5SStefano Zampini . pc - the preconditioner context 124750efa1b5SStefano Zampini . r - input vector (global) 124850efa1b5SStefano Zampini 124950efa1b5SStefano Zampini Output Parameter: 125050efa1b5SStefano Zampini . z - output vector (global) 125150efa1b5SStefano Zampini 125250efa1b5SStefano Zampini Application Interface Routine: PCApplyTranspose() 125350efa1b5SStefano Zampini */ 125450efa1b5SStefano Zampini #undef __FUNCT__ 125550efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC" 125650efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z) 125750efa1b5SStefano Zampini { 125850efa1b5SStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 125950efa1b5SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 126050efa1b5SStefano Zampini PetscErrorCode ierr; 126150efa1b5SStefano Zampini const PetscScalar one = 1.0; 126250efa1b5SStefano Zampini const PetscScalar m_one = -1.0; 126350efa1b5SStefano Zampini const PetscScalar zero = 0.0; 126450efa1b5SStefano Zampini 126550efa1b5SStefano Zampini PetscFunctionBegin; 126650efa1b5SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 126750efa1b5SStefano Zampini /* First Dirichlet solve */ 126850efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 126950efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 127050efa1b5SStefano Zampini ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 127150efa1b5SStefano Zampini /* 127250efa1b5SStefano Zampini Assembling right hand side for BDDC operator 127350efa1b5SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 127450efa1b5SStefano Zampini - pcis->vec1_B the interface part of the global vector z 127550efa1b5SStefano Zampini */ 127650efa1b5SStefano Zampini ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 127750efa1b5SStefano Zampini ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 127850efa1b5SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 127950efa1b5SStefano Zampini ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 128050efa1b5SStefano Zampini ierr = VecCopy(r,z);CHKERRQ(ierr); 128150efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 128250efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 128350efa1b5SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 128450efa1b5SStefano Zampini } else { 128550efa1b5SStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 128650efa1b5SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 128750efa1b5SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 128850efa1b5SStefano Zampini } 128950efa1b5SStefano Zampini 129050efa1b5SStefano Zampini /* Apply interface preconditioner 129150efa1b5SStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1292dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr); 129350efa1b5SStefano Zampini 129450efa1b5SStefano Zampini /* Apply transpose of partition of unity operator */ 129550efa1b5SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 129650efa1b5SStefano Zampini 129750efa1b5SStefano Zampini /* Second Dirichlet solve and assembling of output */ 129850efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 129950efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 130050efa1b5SStefano Zampini ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 130150efa1b5SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1302b0147a47SStefano Zampini ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1303b0147a47SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1304b0147a47SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1305b0147a47SStefano Zampini ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 130650efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 130750efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 130850efa1b5SStefano Zampini PetscFunctionReturn(0); 130950efa1b5SStefano Zampini } 1310da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 1311674ae819SStefano Zampini 1312da1bb401SStefano Zampini #undef __FUNCT__ 1313da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC" 1314da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc) 1315da1bb401SStefano Zampini { 1316da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1317da1bb401SStefano Zampini PetscErrorCode ierr; 1318da1bb401SStefano Zampini 1319da1bb401SStefano Zampini PetscFunctionBegin; 1320da1bb401SStefano Zampini /* free data created by PCIS */ 1321da1bb401SStefano Zampini ierr = PCISDestroy(pc);CHKERRQ(ierr); 1322674ae819SStefano Zampini /* free BDDC custom data */ 1323674ae819SStefano Zampini ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 1324674ae819SStefano Zampini /* destroy objects related to topography */ 1325674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 1326674ae819SStefano Zampini /* free allocated graph structure */ 1327da1bb401SStefano Zampini ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 1328674ae819SStefano Zampini /* free data for scaling operator */ 1329674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 1330674ae819SStefano Zampini /* free solvers stuff */ 1331674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 133262a6ff1dSStefano Zampini /* free global vectors needed in presolve */ 133362a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 133462a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 133539e2fb2aSStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr); 13363425bc38SStefano Zampini /* remove functions */ 1337674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 1338bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 13392b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 1340b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 13412b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 1342bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 1343bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 134482ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1345bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 134682ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1347bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 134882ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1349bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1350785d1243SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1351bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 135263602bcaSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr); 1353bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 1354bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 1355bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 1356bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 1357674ae819SStefano Zampini /* Free the private data structure */ 1358674ae819SStefano Zampini ierr = PetscFree(pc->data);CHKERRQ(ierr); 1359da1bb401SStefano Zampini PetscFunctionReturn(0); 1360da1bb401SStefano Zampini } 13613425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 13621e6b0712SBarry Smith 13633425bc38SStefano Zampini #undef __FUNCT__ 13643425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 13653425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 13663425bc38SStefano Zampini { 1367674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 13683425bc38SStefano Zampini PC_IS* pcis; 13693425bc38SStefano Zampini PC_BDDC* pcbddc; 13703425bc38SStefano Zampini PetscErrorCode ierr; 13710c7d97c5SJed Brown 13723425bc38SStefano Zampini PetscFunctionBegin; 13733425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 13743425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 13753425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 13763425bc38SStefano Zampini 13773425bc38SStefano Zampini /* change of basis for physical rhs if needed 13783425bc38SStefano Zampini It also changes the rhs in case of dirichlet boundaries */ 13793308cffdSStefano Zampini ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 13803425bc38SStefano Zampini /* store vectors for computation of fetidp final solution */ 13813425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 13823425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1383fb223d50SStefano Zampini /* scale rhs since it should be unassembled */ 1384fb223d50SStefano Zampini /* TODO use counter scaling? (also below) */ 13853425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 13863425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1387674ae819SStefano Zampini /* Apply partition of unity */ 13883425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1389674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 13908eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 13913425bc38SStefano Zampini /* compute partially subassembled Schur complement right-hand side */ 13923425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 13933425bc38SStefano Zampini ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 13943425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 13953425bc38SStefano Zampini ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 13963425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 13973425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1398674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 13993425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14003425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14013425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 14023425bc38SStefano Zampini } 14033425bc38SStefano Zampini /* BDDC rhs */ 14043425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 14058eeda7d8SStefano Zampini if (pcbddc->switch_static) { 14063425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 14073425bc38SStefano Zampini } 14083425bc38SStefano Zampini /* apply BDDC */ 1409dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 14103425bc38SStefano Zampini /* Application of B_delta and assembling of rhs for fetidp fluxes */ 14113425bc38SStefano Zampini ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 14123425bc38SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 14133425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14143425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14153425bc38SStefano Zampini /* restore original rhs */ 14163425bc38SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 14173425bc38SStefano Zampini PetscFunctionReturn(0); 14183425bc38SStefano Zampini } 14191e6b0712SBarry Smith 14203425bc38SStefano Zampini #undef __FUNCT__ 14213425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 14223425bc38SStefano Zampini /*@ 142328509bceSStefano Zampini PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system 14243425bc38SStefano Zampini 14253425bc38SStefano Zampini Collective 14263425bc38SStefano Zampini 14273425bc38SStefano Zampini Input Parameters: 142828509bceSStefano Zampini + fetidp_mat - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators 142928509bceSStefano Zampini . standard_rhs - the right-hand side for your linear system 14303425bc38SStefano Zampini 14313425bc38SStefano Zampini Output Parameters: 143228509bceSStefano Zampini - fetidp_flux_rhs - the right-hand side for the FETIDP linear system 14333425bc38SStefano Zampini 14343425bc38SStefano Zampini Level: developer 14353425bc38SStefano Zampini 14363425bc38SStefano Zampini Notes: 14373425bc38SStefano Zampini 143828509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 14393425bc38SStefano Zampini @*/ 14403425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 14413425bc38SStefano Zampini { 1442674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 14433425bc38SStefano Zampini PetscErrorCode ierr; 14443425bc38SStefano Zampini 14453425bc38SStefano Zampini PetscFunctionBegin; 14463425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 14473425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 14483425bc38SStefano Zampini PetscFunctionReturn(0); 14493425bc38SStefano Zampini } 14503425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 14511e6b0712SBarry Smith 14523425bc38SStefano Zampini #undef __FUNCT__ 14533425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 14543425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 14553425bc38SStefano Zampini { 1456674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 14573425bc38SStefano Zampini PC_IS* pcis; 14583425bc38SStefano Zampini PC_BDDC* pcbddc; 14593425bc38SStefano Zampini PetscErrorCode ierr; 14603425bc38SStefano Zampini 14613425bc38SStefano Zampini PetscFunctionBegin; 14623425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 14633425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 14643425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 14653425bc38SStefano Zampini 14663425bc38SStefano Zampini /* apply B_delta^T */ 14673425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14683425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14693425bc38SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 14703425bc38SStefano Zampini /* compute rhs for BDDC application */ 14713425bc38SStefano Zampini ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 14728eeda7d8SStefano Zampini if (pcbddc->switch_static) { 14733425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 14743425bc38SStefano Zampini } 14753425bc38SStefano Zampini /* apply BDDC */ 1476dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 14773425bc38SStefano Zampini /* put values into standard global vector */ 14783425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14793425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14808eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 14813425bc38SStefano Zampini /* compute values into the interior if solved for the partially subassembled Schur complement */ 14823425bc38SStefano Zampini ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 14833425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 14843425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 14853425bc38SStefano Zampini } 14863425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14873425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14883425bc38SStefano Zampini /* final change of basis if needed 14893425bc38SStefano Zampini Is also sums the dirichlet part removed during RHS assembling */ 14903308cffdSStefano Zampini ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 14913425bc38SStefano Zampini PetscFunctionReturn(0); 14923425bc38SStefano Zampini } 14931e6b0712SBarry Smith 14943425bc38SStefano Zampini #undef __FUNCT__ 14953425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 14963425bc38SStefano Zampini /*@ 149728509bceSStefano Zampini PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system 14983425bc38SStefano Zampini 14993425bc38SStefano Zampini Collective 15003425bc38SStefano Zampini 15013425bc38SStefano Zampini Input Parameters: 150228509bceSStefano Zampini + fetidp_mat - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators 150328509bceSStefano Zampini . fetidp_flux_sol - the solution of the FETIDP linear system 15043425bc38SStefano Zampini 15053425bc38SStefano Zampini Output Parameters: 150628509bceSStefano Zampini - standard_sol - the solution defined on the physical domain 15073425bc38SStefano Zampini 15083425bc38SStefano Zampini Level: developer 15093425bc38SStefano Zampini 15103425bc38SStefano Zampini Notes: 15113425bc38SStefano Zampini 151228509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 15133425bc38SStefano Zampini @*/ 15143425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 15153425bc38SStefano Zampini { 1516674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 15173425bc38SStefano Zampini PetscErrorCode ierr; 15183425bc38SStefano Zampini 15193425bc38SStefano Zampini PetscFunctionBegin; 15203425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15213425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 15223425bc38SStefano Zampini PetscFunctionReturn(0); 15233425bc38SStefano Zampini } 15243425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 15251e6b0712SBarry Smith 1526f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1527edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec); 1528f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1529f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1530edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec); 1531f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1532674ae819SStefano Zampini 15333425bc38SStefano Zampini #undef __FUNCT__ 15343425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 15353425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 15363425bc38SStefano Zampini { 1537674ae819SStefano Zampini 1538674ae819SStefano Zampini FETIDPMat_ctx fetidpmat_ctx; 15393425bc38SStefano Zampini Mat newmat; 1540674ae819SStefano Zampini FETIDPPC_ctx fetidppc_ctx; 15413425bc38SStefano Zampini PC newpc; 1542ce94432eSBarry Smith MPI_Comm comm; 15433425bc38SStefano Zampini PetscErrorCode ierr; 15443425bc38SStefano Zampini 15453425bc38SStefano Zampini PetscFunctionBegin; 1546ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 15473425bc38SStefano Zampini /* FETIDP linear matrix */ 15483425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 15493425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 15503425bc38SStefano Zampini ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 15513425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 1552edf7251bSStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr); 15533425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 15543425bc38SStefano Zampini ierr = MatSetUp(newmat);CHKERRQ(ierr); 15553425bc38SStefano Zampini /* FETIDP preconditioner */ 15563425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 15573425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 15583425bc38SStefano Zampini ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 15593425bc38SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 15603425bc38SStefano Zampini ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 15613425bc38SStefano Zampini ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 1562edf7251bSStefano Zampini ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr); 15633425bc38SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 156423ee1639SBarry Smith ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr); 15653425bc38SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 15663425bc38SStefano Zampini /* return pointers for objects created */ 15673425bc38SStefano Zampini *fetidp_mat=newmat; 15683425bc38SStefano Zampini *fetidp_pc=newpc; 15693425bc38SStefano Zampini PetscFunctionReturn(0); 15703425bc38SStefano Zampini } 15711e6b0712SBarry Smith 15723425bc38SStefano Zampini #undef __FUNCT__ 15733425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 15743425bc38SStefano Zampini /*@ 157528509bceSStefano Zampini PCBDDCCreateFETIDPOperators - Create operators for FETIDP 15763425bc38SStefano Zampini 15773425bc38SStefano Zampini Collective 15783425bc38SStefano Zampini 15793425bc38SStefano Zampini Input Parameters: 158028509bceSStefano Zampini + pc - the BDDC preconditioning context already setup 158128509bceSStefano Zampini 158228509bceSStefano Zampini Output Parameters: 158328509bceSStefano Zampini . fetidp_mat - shell FETIDP matrix object 158428509bceSStefano Zampini . fetidp_pc - shell Dirichlet preconditioner for FETIDP matrix 158528509bceSStefano Zampini 158628509bceSStefano Zampini Options Database Keys: 158728509bceSStefano Zampini - -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers 15883425bc38SStefano Zampini 15893425bc38SStefano Zampini Level: developer 15903425bc38SStefano Zampini 15913425bc38SStefano Zampini Notes: 159228509bceSStefano Zampini Currently the only operation provided for FETIDP matrix is MatMult 15933425bc38SStefano Zampini 15943425bc38SStefano Zampini .seealso: PCBDDC 15953425bc38SStefano Zampini @*/ 15963425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 15973425bc38SStefano Zampini { 15983425bc38SStefano Zampini PetscErrorCode ierr; 15993425bc38SStefano Zampini 16003425bc38SStefano Zampini PetscFunctionBegin; 16013425bc38SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16023425bc38SStefano Zampini if (pc->setupcalled) { 1603516d51a7SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1604f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 16053425bc38SStefano Zampini PetscFunctionReturn(0); 16063425bc38SStefano Zampini } 16070c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 1608da1bb401SStefano Zampini /*MC 1609da1bb401SStefano Zampini PCBDDC - Balancing Domain Decomposition by Constraints. 16100c7d97c5SJed Brown 161128509bceSStefano Zampini An implementation of the BDDC preconditioner based on 161228509bceSStefano Zampini 161328509bceSStefano Zampini .vb 161428509bceSStefano Zampini [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 161528509bceSStefano 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 161628509bceSStefano Zampini [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 161728509bceSStefano Zampini .ve 161828509bceSStefano Zampini 161928509bceSStefano Zampini The matrix to be preconditioned (Pmat) must be of type MATIS. 162028509bceSStefano Zampini 1621b6fdb6dfSStefano Zampini Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 162228509bceSStefano Zampini 162328509bceSStefano Zampini It also works with unsymmetric and indefinite problems. 162428509bceSStefano Zampini 1625b6fdb6dfSStefano 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. 1626b6fdb6dfSStefano Zampini 162728509bceSStefano Zampini Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace 162828509bceSStefano Zampini 162928509bceSStefano Zampini Boundary nodes are split in vertices, edges and faces using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph 163028509bceSStefano Zampini 163128509bceSStefano Zampini Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace. 163228509bceSStefano Zampini 1633b6fdb6dfSStefano 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. 163428509bceSStefano Zampini 163528509bceSStefano Zampini The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object. 163628509bceSStefano Zampini 1637da1bb401SStefano Zampini Options Database Keys: 163828509bceSStefano Zampini 163928509bceSStefano Zampini . -pc_bddc_use_vertices <1> - use or not vertices in primal space 164028509bceSStefano Zampini . -pc_bddc_use_edges <1> - use or not edges in primal space 1641b6fdb6dfSStefano Zampini . -pc_bddc_use_faces <0> - use or not faces in primal space 164228509bceSStefano Zampini . -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only) 164328509bceSStefano Zampini . -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested 164428509bceSStefano Zampini . -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1]) 164528509bceSStefano Zampini . -pc_bddc_levels <0> - maximum number of levels for multilevel 164628509bceSStefano Zampini . -pc_bddc_coarsening_ratio - H/h ratio at the coarser level 164728509bceSStefano Zampini - -pc_bddc_check_level <0> - set verbosity level of debugging output 164828509bceSStefano Zampini 164928509bceSStefano Zampini Options for Dirichlet, Neumann or coarse solver can be set with 165028509bceSStefano Zampini .vb 165128509bceSStefano Zampini -pc_bddc_dirichlet_ 165228509bceSStefano Zampini -pc_bddc_neumann_ 165328509bceSStefano Zampini -pc_bddc_coarse_ 165428509bceSStefano Zampini .ve 165528509bceSStefano Zampini e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg 165628509bceSStefano Zampini 165728509bceSStefano Zampini When using a multilevel approach, solvers' options at the N-th level can be specified as 165828509bceSStefano Zampini .vb 1659312be037SStefano Zampini -pc_bddc_dirichlet_lN_ 1660312be037SStefano Zampini -pc_bddc_neumann_lN_ 1661312be037SStefano Zampini -pc_bddc_coarse_lN_ 166228509bceSStefano Zampini .ve 1663312be037SStefano Zampini Note that level number ranges from the finest 0 to the coarsest N. 1664da1bb401SStefano Zampini 1665da1bb401SStefano Zampini Level: intermediate 1666da1bb401SStefano Zampini 1667b6fdb6dfSStefano Zampini Developer notes: 1668da1bb401SStefano Zampini 166928509bceSStefano Zampini New deluxe scaling operator will be available soon. 1670da1bb401SStefano Zampini 1671da1bb401SStefano Zampini Contributed by Stefano Zampini 1672da1bb401SStefano Zampini 1673da1bb401SStefano Zampini .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1674da1bb401SStefano Zampini M*/ 1675b2573a8aSBarry Smith 1676da1bb401SStefano Zampini #undef __FUNCT__ 1677da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC" 16788cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1679da1bb401SStefano Zampini { 1680da1bb401SStefano Zampini PetscErrorCode ierr; 1681da1bb401SStefano Zampini PC_BDDC *pcbddc; 1682da1bb401SStefano Zampini 1683da1bb401SStefano Zampini PetscFunctionBegin; 1684da1bb401SStefano Zampini /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1685b00a9115SJed Brown ierr = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr); 1686da1bb401SStefano Zampini pc->data = (void*)pcbddc; 1687da1bb401SStefano Zampini 1688da1bb401SStefano Zampini /* create PCIS data structure */ 1689da1bb401SStefano Zampini ierr = PCISCreate(pc);CHKERRQ(ierr); 1690da1bb401SStefano Zampini 169147d04d0dSStefano Zampini /* BDDC customization */ 1692*08a5cf49SStefano Zampini pcbddc->use_local_adj = PETSC_TRUE; 169347d04d0dSStefano Zampini pcbddc->use_vertices = PETSC_TRUE; 169447d04d0dSStefano Zampini pcbddc->use_edges = PETSC_TRUE; 169547d04d0dSStefano Zampini pcbddc->use_faces = PETSC_FALSE; 169647d04d0dSStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 169747d04d0dSStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 169847d04d0dSStefano Zampini pcbddc->switch_static = PETSC_FALSE; 169947d04d0dSStefano Zampini pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 170047d04d0dSStefano Zampini pcbddc->dbg_flag = 0; 170147d04d0dSStefano Zampini 170239e2fb2aSStefano Zampini pcbddc->BtoNmap = 0; 1703727cdba6SStefano Zampini pcbddc->local_primal_size = 0; 1704e9189074SStefano Zampini pcbddc->n_vertices = 0; 1705e9189074SStefano Zampini pcbddc->n_actual_vertices = 0; 1706e9189074SStefano Zampini pcbddc->n_constraints = 0; 1707727cdba6SStefano Zampini pcbddc->primal_indices_local_idxs = 0; 1708fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_FALSE; 170968457ee5SStefano Zampini pcbddc->coarse_size = -1; 1710f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1711727cdba6SStefano Zampini pcbddc->new_primal_space_local = PETSC_FALSE; 1712f4ddd8eeSStefano Zampini pcbddc->global_primal_indices = 0; 1713f4ddd8eeSStefano Zampini pcbddc->onearnullspace = 0; 1714f4ddd8eeSStefano Zampini pcbddc->onearnullvecs_state = 0; 1715674ae819SStefano Zampini pcbddc->user_primal_vertices = 0; 17160bdf917eSStefano Zampini pcbddc->NullSpace = 0; 17173972b0daSStefano Zampini pcbddc->temp_solution = 0; 1718534831adSStefano Zampini pcbddc->original_rhs = 0; 1719534831adSStefano Zampini pcbddc->local_mat = 0; 1720534831adSStefano Zampini pcbddc->ChangeOfBasisMatrix = 0; 1721da1bb401SStefano Zampini pcbddc->coarse_vec = 0; 1722da1bb401SStefano Zampini pcbddc->coarse_rhs = 0; 1723da1bb401SStefano Zampini pcbddc->coarse_ksp = 0; 1724da1bb401SStefano Zampini pcbddc->coarse_phi_B = 0; 1725da1bb401SStefano Zampini pcbddc->coarse_phi_D = 0; 172615aaf578SStefano Zampini pcbddc->coarse_psi_B = 0; 172715aaf578SStefano Zampini pcbddc->coarse_psi_D = 0; 1728da1bb401SStefano Zampini pcbddc->vec1_P = 0; 1729da1bb401SStefano Zampini pcbddc->vec1_R = 0; 1730da1bb401SStefano Zampini pcbddc->vec2_R = 0; 1731da1bb401SStefano Zampini pcbddc->local_auxmat1 = 0; 1732da1bb401SStefano Zampini pcbddc->local_auxmat2 = 0; 1733da1bb401SStefano Zampini pcbddc->R_to_B = 0; 1734da1bb401SStefano Zampini pcbddc->R_to_D = 0; 1735da1bb401SStefano Zampini pcbddc->ksp_D = 0; 1736da1bb401SStefano Zampini pcbddc->ksp_R = 0; 1737da1bb401SStefano Zampini pcbddc->NeumannBoundaries = 0; 1738785d1243SStefano Zampini pcbddc->NeumannBoundariesLocal = 0; 1739785d1243SStefano Zampini pcbddc->DirichletBoundaries = 0; 1740785d1243SStefano Zampini pcbddc->DirichletBoundariesLocal = 0; 174160d44989SStefano Zampini pcbddc->user_provided_isfordofs = PETSC_FALSE; 174260d44989SStefano Zampini pcbddc->n_ISForDofs = 0; 174363602bcaSStefano Zampini pcbddc->n_ISForDofsLocal = 0; 1744da1bb401SStefano Zampini pcbddc->ISForDofs = 0; 174563602bcaSStefano Zampini pcbddc->ISForDofsLocal = 0; 1746da1bb401SStefano Zampini pcbddc->ConstraintMatrix = 0; 174785c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 174847d04d0dSStefano Zampini pcbddc->coarse_loc_to_glob = 0; 174947d04d0dSStefano Zampini pcbddc->coarsening_ratio = 8; 17504fad6a16SStefano Zampini pcbddc->current_level = 0; 17512b510759SStefano Zampini pcbddc->max_levels = 0; 1752323d395dSStefano Zampini pcbddc->use_coarse_estimates = PETSC_FALSE; 1753f3bde8b3SStefano Zampini pcbddc->coarse_subassembling = 0; 1754323d395dSStefano Zampini pcbddc->coarse_subassembling_init = 0; 1755da1bb401SStefano Zampini 1756674ae819SStefano Zampini /* create local graph structure */ 1757674ae819SStefano Zampini ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1758674ae819SStefano Zampini 1759674ae819SStefano Zampini /* scaling */ 1760674ae819SStefano Zampini pcbddc->use_deluxe_scaling = PETSC_FALSE; 1761*08a5cf49SStefano Zampini pcbddc->deluxe_ctx = 0; 1762674ae819SStefano Zampini pcbddc->work_scaling = 0; 1763da1bb401SStefano Zampini 1764da1bb401SStefano Zampini /* function pointers */ 1765da1bb401SStefano Zampini pc->ops->apply = PCApply_BDDC; 176693bd9ae7SStefano Zampini pc->ops->applytranspose = PCApplyTranspose_BDDC; 1767da1bb401SStefano Zampini pc->ops->setup = PCSetUp_BDDC; 1768da1bb401SStefano Zampini pc->ops->destroy = PCDestroy_BDDC; 1769da1bb401SStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1770da1bb401SStefano Zampini pc->ops->view = 0; 1771da1bb401SStefano Zampini pc->ops->applyrichardson = 0; 1772da1bb401SStefano Zampini pc->ops->applysymmetricleft = 0; 1773da1bb401SStefano Zampini pc->ops->applysymmetricright = 0; 1774534831adSStefano Zampini pc->ops->presolve = PCPreSolve_BDDC; 1775534831adSStefano Zampini pc->ops->postsolve = PCPostSolve_BDDC; 1776da1bb401SStefano Zampini 1777da1bb401SStefano Zampini /* composing function */ 1778674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1779bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 17802b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1781b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 17822b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1783bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1784bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 178582ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1786bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 178782ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1788bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 178982ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1790bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 179182ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1792bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 179363602bcaSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr); 1794bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1795bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1796bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1797bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1798da1bb401SStefano Zampini PetscFunctionReturn(0); 1799da1bb401SStefano Zampini } 18003425bc38SStefano Zampini 1801