153cdbc3dSStefano Zampini /* TODOLIST 2da1bb401SStefano Zampini DofSplitting and DM attached to pc? 3da1bb401SStefano Zampini Change SetNeumannBoundaries to SetNeumannBoundariesLocal and provide new SetNeumannBoundaries (same Dirichlet) 4a0ba757dSStefano Zampini change how to deal with the coarse problem (PCBDDCSetCoarseEnvironment): 5a0ba757dSStefano Zampini - simplify coarse problem structure -> PCBDDC or PCREDUDANT, nothing else -> same comm for all levels? 6a0ba757dSStefano Zampini - remove coarse enums and allow use of PCBDDCGetCoarseKSP 7674ae819SStefano Zampini - remove metis dependency -> use MatPartitioning for multilevel -> Assemble serial adjacency in PCBDDCAnalyzeInterface? 8a0ba757dSStefano Zampini code refactoring: 9a0ba757dSStefano Zampini - pick up better names for static functions 10a0ba757dSStefano Zampini change options structure: 11a0ba757dSStefano Zampini - insert BDDC into MG framework? 12a0ba757dSStefano Zampini provide other ops? Ask to developers 13a0ba757dSStefano Zampini remove all unused printf 14a0ba757dSStefano Zampini man pages 1553cdbc3dSStefano Zampini */ 160c7d97c5SJed Brown 1753cdbc3dSStefano Zampini /* ---------------------------------------------------------------------------------------------------------------------------------------------- 180c7d97c5SJed Brown Implementation of BDDC preconditioner based on: 190c7d97c5SJed Brown C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 2053cdbc3dSStefano Zampini ---------------------------------------------------------------------------------------------------------------------------------------------- */ 2153cdbc3dSStefano Zampini 22674ae819SStefano Zampini #include "bddc.h" /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 23674ae819SStefano Zampini #include "bddcprivate.h" 243b03a366Sstefano_zampini #include <petscblaslapack.h> 25674ae819SStefano Zampini 260c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 270c7d97c5SJed Brown #undef __FUNCT__ 280c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC" 290c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc) 300c7d97c5SJed Brown { 310c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 320c7d97c5SJed Brown PetscErrorCode ierr; 330c7d97c5SJed Brown 340c7d97c5SJed Brown PetscFunctionBegin; 350c7d97c5SJed Brown ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 368eeda7d8SStefano Zampini /* Verbose debugging */ 378eeda7d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 388eeda7d8SStefano Zampini /* Primal space cumstomization */ 398eeda7d8SStefano 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); 408eeda7d8SStefano 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); 418eeda7d8SStefano 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); 428eeda7d8SStefano Zampini /* Change of basis */ 438eeda7d8SStefano 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); 448eeda7d8SStefano 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); 45674ae819SStefano Zampini if (!pcbddc->use_change_of_basis) { 46674ae819SStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 47674ae819SStefano Zampini } 488eeda7d8SStefano Zampini /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 498eeda7d8SStefano 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); 500298fd71SBarry Smith ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 51*2b510759SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 52674ae819SStefano 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); 530c7d97c5SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 540c7d97c5SJed Brown PetscFunctionReturn(0); 550c7d97c5SJed Brown } 560c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 57674ae819SStefano Zampini #undef __FUNCT__ 58674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 59674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 60674ae819SStefano Zampini { 61674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 62674ae819SStefano Zampini PetscErrorCode ierr; 631e6b0712SBarry Smith 64674ae819SStefano Zampini PetscFunctionBegin; 65674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 66674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 67674ae819SStefano Zampini pcbddc->user_primal_vertices = PrimalVertices; 68674ae819SStefano Zampini PetscFunctionReturn(0); 69674ae819SStefano Zampini } 70674ae819SStefano Zampini #undef __FUNCT__ 71674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 72674ae819SStefano Zampini /*@ 73674ae819SStefano Zampini PCBDDCSetPrimalVerticesLocalIS - Set user defined primal vertices in PCBDDC. 74674ae819SStefano Zampini 75674ae819SStefano Zampini Not collective 76674ae819SStefano Zampini 77674ae819SStefano Zampini Input Parameters: 78674ae819SStefano Zampini + pc - the preconditioning context 79674ae819SStefano Zampini - PrimalVertices - index sets of primal vertices in local numbering 80674ae819SStefano Zampini 81674ae819SStefano Zampini Level: intermediate 82674ae819SStefano Zampini 83674ae819SStefano Zampini Notes: 84674ae819SStefano Zampini 85674ae819SStefano Zampini .seealso: PCBDDC 86674ae819SStefano Zampini @*/ 87674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 88674ae819SStefano Zampini { 89674ae819SStefano Zampini PetscErrorCode ierr; 90674ae819SStefano Zampini 91674ae819SStefano Zampini PetscFunctionBegin; 92674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 93674ae819SStefano Zampini PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 94674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 95674ae819SStefano Zampini PetscFunctionReturn(0); 96674ae819SStefano Zampini } 97674ae819SStefano Zampini /* -------------------------------------------------------------------------- */ 980c7d97c5SJed Brown #undef __FUNCT__ 994fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 1004fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 1014fad6a16SStefano Zampini { 1024fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1034fad6a16SStefano Zampini 1044fad6a16SStefano Zampini PetscFunctionBegin; 1054fad6a16SStefano Zampini pcbddc->coarsening_ratio = k; 1064fad6a16SStefano Zampini PetscFunctionReturn(0); 1074fad6a16SStefano Zampini } 1081e6b0712SBarry Smith 1094fad6a16SStefano Zampini #undef __FUNCT__ 1104fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio" 1114fad6a16SStefano Zampini /*@ 1124fad6a16SStefano Zampini PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel coarsening 1134fad6a16SStefano Zampini 1144fad6a16SStefano Zampini Logically collective on PC 1154fad6a16SStefano Zampini 1164fad6a16SStefano Zampini Input Parameters: 1174fad6a16SStefano Zampini + pc - the preconditioning context 1184fad6a16SStefano Zampini - k - coarsening ratio 1194fad6a16SStefano Zampini 1204fad6a16SStefano Zampini Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level. 1214fad6a16SStefano Zampini 1224fad6a16SStefano Zampini Level: intermediate 1234fad6a16SStefano Zampini 1244fad6a16SStefano Zampini Notes: 1254fad6a16SStefano Zampini 1264fad6a16SStefano Zampini .seealso: PCBDDC 1274fad6a16SStefano Zampini @*/ 1284fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 1294fad6a16SStefano Zampini { 1304fad6a16SStefano Zampini PetscErrorCode ierr; 1314fad6a16SStefano Zampini 1324fad6a16SStefano Zampini PetscFunctionBegin; 1334fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 134*2b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,k,2); 1354fad6a16SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 1364fad6a16SStefano Zampini PetscFunctionReturn(0); 1374fad6a16SStefano Zampini } 138*2b510759SStefano Zampini 139*2b510759SStefano Zampini /* Set level is not a public function */ 140*2b510759SStefano Zampini #undef __FUNCT__ 141*2b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel" 142*2b510759SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 143*2b510759SStefano Zampini { 144*2b510759SStefano Zampini PetscErrorCode ierr; 145*2b510759SStefano Zampini 146*2b510759SStefano Zampini PetscFunctionBegin; 147*2b510759SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 148*2b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,level,2); 149*2b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 150*2b510759SStefano Zampini PetscFunctionReturn(0); 151*2b510759SStefano Zampini } 1521e6b0712SBarry Smith 1534fad6a16SStefano Zampini #undef __FUNCT__ 154*2b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC" 155*2b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 1564fad6a16SStefano Zampini { 1574fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1584fad6a16SStefano Zampini 1594fad6a16SStefano Zampini PetscFunctionBegin; 160*2b510759SStefano Zampini pcbddc->current_level=level; 1614fad6a16SStefano Zampini PetscFunctionReturn(0); 1624fad6a16SStefano Zampini } 1631e6b0712SBarry Smith 1644fad6a16SStefano Zampini #undef __FUNCT__ 165*2b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC" 166*2b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 167*2b510759SStefano Zampini { 168*2b510759SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 169*2b510759SStefano Zampini 170*2b510759SStefano Zampini PetscFunctionBegin; 171*2b510759SStefano Zampini pcbddc->max_levels = levels; 172*2b510759SStefano Zampini PetscFunctionReturn(0); 173*2b510759SStefano Zampini } 174*2b510759SStefano Zampini 175*2b510759SStefano Zampini #undef __FUNCT__ 176*2b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels" 1774fad6a16SStefano Zampini /*@ 178*2b510759SStefano Zampini PCBDDCSetLevels - Sets the maximum number of levels within the multilevel approach. 1794fad6a16SStefano Zampini 1804fad6a16SStefano Zampini Logically collective on PC 1814fad6a16SStefano Zampini 1824fad6a16SStefano Zampini Input Parameters: 1834fad6a16SStefano Zampini + pc - the preconditioning context 184*2b510759SStefano Zampini - levels - the maximum number of levels 1854fad6a16SStefano Zampini 186*2b510759SStefano Zampini Default value is 0, i.e. coarse problem will be solved exactly 1874fad6a16SStefano Zampini 1884fad6a16SStefano Zampini Level: intermediate 1894fad6a16SStefano Zampini 1904fad6a16SStefano Zampini Notes: 1914fad6a16SStefano Zampini 1924fad6a16SStefano Zampini .seealso: PCBDDC 1934fad6a16SStefano Zampini @*/ 194*2b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 1954fad6a16SStefano Zampini { 1964fad6a16SStefano Zampini PetscErrorCode ierr; 1974fad6a16SStefano Zampini 1984fad6a16SStefano Zampini PetscFunctionBegin; 1994fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 200*2b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,levels,2); 201*2b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 2024fad6a16SStefano Zampini PetscFunctionReturn(0); 2034fad6a16SStefano Zampini } 2044fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */ 2051e6b0712SBarry Smith 2064fad6a16SStefano Zampini #undef __FUNCT__ 2070bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 2080bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 2090bdf917eSStefano Zampini { 2100bdf917eSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2110bdf917eSStefano Zampini PetscErrorCode ierr; 2120bdf917eSStefano Zampini 2130bdf917eSStefano Zampini PetscFunctionBegin; 2140bdf917eSStefano Zampini ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 2150bdf917eSStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 2160bdf917eSStefano Zampini pcbddc->NullSpace = NullSpace; 2170bdf917eSStefano Zampini PetscFunctionReturn(0); 2180bdf917eSStefano Zampini } 2191e6b0712SBarry Smith 2200bdf917eSStefano Zampini #undef __FUNCT__ 2210bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace" 2220bdf917eSStefano Zampini /*@ 2230bdf917eSStefano Zampini PCBDDCSetNullSpace - Set NullSpace of global operator of BDDC preconditioned mat. 2240bdf917eSStefano Zampini 2250bdf917eSStefano Zampini Logically collective on PC and MatNullSpace 2260bdf917eSStefano Zampini 2270bdf917eSStefano Zampini Input Parameters: 2280bdf917eSStefano Zampini + pc - the preconditioning context 2290bdf917eSStefano Zampini - NullSpace - Null space of the linear operator to be preconditioned. 2300bdf917eSStefano Zampini 2310bdf917eSStefano Zampini Level: intermediate 2320bdf917eSStefano Zampini 2330bdf917eSStefano Zampini Notes: 2340bdf917eSStefano Zampini 2350bdf917eSStefano Zampini .seealso: PCBDDC 2360bdf917eSStefano Zampini @*/ 2370bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 2380bdf917eSStefano Zampini { 2390bdf917eSStefano Zampini PetscErrorCode ierr; 2400bdf917eSStefano Zampini 2410bdf917eSStefano Zampini PetscFunctionBegin; 2420bdf917eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 243674ae819SStefano Zampini PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 244*2b510759SStefano Zampini PetscCheckSameComm(pc,1,NullSpace,2); 2450bdf917eSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 2460bdf917eSStefano Zampini PetscFunctionReturn(0); 2470bdf917eSStefano Zampini } 2480bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */ 2491e6b0712SBarry Smith 2500bdf917eSStefano Zampini #undef __FUNCT__ 2513b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 2523b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 2533b03a366Sstefano_zampini { 2543b03a366Sstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2553b03a366Sstefano_zampini PetscErrorCode ierr; 2563b03a366Sstefano_zampini 2573b03a366Sstefano_zampini PetscFunctionBegin; 2583b03a366Sstefano_zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 25936e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 26036e030ebSStefano Zampini pcbddc->DirichletBoundaries=DirichletBoundaries; 2613b03a366Sstefano_zampini PetscFunctionReturn(0); 2623b03a366Sstefano_zampini } 2631e6b0712SBarry Smith 2643b03a366Sstefano_zampini #undef __FUNCT__ 2653b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 2663b03a366Sstefano_zampini /*@ 267da1bb401SStefano Zampini PCBDDCSetDirichletBoundaries - Set index set defining subdomain part (in local ordering) 268da1bb401SStefano Zampini of Dirichlet boundaries for the global problem. 2693b03a366Sstefano_zampini 2703b03a366Sstefano_zampini Not collective 2713b03a366Sstefano_zampini 2723b03a366Sstefano_zampini Input Parameters: 2733b03a366Sstefano_zampini + pc - the preconditioning context 2740298fd71SBarry Smith - DirichletBoundaries - sequential index set defining the subdomain part of Dirichlet boundaries (can be NULL) 2753b03a366Sstefano_zampini 2763b03a366Sstefano_zampini Level: intermediate 2773b03a366Sstefano_zampini 2783b03a366Sstefano_zampini Notes: 2793b03a366Sstefano_zampini 2803b03a366Sstefano_zampini .seealso: PCBDDC 2813b03a366Sstefano_zampini @*/ 2823b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 2833b03a366Sstefano_zampini { 2843b03a366Sstefano_zampini PetscErrorCode ierr; 2853b03a366Sstefano_zampini 2863b03a366Sstefano_zampini PetscFunctionBegin; 2873b03a366Sstefano_zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 288674ae819SStefano Zampini PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 2893b03a366Sstefano_zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 2903b03a366Sstefano_zampini PetscFunctionReturn(0); 2913b03a366Sstefano_zampini } 2923b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */ 2931e6b0712SBarry Smith 2943b03a366Sstefano_zampini #undef __FUNCT__ 2950c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 29653cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 2970c7d97c5SJed Brown { 2980c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 29953cdbc3dSStefano Zampini PetscErrorCode ierr; 3000c7d97c5SJed Brown 3010c7d97c5SJed Brown PetscFunctionBegin; 30253cdbc3dSStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 30336e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 30436e030ebSStefano Zampini pcbddc->NeumannBoundaries=NeumannBoundaries; 3050c7d97c5SJed Brown PetscFunctionReturn(0); 3060c7d97c5SJed Brown } 3071e6b0712SBarry Smith 3080c7d97c5SJed Brown #undef __FUNCT__ 3090c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 31057527edcSJed Brown /*@ 311da1bb401SStefano Zampini PCBDDCSetNeumannBoundaries - Set index set defining subdomain part (in local ordering) 312da1bb401SStefano Zampini of Neumann boundaries for the global problem. 31357527edcSJed Brown 3149c0446d6SStefano Zampini Not collective 31557527edcSJed Brown 31657527edcSJed Brown Input Parameters: 31757527edcSJed Brown + pc - the preconditioning context 3180298fd71SBarry Smith - NeumannBoundaries - sequential index set defining the subdomain part of Neumann boundaries (can be NULL) 31957527edcSJed Brown 32057527edcSJed Brown Level: intermediate 32157527edcSJed Brown 32257527edcSJed Brown Notes: 32357527edcSJed Brown 32457527edcSJed Brown .seealso: PCBDDC 32557527edcSJed Brown @*/ 32653cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 3270c7d97c5SJed Brown { 3280c7d97c5SJed Brown PetscErrorCode ierr; 3290c7d97c5SJed Brown 3300c7d97c5SJed Brown PetscFunctionBegin; 3310c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 332674ae819SStefano Zampini PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 33353cdbc3dSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 33453cdbc3dSStefano Zampini PetscFunctionReturn(0); 33553cdbc3dSStefano Zampini } 33653cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 3371e6b0712SBarry Smith 33853cdbc3dSStefano Zampini #undef __FUNCT__ 339da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 340da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 341da1bb401SStefano Zampini { 342da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 343da1bb401SStefano Zampini 344da1bb401SStefano Zampini PetscFunctionBegin; 345da1bb401SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundaries; 346da1bb401SStefano Zampini PetscFunctionReturn(0); 347da1bb401SStefano Zampini } 3481e6b0712SBarry Smith 349da1bb401SStefano Zampini #undef __FUNCT__ 350da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 351da1bb401SStefano Zampini /*@ 352da1bb401SStefano Zampini PCBDDCGetDirichletBoundaries - Get index set defining subdomain part (in local ordering) 353da1bb401SStefano Zampini of Dirichlet boundaries for the global problem. 354da1bb401SStefano Zampini 355da1bb401SStefano Zampini Not collective 356da1bb401SStefano Zampini 357da1bb401SStefano Zampini Input Parameters: 358da1bb401SStefano Zampini + pc - the preconditioning context 359da1bb401SStefano Zampini 360da1bb401SStefano Zampini Output Parameters: 361da1bb401SStefano Zampini + DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 362da1bb401SStefano Zampini 363da1bb401SStefano Zampini Level: intermediate 364da1bb401SStefano Zampini 365da1bb401SStefano Zampini Notes: 366da1bb401SStefano Zampini 367da1bb401SStefano Zampini .seealso: PCBDDC 368da1bb401SStefano Zampini @*/ 369da1bb401SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 370da1bb401SStefano Zampini { 371da1bb401SStefano Zampini PetscErrorCode ierr; 372da1bb401SStefano Zampini 373da1bb401SStefano Zampini PetscFunctionBegin; 374da1bb401SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 375da1bb401SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 376da1bb401SStefano Zampini PetscFunctionReturn(0); 377da1bb401SStefano Zampini } 378da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 3791e6b0712SBarry Smith 380da1bb401SStefano Zampini #undef __FUNCT__ 38153cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 38253cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 38353cdbc3dSStefano Zampini { 38453cdbc3dSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 38553cdbc3dSStefano Zampini 38653cdbc3dSStefano Zampini PetscFunctionBegin; 38753cdbc3dSStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundaries; 38853cdbc3dSStefano Zampini PetscFunctionReturn(0); 38953cdbc3dSStefano Zampini } 3901e6b0712SBarry Smith 39153cdbc3dSStefano Zampini #undef __FUNCT__ 39253cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 39353cdbc3dSStefano Zampini /*@ 394da1bb401SStefano Zampini PCBDDCGetNeumannBoundaries - Get index set defining subdomain part (in local ordering) 395da1bb401SStefano Zampini of Neumann boundaries for the global problem. 39653cdbc3dSStefano Zampini 3979c0446d6SStefano Zampini Not collective 39853cdbc3dSStefano Zampini 39953cdbc3dSStefano Zampini Input Parameters: 40053cdbc3dSStefano Zampini + pc - the preconditioning context 40153cdbc3dSStefano Zampini 40253cdbc3dSStefano Zampini Output Parameters: 40353cdbc3dSStefano Zampini + NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 40453cdbc3dSStefano Zampini 40553cdbc3dSStefano Zampini Level: intermediate 40653cdbc3dSStefano Zampini 40753cdbc3dSStefano Zampini Notes: 40853cdbc3dSStefano Zampini 40953cdbc3dSStefano Zampini .seealso: PCBDDC 41053cdbc3dSStefano Zampini @*/ 41153cdbc3dSStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 41253cdbc3dSStefano Zampini { 41353cdbc3dSStefano Zampini PetscErrorCode ierr; 41453cdbc3dSStefano Zampini 41553cdbc3dSStefano Zampini PetscFunctionBegin; 41653cdbc3dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 41753cdbc3dSStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 4180c7d97c5SJed Brown PetscFunctionReturn(0); 4190c7d97c5SJed Brown } 42036e030ebSStefano Zampini /* -------------------------------------------------------------------------- */ 4211e6b0712SBarry Smith 42236e030ebSStefano Zampini #undef __FUNCT__ 423da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 4241a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 42536e030ebSStefano Zampini { 42636e030ebSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 427da1bb401SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 428da1bb401SStefano Zampini PetscErrorCode ierr; 42936e030ebSStefano Zampini 43036e030ebSStefano Zampini PetscFunctionBegin; 431674ae819SStefano Zampini /* free old CSR */ 432674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 433fb223d50SStefano Zampini /* TODO: PCBDDCGraphSetAdjacency */ 434674ae819SStefano Zampini /* get CSR into graph structure */ 435da1bb401SStefano Zampini if (copymode == PETSC_COPY_VALUES) { 436674ae819SStefano Zampini ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr); 437674ae819SStefano Zampini ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr); 438674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 439674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 440da1bb401SStefano Zampini } else if (copymode == PETSC_OWN_POINTER) { 4411a83f524SJed Brown mat_graph->xadj = (PetscInt*)xadj; 4421a83f524SJed Brown mat_graph->adjncy = (PetscInt*)adjncy; 443674ae819SStefano Zampini } 444575ad6abSStefano Zampini mat_graph->nvtxs_csr = nvtxs; 44536e030ebSStefano Zampini PetscFunctionReturn(0); 44636e030ebSStefano Zampini } 4471e6b0712SBarry Smith 44836e030ebSStefano Zampini #undef __FUNCT__ 449da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 45036e030ebSStefano Zampini /*@ 451da1bb401SStefano Zampini PCBDDCSetLocalAdjacencyGraph - Set CSR graph of local matrix for use of PCBDDC. 45236e030ebSStefano Zampini 45336e030ebSStefano Zampini Not collective 45436e030ebSStefano Zampini 45536e030ebSStefano Zampini Input Parameters: 45636e030ebSStefano Zampini + pc - the preconditioning context 457da1bb401SStefano Zampini - nvtxs - number of local vertices of the graph 458da1bb401SStefano Zampini - xadj, adjncy - the CSR graph 459da1bb401SStefano Zampini - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. In the former case the user must free the array passed in; 460da1bb401SStefano Zampini in the latter case, memory must be obtained with PetscMalloc. 46136e030ebSStefano Zampini 46236e030ebSStefano Zampini Level: intermediate 46336e030ebSStefano Zampini 46436e030ebSStefano Zampini Notes: 46536e030ebSStefano Zampini 46636e030ebSStefano Zampini .seealso: PCBDDC 46736e030ebSStefano Zampini @*/ 4681a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 46936e030ebSStefano Zampini { 470575ad6abSStefano Zampini void (*f)(void) = 0; 47136e030ebSStefano Zampini PetscErrorCode ierr; 47236e030ebSStefano Zampini 47336e030ebSStefano Zampini PetscFunctionBegin; 47436e030ebSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 475674ae819SStefano Zampini PetscValidIntPointer(xadj,3); 476674ae819SStefano Zampini PetscValidIntPointer(xadj,4); 477674ae819SStefano Zampini if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 478674ae819SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 479da1bb401SStefano Zampini } 48036e030ebSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 481575ad6abSStefano Zampini /* free arrays if PCBDDC is not the PC type */ 482575ad6abSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 483575ad6abSStefano Zampini if (!f && copymode == PETSC_OWN_POINTER) { 484575ad6abSStefano Zampini ierr = PetscFree(xadj);CHKERRQ(ierr); 485575ad6abSStefano Zampini ierr = PetscFree(adjncy);CHKERRQ(ierr); 48636e030ebSStefano Zampini } 48736e030ebSStefano Zampini PetscFunctionReturn(0); 48836e030ebSStefano Zampini } 4899c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */ 4901e6b0712SBarry Smith 4919c0446d6SStefano Zampini #undef __FUNCT__ 4929c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 4939c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 4949c0446d6SStefano Zampini { 4959c0446d6SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 4969c0446d6SStefano Zampini PetscInt i; 4979c0446d6SStefano Zampini PetscErrorCode ierr; 4989c0446d6SStefano Zampini 4999c0446d6SStefano Zampini PetscFunctionBegin; 500da1bb401SStefano Zampini /* Destroy ISes if they were already set */ 5019c0446d6SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 5029c0446d6SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 5039c0446d6SStefano Zampini } 504d11ae9bbSstefano_zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 505da1bb401SStefano Zampini /* allocate space then set */ 5069c0446d6SStefano Zampini ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr); 5079c0446d6SStefano Zampini for (i=0;i<n_is;i++) { 508da1bb401SStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 509da1bb401SStefano Zampini pcbddc->ISForDofs[i]=ISForDofs[i]; 5109c0446d6SStefano Zampini } 5119c0446d6SStefano Zampini pcbddc->n_ISForDofs=n_is; 5129c0446d6SStefano Zampini PetscFunctionReturn(0); 5139c0446d6SStefano Zampini } 5141e6b0712SBarry Smith 5159c0446d6SStefano Zampini #undef __FUNCT__ 5169c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting" 5179c0446d6SStefano Zampini /*@ 518da1bb401SStefano Zampini PCBDDCSetDofsSplitting - Set index sets defining fields of local mat. 5199c0446d6SStefano Zampini 5209c0446d6SStefano Zampini Not collective 5219c0446d6SStefano Zampini 5229c0446d6SStefano Zampini Input Parameters: 5239c0446d6SStefano Zampini + pc - the preconditioning context 524da1bb401SStefano Zampini - n - number of index sets defining the fields 525da1bb401SStefano Zampini - IS[] - array of IS describing the fields 5269c0446d6SStefano Zampini 5279c0446d6SStefano Zampini Level: intermediate 5289c0446d6SStefano Zampini 5299c0446d6SStefano Zampini Notes: 5309c0446d6SStefano Zampini 5319c0446d6SStefano Zampini .seealso: PCBDDC 5329c0446d6SStefano Zampini @*/ 5339c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 5349c0446d6SStefano Zampini { 535*2b510759SStefano Zampini PetscInt i; 5369c0446d6SStefano Zampini PetscErrorCode ierr; 5379c0446d6SStefano Zampini 5389c0446d6SStefano Zampini PetscFunctionBegin; 5399c0446d6SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 540*2b510759SStefano Zampini for (i=0;i<n_is;i++) { 541*2b510759SStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2); 542*2b510759SStefano Zampini } 5439c0446d6SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 5449c0446d6SStefano Zampini PetscFunctionReturn(0); 5459c0446d6SStefano Zampini } 546da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 547534831adSStefano Zampini #undef __FUNCT__ 548534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC" 549534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 550534831adSStefano Zampini /* 551534831adSStefano Zampini PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 552534831adSStefano Zampini guess if a transformation of basis approach has been selected. 5539c0446d6SStefano Zampini 554534831adSStefano Zampini Input Parameter: 555534831adSStefano Zampini + pc - the preconditioner contex 556534831adSStefano Zampini 557534831adSStefano Zampini Application Interface Routine: PCPreSolve() 558534831adSStefano Zampini 559534831adSStefano Zampini Notes: 560534831adSStefano Zampini The interface routine PCPreSolve() is not usually called directly by 561534831adSStefano Zampini the user, but instead is called by KSPSolve(). 562534831adSStefano Zampini */ 563534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 564534831adSStefano Zampini { 565534831adSStefano Zampini PetscErrorCode ierr; 566534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 567534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 568534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 569534831adSStefano Zampini Mat temp_mat; 5703972b0daSStefano Zampini IS dirIS; 5713972b0daSStefano Zampini PetscInt dirsize,i,*is_indices; 5723972b0daSStefano Zampini PetscScalar *array_x,*array_diagonal; 5733972b0daSStefano Zampini Vec used_vec; 5743972b0daSStefano Zampini PetscBool guess_nonzero; 575534831adSStefano Zampini 576534831adSStefano Zampini PetscFunctionBegin; 57762a6ff1dSStefano Zampini /* Creates parallel work vectors used in presolve. */ 57862a6ff1dSStefano Zampini if (!pcbddc->original_rhs) { 57962a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 58062a6ff1dSStefano Zampini } 58162a6ff1dSStefano Zampini if (!pcbddc->temp_solution) { 58262a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 58362a6ff1dSStefano Zampini } 5843972b0daSStefano Zampini if (x) { 5853972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 5863972b0daSStefano Zampini used_vec = x; 5873972b0daSStefano Zampini } else { 5883972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 5893972b0daSStefano Zampini used_vec = pcbddc->temp_solution; 5903972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 5913972b0daSStefano Zampini } 5923972b0daSStefano Zampini /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 5933972b0daSStefano Zampini if (ksp) { 5943972b0daSStefano Zampini ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 5953972b0daSStefano Zampini if ( !guess_nonzero ) { 5963972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 5973972b0daSStefano Zampini } 5983972b0daSStefano Zampini } 5993308cffdSStefano Zampini 60062a6ff1dSStefano Zampini if (rhs) { /* TODO: wiser handling of rhs removal, which is only needed in case of zeroed rows */ 6013972b0daSStefano Zampini /* store the original rhs */ 6023972b0daSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 6033972b0daSStefano Zampini 6043972b0daSStefano Zampini /* Take into account zeroed rows -> change rhs and store solution removed */ 6053972b0daSStefano Zampini ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 6063972b0daSStefano Zampini ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 6073972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6083972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6093972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6103972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6113972b0daSStefano Zampini ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr); 6123972b0daSStefano Zampini if (dirIS) { 6133972b0daSStefano Zampini ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr); 6143972b0daSStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 6153972b0daSStefano Zampini ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 6163972b0daSStefano Zampini ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 6172fa5cd67SKarl Rupp for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 6183972b0daSStefano Zampini ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 6193972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 6203972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 6213972b0daSStefano Zampini } 6223972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6233972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 624b76ba322SStefano Zampini 6253972b0daSStefano Zampini /* remove the computed solution from the rhs */ 6263972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 6273972b0daSStefano Zampini ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 6283972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 6293308cffdSStefano Zampini } 630b76ba322SStefano Zampini 631b76ba322SStefano Zampini /* store partially computed solution and set initial guess */ 6323972b0daSStefano Zampini if (x) { 6333972b0daSStefano Zampini ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 6343972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 63515aaf578SStefano Zampini if (pcbddc->use_exact_dirichlet && !pcbddc->coarse_psi_B) { 636b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 637b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 638b76ba322SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 639b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 640b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 641b76ba322SStefano Zampini if (ksp) { 642b76ba322SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 643b76ba322SStefano Zampini } 644b76ba322SStefano Zampini } 6453972b0daSStefano Zampini } 646b76ba322SStefano Zampini 647674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 648b76ba322SStefano Zampini /* swap pointers for local matrices */ 649b76ba322SStefano Zampini temp_mat = matis->A; 650b76ba322SStefano Zampini matis->A = pcbddc->local_mat; 651b76ba322SStefano Zampini pcbddc->local_mat = temp_mat; 6523308cffdSStefano Zampini } 6533308cffdSStefano Zampini if (pcbddc->use_change_of_basis && rhs) { 654b76ba322SStefano Zampini /* Get local rhs and apply transformation of basis */ 655b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 656b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 657b76ba322SStefano Zampini /* from original basis to modified basis */ 658b76ba322SStefano Zampini ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 659b76ba322SStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 660b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 661b76ba322SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 662674ae819SStefano Zampini } 6630bdf917eSStefano Zampini if (ksp && pcbddc->NullSpace) { 664d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 665d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 666b76ba322SStefano Zampini } 6670bdf917eSStefano Zampini ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 668534831adSStefano Zampini PetscFunctionReturn(0); 669534831adSStefano Zampini } 670534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 671534831adSStefano Zampini #undef __FUNCT__ 672534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC" 673534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 674534831adSStefano Zampini /* 675534831adSStefano Zampini PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 676534831adSStefano Zampini approach has been selected. Also, restores rhs to its original state. 677534831adSStefano Zampini 678534831adSStefano Zampini Input Parameter: 679534831adSStefano Zampini + pc - the preconditioner contex 680534831adSStefano Zampini 681534831adSStefano Zampini Application Interface Routine: PCPostSolve() 682534831adSStefano Zampini 683534831adSStefano Zampini Notes: 684534831adSStefano Zampini The interface routine PCPostSolve() is not usually called directly by 685534831adSStefano Zampini the user, but instead is called by KSPSolve(). 686534831adSStefano Zampini */ 687534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 688534831adSStefano Zampini { 689534831adSStefano Zampini PetscErrorCode ierr; 690534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 691534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 692534831adSStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 693534831adSStefano Zampini Mat temp_mat; 694534831adSStefano Zampini 695534831adSStefano Zampini PetscFunctionBegin; 696674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 697534831adSStefano Zampini /* swap pointers for local matrices */ 698534831adSStefano Zampini temp_mat = matis->A; 699534831adSStefano Zampini matis->A = pcbddc->local_mat; 700534831adSStefano Zampini pcbddc->local_mat = temp_mat; 7013425bc38SStefano Zampini } 7023308cffdSStefano Zampini if (pcbddc->use_change_of_basis && x) { 703534831adSStefano Zampini /* Get Local boundary and apply transformation of basis to solution vector */ 704534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 705534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 706534831adSStefano Zampini /* from modified basis to original basis */ 707534831adSStefano Zampini ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 708534831adSStefano Zampini /* put back modified values into the global vec using INSERT_VALUES copy mode */ 709534831adSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 710534831adSStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 711534831adSStefano Zampini } 7123972b0daSStefano Zampini /* add solution removed in presolve */ 7133425bc38SStefano Zampini if (x) { 7143425bc38SStefano Zampini ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 7153425bc38SStefano Zampini } 716fb223d50SStefano Zampini /* restore rhs to its original state */ 717fb223d50SStefano Zampini if (rhs) { 718fb223d50SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 719fb223d50SStefano Zampini } 720534831adSStefano Zampini PetscFunctionReturn(0); 721534831adSStefano Zampini } 722534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 72353cdbc3dSStefano Zampini #undef __FUNCT__ 72453cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC" 7250c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 7260c7d97c5SJed Brown /* 7270c7d97c5SJed Brown PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 7280c7d97c5SJed Brown by setting data structures and options. 7290c7d97c5SJed Brown 7300c7d97c5SJed Brown Input Parameter: 73153cdbc3dSStefano Zampini + pc - the preconditioner context 7320c7d97c5SJed Brown 7330c7d97c5SJed Brown Application Interface Routine: PCSetUp() 7340c7d97c5SJed Brown 7350c7d97c5SJed Brown Notes: 7360c7d97c5SJed Brown The interface routine PCSetUp() is not usually called directly by 7370c7d97c5SJed Brown the user, but instead is called by PCApply() if necessary. 7380c7d97c5SJed Brown */ 73953cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc) 7400c7d97c5SJed Brown { 7410c7d97c5SJed Brown PetscErrorCode ierr; 7420c7d97c5SJed Brown PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 743674ae819SStefano Zampini MatStructure flag; 744674ae819SStefano Zampini PetscBool computeis,computetopography,computesolvers; 7450c7d97c5SJed Brown 7460c7d97c5SJed Brown PetscFunctionBegin; 747674ae819SStefano Zampini /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */ 7483b03a366Sstefano_zampini /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 7499c0446d6SStefano Zampini So, we set to pcnone the Neumann problem of pcis in order to avoid unneeded computation 7500c7d97c5SJed Brown Also, we decide to directly build the (same) Dirichlet problem */ 7510c7d97c5SJed Brown ierr = PetscOptionsSetValue("-is_localN_pc_type","none");CHKERRQ(ierr); 7520c7d97c5SJed Brown ierr = PetscOptionsSetValue("-is_localD_pc_type","none");CHKERRQ(ierr); 7533b03a366Sstefano_zampini /* Get stdout for dbg */ 754674ae819SStefano Zampini if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) { 755ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr); 756e269702eSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 757*2b510759SStefano Zampini if (pcbddc->current_level) { 758*2b510759SStefano Zampini ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 759*2b510759SStefano Zampini } 760e269702eSStefano Zampini } 761674ae819SStefano Zampini /* first attempt to split work */ 762674ae819SStefano Zampini if (pc->setupcalled) { 763674ae819SStefano Zampini computeis = PETSC_FALSE; 764674ae819SStefano Zampini ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr); 765674ae819SStefano Zampini if (flag == SAME_PRECONDITIONER) { 766674ae819SStefano Zampini computetopography = PETSC_FALSE; 767674ae819SStefano Zampini computesolvers = PETSC_FALSE; 768674ae819SStefano Zampini } else if (flag == SAME_NONZERO_PATTERN) { 769674ae819SStefano Zampini computetopography = PETSC_FALSE; 770674ae819SStefano Zampini computesolvers = PETSC_TRUE; 771674ae819SStefano Zampini } else { /* DIFFERENT_NONZERO_PATTERN */ 772674ae819SStefano Zampini computetopography = PETSC_TRUE; 773674ae819SStefano Zampini computesolvers = PETSC_TRUE; 774674ae819SStefano Zampini } 775674ae819SStefano Zampini } else { 776674ae819SStefano Zampini computeis = PETSC_TRUE; 777674ae819SStefano Zampini computetopography = PETSC_TRUE; 778674ae819SStefano Zampini computesolvers = PETSC_TRUE; 779674ae819SStefano Zampini } 780674ae819SStefano Zampini /* Set up all the "iterative substructuring" common block */ 781674ae819SStefano Zampini if (computeis) { 782674ae819SStefano Zampini ierr = PCISSetUp(pc);CHKERRQ(ierr); 783674ae819SStefano Zampini } 784674ae819SStefano Zampini /* Analyze interface and set up local constraint and change of basis matrices */ 785674ae819SStefano Zampini if (computetopography) { 786674ae819SStefano Zampini /* reset data */ 787674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 788674ae819SStefano Zampini ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 789674ae819SStefano Zampini ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 790674ae819SStefano Zampini } 791674ae819SStefano Zampini if (computesolvers) { 792674ae819SStefano Zampini /* reset data */ 793674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 794674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 79599cc7994SStefano Zampini /* Create coarse and local stuffs */ 79699cc7994SStefano Zampini ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 797674ae819SStefano Zampini ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 7980c7d97c5SJed Brown } 799*2b510759SStefano Zampini if (pcbddc->dbg_flag && pcbddc->current_level) { 800*2b510759SStefano Zampini ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 801*2b510759SStefano Zampini } 8020c7d97c5SJed Brown PetscFunctionReturn(0); 8030c7d97c5SJed Brown } 8040c7d97c5SJed Brown 8050c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 8060c7d97c5SJed Brown /* 8070c7d97c5SJed Brown PCApply_BDDC - Applies the BDDC preconditioner to a vector. 8080c7d97c5SJed Brown 8090c7d97c5SJed Brown Input Parameters: 8100c7d97c5SJed Brown . pc - the preconditioner context 8110c7d97c5SJed Brown . r - input vector (global) 8120c7d97c5SJed Brown 8130c7d97c5SJed Brown Output Parameter: 8140c7d97c5SJed Brown . z - output vector (global) 8150c7d97c5SJed Brown 8160c7d97c5SJed Brown Application Interface Routine: PCApply() 8170c7d97c5SJed Brown */ 8180c7d97c5SJed Brown #undef __FUNCT__ 8190c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC" 82053cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 8210c7d97c5SJed Brown { 8220c7d97c5SJed Brown PC_IS *pcis = (PC_IS*)(pc->data); 8230c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 8240c7d97c5SJed Brown PetscErrorCode ierr; 8253b03a366Sstefano_zampini const PetscScalar one = 1.0; 8263b03a366Sstefano_zampini const PetscScalar m_one = -1.0; 8272617d88aSStefano Zampini const PetscScalar zero = 0.0; 8280c7d97c5SJed Brown 8290c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN 8300c7d97c5SJed Brown NN interface preconditioner changed to BDDC 8318eeda7d8SStefano Zampini Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 8320c7d97c5SJed Brown 8330c7d97c5SJed Brown PetscFunctionBegin; 83415aaf578SStefano Zampini if (!pcbddc->use_exact_dirichlet || pcbddc->coarse_psi_B) { 8350c7d97c5SJed Brown /* First Dirichlet solve */ 8360c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8370c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 83853cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 8390c7d97c5SJed Brown /* 8400c7d97c5SJed Brown Assembling right hand side for BDDC operator 841674ae819SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 842674ae819SStefano Zampini - pcis->vec1_B the interface part of the global vector z 8430c7d97c5SJed Brown */ 8440c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 8450c7d97c5SJed Brown ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 8468eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 8470c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 8480c7d97c5SJed Brown ierr = VecCopy(r,z);CHKERRQ(ierr); 8490c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8500c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 851674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 852b76ba322SStefano Zampini } else { 8530bdf917eSStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 854b76ba322SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 855674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 856b76ba322SStefano Zampini } 857b76ba322SStefano Zampini 8582617d88aSStefano Zampini /* Apply interface preconditioner 8592617d88aSStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 8602617d88aSStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr); 8612617d88aSStefano Zampini 862674ae819SStefano Zampini /* Apply transpose of partition of unity operator */ 863674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 8640c7d97c5SJed Brown 8653b03a366Sstefano_zampini /* Second Dirichlet solve and assembling of output */ 8660c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8670c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8680c7d97c5SJed Brown ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 8698eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 87053cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr); 8710c7d97c5SJed Brown ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr); 8728eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 8730c7d97c5SJed Brown ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr); 8740c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8750c7d97c5SJed Brown ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8760c7d97c5SJed Brown PetscFunctionReturn(0); 8770c7d97c5SJed Brown } 878da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 879674ae819SStefano Zampini 880da1bb401SStefano Zampini #undef __FUNCT__ 881da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC" 882da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc) 883da1bb401SStefano Zampini { 884da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 885da1bb401SStefano Zampini PetscErrorCode ierr; 886da1bb401SStefano Zampini 887da1bb401SStefano Zampini PetscFunctionBegin; 888da1bb401SStefano Zampini /* free data created by PCIS */ 889da1bb401SStefano Zampini ierr = PCISDestroy(pc);CHKERRQ(ierr); 890674ae819SStefano Zampini /* free BDDC custom data */ 891674ae819SStefano Zampini ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 892674ae819SStefano Zampini /* destroy objects related to topography */ 893674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 894674ae819SStefano Zampini /* free allocated graph structure */ 895da1bb401SStefano Zampini ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 896674ae819SStefano Zampini /* free data for scaling operator */ 897674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 898674ae819SStefano Zampini /* free solvers stuff */ 899674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 90033bc96a4SStefano Zampini ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr); 90133bc96a4SStefano Zampini ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr); 9029881197aSStefano Zampini ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr); 903ac78edfcSStefano Zampini ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 90462a6ff1dSStefano Zampini /* free global vectors needed in presolve */ 90562a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 90662a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 9073425bc38SStefano Zampini /* remove functions */ 908674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 909bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 910*2b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 911*2b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 912bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 913bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 914bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 915bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 916bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 917bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 918bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 919bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 920bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 921bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 922674ae819SStefano Zampini /* Free the private data structure */ 923674ae819SStefano Zampini ierr = PetscFree(pc->data);CHKERRQ(ierr); 924da1bb401SStefano Zampini PetscFunctionReturn(0); 925da1bb401SStefano Zampini } 9263425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 9271e6b0712SBarry Smith 9283425bc38SStefano Zampini #undef __FUNCT__ 9293425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 9303425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 9313425bc38SStefano Zampini { 932674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 9333425bc38SStefano Zampini PC_IS* pcis; 9343425bc38SStefano Zampini PC_BDDC* pcbddc; 9353425bc38SStefano Zampini PetscErrorCode ierr; 9360c7d97c5SJed Brown 9373425bc38SStefano Zampini PetscFunctionBegin; 9383425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 9393425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 9403425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 9413425bc38SStefano Zampini 9423425bc38SStefano Zampini /* change of basis for physical rhs if needed 9433425bc38SStefano Zampini It also changes the rhs in case of dirichlet boundaries */ 9443308cffdSStefano Zampini ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 9453425bc38SStefano Zampini /* store vectors for computation of fetidp final solution */ 9463425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9473425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 948fb223d50SStefano Zampini /* scale rhs since it should be unassembled */ 949fb223d50SStefano Zampini /* TODO use counter scaling? (also below) */ 9503425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9513425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 952674ae819SStefano Zampini /* Apply partition of unity */ 9533425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 954674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 9558eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 9563425bc38SStefano Zampini /* compute partially subassembled Schur complement right-hand side */ 9573425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 9583425bc38SStefano Zampini ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 9593425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 9603425bc38SStefano Zampini ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 9613425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9623425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 963674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 9643425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9653425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9663425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 9673425bc38SStefano Zampini } 9683425bc38SStefano Zampini /* BDDC rhs */ 9693425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 9708eeda7d8SStefano Zampini if (pcbddc->switch_static) { 9713425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 9723425bc38SStefano Zampini } 9733425bc38SStefano Zampini /* apply BDDC */ 9743425bc38SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 9753425bc38SStefano Zampini /* Application of B_delta and assembling of rhs for fetidp fluxes */ 9763425bc38SStefano Zampini ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 9773425bc38SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 9783425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9793425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9803425bc38SStefano Zampini /* restore original rhs */ 9813425bc38SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 9823425bc38SStefano Zampini PetscFunctionReturn(0); 9833425bc38SStefano Zampini } 9841e6b0712SBarry Smith 9853425bc38SStefano Zampini #undef __FUNCT__ 9863425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 9873425bc38SStefano Zampini /*@ 9883425bc38SStefano Zampini PCBDDCMatFETIDPGetRHS - Get rhs for FETIDP linear system. 9893425bc38SStefano Zampini 9903425bc38SStefano Zampini Collective 9913425bc38SStefano Zampini 9923425bc38SStefano Zampini Input Parameters: 9933425bc38SStefano Zampini + fetidp_mat - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators 9943425bc38SStefano Zampini + standard_rhs - the rhs of your linear system 9953425bc38SStefano Zampini 9963425bc38SStefano Zampini Output Parameters: 9973425bc38SStefano Zampini + fetidp_flux_rhs - the rhs of the FETIDP linear system 9983425bc38SStefano Zampini 9993425bc38SStefano Zampini Level: developer 10003425bc38SStefano Zampini 10013425bc38SStefano Zampini Notes: 10023425bc38SStefano Zampini 10033425bc38SStefano Zampini .seealso: PCBDDC 10043425bc38SStefano Zampini @*/ 10053425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 10063425bc38SStefano Zampini { 1007674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10083425bc38SStefano Zampini PetscErrorCode ierr; 10093425bc38SStefano Zampini 10103425bc38SStefano Zampini PetscFunctionBegin; 10113425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10123425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 10133425bc38SStefano Zampini PetscFunctionReturn(0); 10143425bc38SStefano Zampini } 10153425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 10161e6b0712SBarry Smith 10173425bc38SStefano Zampini #undef __FUNCT__ 10183425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 10193425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 10203425bc38SStefano Zampini { 1021674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10223425bc38SStefano Zampini PC_IS* pcis; 10233425bc38SStefano Zampini PC_BDDC* pcbddc; 10243425bc38SStefano Zampini PetscErrorCode ierr; 10253425bc38SStefano Zampini 10263425bc38SStefano Zampini PetscFunctionBegin; 10273425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10283425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 10293425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 10303425bc38SStefano Zampini 10313425bc38SStefano Zampini /* apply B_delta^T */ 10323425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10333425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10343425bc38SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 10353425bc38SStefano Zampini /* compute rhs for BDDC application */ 10363425bc38SStefano Zampini ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 10378eeda7d8SStefano Zampini if (pcbddc->switch_static) { 10383425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10393425bc38SStefano Zampini } 10403425bc38SStefano Zampini /* apply BDDC */ 10413425bc38SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 10423425bc38SStefano Zampini /* put values into standard global vector */ 10433425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10443425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10458eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 10463425bc38SStefano Zampini /* compute values into the interior if solved for the partially subassembled Schur complement */ 10473425bc38SStefano Zampini ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 10483425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 10493425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 10503425bc38SStefano Zampini } 10513425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10523425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10533425bc38SStefano Zampini /* final change of basis if needed 10543425bc38SStefano Zampini Is also sums the dirichlet part removed during RHS assembling */ 10553308cffdSStefano Zampini ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 10563425bc38SStefano Zampini PetscFunctionReturn(0); 10573425bc38SStefano Zampini } 10581e6b0712SBarry Smith 10593425bc38SStefano Zampini #undef __FUNCT__ 10603425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 10613425bc38SStefano Zampini /*@ 10623425bc38SStefano Zampini PCBDDCMatFETIDPGetSolution - Get Solution for FETIDP linear system. 10633425bc38SStefano Zampini 10643425bc38SStefano Zampini Collective 10653425bc38SStefano Zampini 10663425bc38SStefano Zampini Input Parameters: 10673425bc38SStefano Zampini + fetidp_mat - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators 10683425bc38SStefano Zampini + fetidp_flux_sol - the solution of the FETIDP linear system 10693425bc38SStefano Zampini 10703425bc38SStefano Zampini Output Parameters: 10713425bc38SStefano Zampini + standard_sol - the solution on the global domain 10723425bc38SStefano Zampini 10733425bc38SStefano Zampini Level: developer 10743425bc38SStefano Zampini 10753425bc38SStefano Zampini Notes: 10763425bc38SStefano Zampini 10773425bc38SStefano Zampini .seealso: PCBDDC 10783425bc38SStefano Zampini @*/ 10793425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 10803425bc38SStefano Zampini { 1081674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 10823425bc38SStefano Zampini PetscErrorCode ierr; 10833425bc38SStefano Zampini 10843425bc38SStefano Zampini PetscFunctionBegin; 10853425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 10863425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 10873425bc38SStefano Zampini PetscFunctionReturn(0); 10883425bc38SStefano Zampini } 10893425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 10901e6b0712SBarry Smith 1091f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1092f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1093f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1094f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1095674ae819SStefano Zampini 10963425bc38SStefano Zampini #undef __FUNCT__ 10973425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 10983425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 10993425bc38SStefano Zampini { 1100674ae819SStefano Zampini 1101674ae819SStefano Zampini FETIDPMat_ctx fetidpmat_ctx; 11023425bc38SStefano Zampini Mat newmat; 1103674ae819SStefano Zampini FETIDPPC_ctx fetidppc_ctx; 11043425bc38SStefano Zampini PC newpc; 1105ce94432eSBarry Smith MPI_Comm comm; 11063425bc38SStefano Zampini PetscErrorCode ierr; 11073425bc38SStefano Zampini 11083425bc38SStefano Zampini PetscFunctionBegin; 1109ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 11103425bc38SStefano Zampini /* FETIDP linear matrix */ 11113425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 11123425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 11133425bc38SStefano Zampini ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 11143425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 11153425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 11163425bc38SStefano Zampini ierr = MatSetUp(newmat);CHKERRQ(ierr); 11173425bc38SStefano Zampini /* FETIDP preconditioner */ 11183425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 11193425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 11203425bc38SStefano Zampini ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 11213425bc38SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 11223425bc38SStefano Zampini ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 11233425bc38SStefano Zampini ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 11243425bc38SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 11253425bc38SStefano Zampini ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr); 11263425bc38SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 11273425bc38SStefano Zampini /* return pointers for objects created */ 11283425bc38SStefano Zampini *fetidp_mat=newmat; 11293425bc38SStefano Zampini *fetidp_pc=newpc; 11303425bc38SStefano Zampini PetscFunctionReturn(0); 11313425bc38SStefano Zampini } 11321e6b0712SBarry Smith 11333425bc38SStefano Zampini #undef __FUNCT__ 11343425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 11353425bc38SStefano Zampini /*@ 11363425bc38SStefano Zampini PCBDDCCreateFETIDPOperators - Create operators for FETIDP. 11373425bc38SStefano Zampini 11383425bc38SStefano Zampini Collective 11393425bc38SStefano Zampini 11403425bc38SStefano Zampini Input Parameters: 11413425bc38SStefano Zampini + pc - the BDDC preconditioning context (setup must be already called) 11423425bc38SStefano Zampini 11433425bc38SStefano Zampini Level: developer 11443425bc38SStefano Zampini 11453425bc38SStefano Zampini Notes: 11463425bc38SStefano Zampini 11473425bc38SStefano Zampini .seealso: PCBDDC 11483425bc38SStefano Zampini @*/ 11493425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 11503425bc38SStefano Zampini { 11513425bc38SStefano Zampini PetscErrorCode ierr; 11523425bc38SStefano Zampini 11533425bc38SStefano Zampini PetscFunctionBegin; 11543425bc38SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11553425bc38SStefano Zampini if (pc->setupcalled) { 1156516d51a7SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1157f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 11583425bc38SStefano Zampini PetscFunctionReturn(0); 11593425bc38SStefano Zampini } 11600c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 1161da1bb401SStefano Zampini /*MC 1162da1bb401SStefano Zampini PCBDDC - Balancing Domain Decomposition by Constraints. 11630c7d97c5SJed Brown 1164da1bb401SStefano Zampini Options Database Keys: 1165da1bb401SStefano Zampini . -pcbddc ??? - 1166da1bb401SStefano Zampini 1167da1bb401SStefano Zampini Level: intermediate 1168da1bb401SStefano Zampini 1169da1bb401SStefano Zampini Notes: The matrix used with this preconditioner must be of type MATIS 1170da1bb401SStefano Zampini 1171da1bb401SStefano Zampini Unlike more 'conventional' interface preconditioners, this iterates over ALL the 1172da1bb401SStefano Zampini degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers 1173da1bb401SStefano Zampini on the subdomains). 1174da1bb401SStefano Zampini 1175da1bb401SStefano Zampini Options for the coarse grid preconditioner can be set with - 1176da1bb401SStefano Zampini Options for the Dirichlet subproblem can be set with - 1177da1bb401SStefano Zampini Options for the Neumann subproblem can be set with - 1178da1bb401SStefano Zampini 1179da1bb401SStefano Zampini Contributed by Stefano Zampini 1180da1bb401SStefano Zampini 1181da1bb401SStefano Zampini .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1182da1bb401SStefano Zampini M*/ 1183b2573a8aSBarry Smith 1184da1bb401SStefano Zampini #undef __FUNCT__ 1185da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC" 11868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1187da1bb401SStefano Zampini { 1188da1bb401SStefano Zampini PetscErrorCode ierr; 1189da1bb401SStefano Zampini PC_BDDC *pcbddc; 1190da1bb401SStefano Zampini 1191da1bb401SStefano Zampini PetscFunctionBegin; 1192da1bb401SStefano Zampini /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1193da1bb401SStefano Zampini ierr = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr); 1194da1bb401SStefano Zampini pc->data = (void*)pcbddc; 1195da1bb401SStefano Zampini 1196da1bb401SStefano Zampini /* create PCIS data structure */ 1197da1bb401SStefano Zampini ierr = PCISCreate(pc);CHKERRQ(ierr); 1198da1bb401SStefano Zampini 119947d04d0dSStefano Zampini /* BDDC customization */ 120047d04d0dSStefano Zampini pcbddc->use_vertices = PETSC_TRUE; 120147d04d0dSStefano Zampini pcbddc->use_edges = PETSC_TRUE; 120247d04d0dSStefano Zampini pcbddc->use_faces = PETSC_FALSE; 120347d04d0dSStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 120447d04d0dSStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 120547d04d0dSStefano Zampini pcbddc->switch_static = PETSC_FALSE; 120647d04d0dSStefano Zampini pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 120747d04d0dSStefano Zampini pcbddc->dbg_flag = 0; 120847d04d0dSStefano Zampini 1209674ae819SStefano Zampini pcbddc->user_primal_vertices = 0; 12100bdf917eSStefano Zampini pcbddc->NullSpace = 0; 12113972b0daSStefano Zampini pcbddc->temp_solution = 0; 1212534831adSStefano Zampini pcbddc->original_rhs = 0; 1213534831adSStefano Zampini pcbddc->local_mat = 0; 1214534831adSStefano Zampini pcbddc->ChangeOfBasisMatrix = 0; 1215da1bb401SStefano Zampini pcbddc->coarse_vec = 0; 1216da1bb401SStefano Zampini pcbddc->coarse_rhs = 0; 1217da1bb401SStefano Zampini pcbddc->coarse_ksp = 0; 1218da1bb401SStefano Zampini pcbddc->coarse_phi_B = 0; 1219da1bb401SStefano Zampini pcbddc->coarse_phi_D = 0; 122015aaf578SStefano Zampini pcbddc->coarse_psi_B = 0; 122115aaf578SStefano Zampini pcbddc->coarse_psi_D = 0; 1222da1bb401SStefano Zampini pcbddc->vec1_P = 0; 1223da1bb401SStefano Zampini pcbddc->vec1_R = 0; 1224da1bb401SStefano Zampini pcbddc->vec2_R = 0; 1225da1bb401SStefano Zampini pcbddc->local_auxmat1 = 0; 1226da1bb401SStefano Zampini pcbddc->local_auxmat2 = 0; 1227da1bb401SStefano Zampini pcbddc->R_to_B = 0; 1228da1bb401SStefano Zampini pcbddc->R_to_D = 0; 1229da1bb401SStefano Zampini pcbddc->ksp_D = 0; 1230da1bb401SStefano Zampini pcbddc->ksp_R = 0; 1231da1bb401SStefano Zampini pcbddc->NeumannBoundaries = 0; 1232da1bb401SStefano Zampini pcbddc->ISForDofs = 0; 1233da1bb401SStefano Zampini pcbddc->ConstraintMatrix = 0; 1234b76ba322SStefano Zampini pcbddc->use_exact_dirichlet = PETSC_TRUE; 123547d04d0dSStefano Zampini pcbddc->coarse_loc_to_glob = 0; 123647d04d0dSStefano Zampini pcbddc->coarsening_ratio = 8; 12374fad6a16SStefano Zampini pcbddc->current_level = 0; 1238*2b510759SStefano Zampini pcbddc->max_levels = 0; 1239da1bb401SStefano Zampini 1240674ae819SStefano Zampini /* create local graph structure */ 1241674ae819SStefano Zampini ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1242674ae819SStefano Zampini 1243674ae819SStefano Zampini /* scaling */ 1244674ae819SStefano Zampini pcbddc->use_deluxe_scaling = PETSC_FALSE; 1245674ae819SStefano Zampini pcbddc->work_scaling = 0; 1246da1bb401SStefano Zampini 1247da1bb401SStefano Zampini /* function pointers */ 1248da1bb401SStefano Zampini pc->ops->apply = PCApply_BDDC; 1249da1bb401SStefano Zampini pc->ops->applytranspose = 0; 1250da1bb401SStefano Zampini pc->ops->setup = PCSetUp_BDDC; 1251da1bb401SStefano Zampini pc->ops->destroy = PCDestroy_BDDC; 1252da1bb401SStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1253da1bb401SStefano Zampini pc->ops->view = 0; 1254da1bb401SStefano Zampini pc->ops->applyrichardson = 0; 1255da1bb401SStefano Zampini pc->ops->applysymmetricleft = 0; 1256da1bb401SStefano Zampini pc->ops->applysymmetricright = 0; 1257534831adSStefano Zampini pc->ops->presolve = PCPreSolve_BDDC; 1258534831adSStefano Zampini pc->ops->postsolve = PCPostSolve_BDDC; 1259da1bb401SStefano Zampini 1260da1bb401SStefano Zampini /* composing function */ 1261674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1262bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 1263*2b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1264*2b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1265bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1266bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1267bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1268bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1269bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1270bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 1271bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1272bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1273bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1274bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1275da1bb401SStefano Zampini PetscFunctionReturn(0); 1276da1bb401SStefano Zampini } 12773425bc38SStefano Zampini 1278