153cdbc3dSStefano Zampini /* TODOLIST 2eb97c9d2SStefano Zampini 3eb97c9d2SStefano Zampini ConstraintsSetup 4eb97c9d2SStefano Zampini - tolerances for constraints as an option (take care of single precision!) 5eb97c9d2SStefano Zampini 6eb97c9d2SStefano Zampini Solvers 7a0d3c3abSStefano Zampini - Add support for cholesky for coarse solver (similar to local solvers) 8eb97c9d2SStefano Zampini - Propagate ksp prefixes for solvers to mat objects? 9eb97c9d2SStefano Zampini - Propagate nearnullspace info among levels 10eb97c9d2SStefano Zampini 11eb97c9d2SStefano Zampini User interface 12b9b85e73SStefano Zampini - ** DofSplitting and DM attached to pc? 13eb97c9d2SStefano Zampini 14eb97c9d2SStefano Zampini Debugging output 15b9b85e73SStefano Zampini - * Better management of verbosity levels of debugging output 16eb97c9d2SStefano Zampini 17eb97c9d2SStefano Zampini Build 18eb97c9d2SStefano Zampini - make runexe59 19eb97c9d2SStefano Zampini 20eb97c9d2SStefano Zampini Extra 21b9b85e73SStefano Zampini - ** GetRid of PCBDDCApplySchur, use MatSchur instead 22b9b85e73SStefano Zampini - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)? 23eb97c9d2SStefano Zampini - add support for computing h,H and related using coordinates? 24c0b83709SStefano Zampini - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap) 25eb97c9d2SStefano Zampini - Better management in PCIS code 26eb97c9d2SStefano Zampini - BDDC with MG framework? 27eb97c9d2SStefano Zampini 28eb97c9d2SStefano Zampini FETIDP 29eb97c9d2SStefano Zampini - Move FETIDP code to its own classes 30eb97c9d2SStefano Zampini 31eb97c9d2SStefano Zampini MATIS related operations contained in BDDC code 32eb97c9d2SStefano Zampini - Provide general case for subassembling 33b9b85e73SStefano Zampini - *** Preallocation routines in MatISGetMPIAXAIJ 34eb97c9d2SStefano Zampini 3553cdbc3dSStefano Zampini */ 360c7d97c5SJed Brown 3753cdbc3dSStefano Zampini /* ---------------------------------------------------------------------------------------------------------------------------------------------- 380c7d97c5SJed Brown Implementation of BDDC preconditioner based on: 390c7d97c5SJed Brown C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 4053cdbc3dSStefano Zampini ---------------------------------------------------------------------------------------------------------------------------------------------- */ 4153cdbc3dSStefano Zampini 42ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 43ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 443b03a366Sstefano_zampini #include <petscblaslapack.h> 45674ae819SStefano Zampini 460c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 470c7d97c5SJed Brown #undef __FUNCT__ 480c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC" 490c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc) 500c7d97c5SJed Brown { 510c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 520c7d97c5SJed Brown PetscErrorCode ierr; 530c7d97c5SJed Brown 540c7d97c5SJed Brown PetscFunctionBegin; 550c7d97c5SJed Brown ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 568eeda7d8SStefano Zampini /* Verbose debugging */ 578eeda7d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 588eeda7d8SStefano Zampini /* Primal space cumstomization */ 5908a5cf49SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);CHKERRQ(ierr); 608eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr); 618eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr); 628eeda7d8SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr); 63*6661aa1dSStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_true_nnsp","Use near null space attached to the matrix without modifications","none",pcbddc->use_nnsp_true,&pcbddc->use_nnsp_true,NULL);CHKERRQ(ierr); 648eeda7d8SStefano Zampini /* Change of basis */ 65b9b85e73SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not internal change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr); 66b9b85e73SStefano Zampini ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not internal change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr); 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); 7274e2c79eSStefano Zampini ierr = PetscOptionsInt("-pc_bddc_coarse_redistribute","Number of procs where to redistribute coarse problem","none",pcbddc->redistribute_coarse,&pcbddc->redistribute_coarse,NULL);CHKERRQ(ierr); 730298fd71SBarry Smith ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 742b510759SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 75323d395dSStefano 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); 76674ae819SStefano 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); 771cef56d8SStefano Zampini ierr = PetscOptionsInt("-pc_bddc_deluxe_threshold","Deluxe subproblems (Schur principal minors) smaller than this value are explicilty computed (-1 computes all)","none",pcbddc->deluxe_threshold,&pcbddc->deluxe_threshold,NULL);CHKERRQ(ierr); 7834a97f8cSStefano Zampini ierr = PetscOptionsBool("-pc_bddc_deluxe_rebuild","Whether or not the interface graph for deluxe has to be rebuilt (i.e. use a standard definition of the interface)","none",pcbddc->deluxe_rebuild,&pcbddc->deluxe_rebuild,NULL);CHKERRQ(ierr); 7934a97f8cSStefano Zampini ierr = PetscOptionsInt("-pc_bddc_deluxe_layers","Number of dofs' layers for the application of deluxe cheap version (i.e. -1 uses all dofs)","none",pcbddc->deluxe_layers,&pcbddc->deluxe_layers,NULL);CHKERRQ(ierr); 8034a97f8cSStefano Zampini ierr = PetscOptionsBool("-pc_bddc_deluxe_use_useradj","Whether or not the CSR graph specified by the user should be used for computing layers (default is to use adj of local mat)","none",pcbddc->deluxe_use_useradj,&pcbddc->deluxe_use_useradj,NULL);CHKERRQ(ierr); 810c7d97c5SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 820c7d97c5SJed Brown PetscFunctionReturn(0); 830c7d97c5SJed Brown } 840c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 85674ae819SStefano Zampini #undef __FUNCT__ 86906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC" 87906d46d4SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change) 88b9b85e73SStefano Zampini { 89b9b85e73SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 90b9b85e73SStefano Zampini PetscErrorCode ierr; 91b9b85e73SStefano Zampini 92b9b85e73SStefano Zampini PetscFunctionBegin; 93b9b85e73SStefano Zampini ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr); 94b9b85e73SStefano Zampini ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr); 95b9b85e73SStefano Zampini pcbddc->user_ChangeOfBasisMatrix = change; 96b9b85e73SStefano Zampini PetscFunctionReturn(0); 97b9b85e73SStefano Zampini } 98b9b85e73SStefano Zampini #undef __FUNCT__ 99906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat" 100b9b85e73SStefano Zampini /*@ 101906d46d4SStefano Zampini PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs 102b9b85e73SStefano Zampini 103b9b85e73SStefano Zampini Collective on PC 104b9b85e73SStefano Zampini 105b9b85e73SStefano Zampini Input Parameters: 106b9b85e73SStefano Zampini + pc - the preconditioning context 107906d46d4SStefano Zampini - change - the change of basis matrix 108b9b85e73SStefano Zampini 109b9b85e73SStefano Zampini Level: intermediate 110b9b85e73SStefano Zampini 111b9b85e73SStefano Zampini Notes: 112b9b85e73SStefano Zampini 113b9b85e73SStefano Zampini .seealso: PCBDDC 114b9b85e73SStefano Zampini @*/ 115906d46d4SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change) 116b9b85e73SStefano Zampini { 117b9b85e73SStefano Zampini PetscErrorCode ierr; 118b9b85e73SStefano Zampini 119b9b85e73SStefano Zampini PetscFunctionBegin; 120b9b85e73SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 121b9b85e73SStefano Zampini PetscValidHeaderSpecific(change,MAT_CLASSID,2); 122906d46d4SStefano Zampini PetscCheckSameComm(pc,1,change,2); 123906d46d4SStefano Zampini if (pc->mat) { 124906d46d4SStefano Zampini PetscInt rows_c,cols_c,rows,cols; 125906d46d4SStefano Zampini ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr); 126906d46d4SStefano Zampini ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr); 127906d46d4SStefano Zampini if (rows_c != rows) { 128906d46d4SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows); 129906d46d4SStefano Zampini } 130906d46d4SStefano Zampini if (cols_c != cols) { 131906d46d4SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols); 132906d46d4SStefano Zampini } 133906d46d4SStefano Zampini ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr); 134906d46d4SStefano Zampini ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr); 135906d46d4SStefano Zampini if (rows_c != rows) { 136906d46d4SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows); 137906d46d4SStefano Zampini } 138906d46d4SStefano Zampini if (cols_c != cols) { 139906d46d4SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols); 140906d46d4SStefano Zampini } 141906d46d4SStefano Zampini } 142906d46d4SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr); 143b9b85e73SStefano Zampini PetscFunctionReturn(0); 144b9b85e73SStefano Zampini } 145b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */ 146b9b85e73SStefano Zampini #undef __FUNCT__ 147674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 148674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 149674ae819SStefano Zampini { 150674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 151674ae819SStefano Zampini PetscErrorCode ierr; 1521e6b0712SBarry Smith 153674ae819SStefano Zampini PetscFunctionBegin; 154674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 155674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 156674ae819SStefano Zampini pcbddc->user_primal_vertices = PrimalVertices; 157674ae819SStefano Zampini PetscFunctionReturn(0); 158674ae819SStefano Zampini } 159674ae819SStefano Zampini #undef __FUNCT__ 160674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 161674ae819SStefano Zampini /*@ 16228509bceSStefano Zampini PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 163674ae819SStefano Zampini 164674ae819SStefano Zampini Not collective 165674ae819SStefano Zampini 166674ae819SStefano Zampini Input Parameters: 167674ae819SStefano Zampini + pc - the preconditioning context 16828509bceSStefano Zampini - PrimalVertices - index set of primal vertices in local numbering 169674ae819SStefano Zampini 170674ae819SStefano Zampini Level: intermediate 171674ae819SStefano Zampini 172674ae819SStefano Zampini Notes: 173674ae819SStefano Zampini 174674ae819SStefano Zampini .seealso: PCBDDC 175674ae819SStefano Zampini @*/ 176674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 177674ae819SStefano Zampini { 178674ae819SStefano Zampini PetscErrorCode ierr; 179674ae819SStefano Zampini 180674ae819SStefano Zampini PetscFunctionBegin; 181674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 182674ae819SStefano Zampini PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 183674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 184674ae819SStefano Zampini PetscFunctionReturn(0); 185674ae819SStefano Zampini } 186674ae819SStefano Zampini /* -------------------------------------------------------------------------- */ 1870c7d97c5SJed Brown #undef __FUNCT__ 1884fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 1894fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 1904fad6a16SStefano Zampini { 1914fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1924fad6a16SStefano Zampini 1934fad6a16SStefano Zampini PetscFunctionBegin; 1944fad6a16SStefano Zampini pcbddc->coarsening_ratio = k; 1954fad6a16SStefano Zampini PetscFunctionReturn(0); 1964fad6a16SStefano Zampini } 1971e6b0712SBarry Smith 1984fad6a16SStefano Zampini #undef __FUNCT__ 1994fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio" 2004fad6a16SStefano Zampini /*@ 20128509bceSStefano Zampini PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 2024fad6a16SStefano Zampini 2034fad6a16SStefano Zampini Logically collective on PC 2044fad6a16SStefano Zampini 2054fad6a16SStefano Zampini Input Parameters: 2064fad6a16SStefano Zampini + pc - the preconditioning context 20728509bceSStefano Zampini - k - coarsening ratio (H/h at the coarser level) 2084fad6a16SStefano Zampini 20928509bceSStefano Zampini Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 2104fad6a16SStefano Zampini 2114fad6a16SStefano Zampini Level: intermediate 2124fad6a16SStefano Zampini 2134fad6a16SStefano Zampini Notes: 2144fad6a16SStefano Zampini 2154fad6a16SStefano Zampini .seealso: PCBDDC 2164fad6a16SStefano Zampini @*/ 2174fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 2184fad6a16SStefano Zampini { 2194fad6a16SStefano Zampini PetscErrorCode ierr; 2204fad6a16SStefano Zampini 2214fad6a16SStefano Zampini PetscFunctionBegin; 2224fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2232b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,k,2); 2244fad6a16SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 2254fad6a16SStefano Zampini PetscFunctionReturn(0); 2264fad6a16SStefano Zampini } 2272b510759SStefano Zampini 228b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 2292b510759SStefano Zampini #undef __FUNCT__ 230b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 231b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 232b8ffe317SStefano Zampini { 233b8ffe317SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 234b8ffe317SStefano Zampini 235b8ffe317SStefano Zampini PetscFunctionBegin; 23685c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = flg; 237b8ffe317SStefano Zampini PetscFunctionReturn(0); 238b8ffe317SStefano Zampini } 239b8ffe317SStefano Zampini 240b8ffe317SStefano Zampini #undef __FUNCT__ 241b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 242b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 2432b510759SStefano Zampini { 2442b510759SStefano Zampini PetscErrorCode ierr; 2452b510759SStefano Zampini 2462b510759SStefano Zampini PetscFunctionBegin; 2472b510759SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 248b8ffe317SStefano Zampini PetscValidLogicalCollectiveBool(pc,flg,2); 249b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 2502b510759SStefano Zampini PetscFunctionReturn(0); 2512b510759SStefano Zampini } 2521e6b0712SBarry Smith 2534fad6a16SStefano Zampini #undef __FUNCT__ 2542b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC" 2552b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 2564fad6a16SStefano Zampini { 2574fad6a16SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2584fad6a16SStefano Zampini 2594fad6a16SStefano Zampini PetscFunctionBegin; 2602b510759SStefano Zampini pcbddc->current_level = level; 2614fad6a16SStefano Zampini PetscFunctionReturn(0); 2624fad6a16SStefano Zampini } 2631e6b0712SBarry Smith 2644fad6a16SStefano Zampini #undef __FUNCT__ 265b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel" 266b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 267b8ffe317SStefano Zampini { 268b8ffe317SStefano Zampini PetscErrorCode ierr; 269b8ffe317SStefano Zampini 270b8ffe317SStefano Zampini PetscFunctionBegin; 271b8ffe317SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 272b8ffe317SStefano Zampini PetscValidLogicalCollectiveInt(pc,level,2); 273b8ffe317SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 274b8ffe317SStefano Zampini PetscFunctionReturn(0); 275b8ffe317SStefano Zampini } 276b8ffe317SStefano Zampini 277b8ffe317SStefano Zampini #undef __FUNCT__ 2782b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC" 2792b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 2802b510759SStefano Zampini { 2812b510759SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2822b510759SStefano Zampini 2832b510759SStefano Zampini PetscFunctionBegin; 2842b510759SStefano Zampini pcbddc->max_levels = levels; 2852b510759SStefano Zampini PetscFunctionReturn(0); 2862b510759SStefano Zampini } 2872b510759SStefano Zampini 2882b510759SStefano Zampini #undef __FUNCT__ 2892b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels" 2904fad6a16SStefano Zampini /*@ 29128509bceSStefano Zampini PCBDDCSetLevels - Sets the maximum number of levels for multilevel 2924fad6a16SStefano Zampini 2934fad6a16SStefano Zampini Logically collective on PC 2944fad6a16SStefano Zampini 2954fad6a16SStefano Zampini Input Parameters: 2964fad6a16SStefano Zampini + pc - the preconditioning context 29728509bceSStefano Zampini - levels - the maximum number of levels (max 9) 2984fad6a16SStefano Zampini 29928509bceSStefano Zampini Default value is 0, i.e. traditional one-level BDDC 3004fad6a16SStefano Zampini 3014fad6a16SStefano Zampini Level: intermediate 3024fad6a16SStefano Zampini 3034fad6a16SStefano Zampini Notes: 3044fad6a16SStefano Zampini 3054fad6a16SStefano Zampini .seealso: PCBDDC 3064fad6a16SStefano Zampini @*/ 3072b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 3084fad6a16SStefano Zampini { 3094fad6a16SStefano Zampini PetscErrorCode ierr; 3104fad6a16SStefano Zampini 3114fad6a16SStefano Zampini PetscFunctionBegin; 3124fad6a16SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3132b510759SStefano Zampini PetscValidLogicalCollectiveInt(pc,levels,2); 314312be037SStefano Zampini if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n"); 3152b510759SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 3164fad6a16SStefano Zampini PetscFunctionReturn(0); 3174fad6a16SStefano Zampini } 3184fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */ 3191e6b0712SBarry Smith 3204fad6a16SStefano Zampini #undef __FUNCT__ 3210bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 3220bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 3230bdf917eSStefano Zampini { 3240bdf917eSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3250bdf917eSStefano Zampini PetscErrorCode ierr; 3260bdf917eSStefano Zampini 3270bdf917eSStefano Zampini PetscFunctionBegin; 3280bdf917eSStefano Zampini ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 3290bdf917eSStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 3300bdf917eSStefano Zampini pcbddc->NullSpace = NullSpace; 3310bdf917eSStefano Zampini PetscFunctionReturn(0); 3320bdf917eSStefano Zampini } 3331e6b0712SBarry Smith 3340bdf917eSStefano Zampini #undef __FUNCT__ 3350bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace" 3360bdf917eSStefano Zampini /*@ 33728509bceSStefano Zampini PCBDDCSetNullSpace - Set nullspace for BDDC operator 3380bdf917eSStefano Zampini 3390bdf917eSStefano Zampini Logically collective on PC and MatNullSpace 3400bdf917eSStefano Zampini 3410bdf917eSStefano Zampini Input Parameters: 3420bdf917eSStefano Zampini + pc - the preconditioning context 34328509bceSStefano Zampini - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 3440bdf917eSStefano Zampini 3450bdf917eSStefano Zampini Level: intermediate 3460bdf917eSStefano Zampini 3470bdf917eSStefano Zampini Notes: 3480bdf917eSStefano Zampini 3490bdf917eSStefano Zampini .seealso: PCBDDC 3500bdf917eSStefano Zampini @*/ 3510bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 3520bdf917eSStefano Zampini { 3530bdf917eSStefano Zampini PetscErrorCode ierr; 3540bdf917eSStefano Zampini 3550bdf917eSStefano Zampini PetscFunctionBegin; 3560bdf917eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 357674ae819SStefano Zampini PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 3582b510759SStefano Zampini PetscCheckSameComm(pc,1,NullSpace,2); 3590bdf917eSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 3600bdf917eSStefano Zampini PetscFunctionReturn(0); 3610bdf917eSStefano Zampini } 3620bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */ 3631e6b0712SBarry Smith 3640bdf917eSStefano Zampini #undef __FUNCT__ 3653b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 3663b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 3673b03a366Sstefano_zampini { 3683b03a366Sstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 3693b03a366Sstefano_zampini PetscErrorCode ierr; 3703b03a366Sstefano_zampini 3713b03a366Sstefano_zampini PetscFunctionBegin; 372785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 373785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 3743b03a366Sstefano_zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 37536e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 37636e030ebSStefano Zampini pcbddc->DirichletBoundaries = DirichletBoundaries; 377fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 3783b03a366Sstefano_zampini PetscFunctionReturn(0); 3793b03a366Sstefano_zampini } 3801e6b0712SBarry Smith 3813b03a366Sstefano_zampini #undef __FUNCT__ 3823b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 3833b03a366Sstefano_zampini /*@ 38428509bceSStefano Zampini PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 3853b03a366Sstefano_zampini 386785d1243SStefano Zampini Collective 3873b03a366Sstefano_zampini 3883b03a366Sstefano_zampini Input Parameters: 3893b03a366Sstefano_zampini + pc - the preconditioning context 390785d1243SStefano Zampini - DirichletBoundaries - parallel IS defining the Dirichlet boundaries 3913b03a366Sstefano_zampini 3923b03a366Sstefano_zampini Level: intermediate 3933b03a366Sstefano_zampini 394785d1243SStefano Zampini Notes: Any process can list any global node 3953b03a366Sstefano_zampini 3963b03a366Sstefano_zampini .seealso: PCBDDC 3973b03a366Sstefano_zampini @*/ 3983b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 3993b03a366Sstefano_zampini { 4003b03a366Sstefano_zampini PetscErrorCode ierr; 4013b03a366Sstefano_zampini 4023b03a366Sstefano_zampini PetscFunctionBegin; 4033b03a366Sstefano_zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 404674ae819SStefano Zampini PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 405785d1243SStefano Zampini PetscCheckSameComm(pc,1,DirichletBoundaries,2); 4063b03a366Sstefano_zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 4073b03a366Sstefano_zampini PetscFunctionReturn(0); 4083b03a366Sstefano_zampini } 4093b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */ 4101e6b0712SBarry Smith 4113b03a366Sstefano_zampini #undef __FUNCT__ 41282ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC" 41382ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries) 4140c7d97c5SJed Brown { 4150c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 4160c7d97c5SJed Brown PetscErrorCode ierr; 4170c7d97c5SJed Brown 4180c7d97c5SJed Brown PetscFunctionBegin; 419785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 420785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 4210c7d97c5SJed Brown ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 4220c7d97c5SJed Brown ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 423785d1243SStefano Zampini pcbddc->DirichletBoundariesLocal = DirichletBoundaries; 4247d24d155SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4250c7d97c5SJed Brown PetscFunctionReturn(0); 4260c7d97c5SJed Brown } 4270c7d97c5SJed Brown 4280c7d97c5SJed Brown #undef __FUNCT__ 42982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal" 4309c0446d6SStefano Zampini /*@ 43182ba6b80SStefano Zampini PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering. 4329c0446d6SStefano Zampini 433785d1243SStefano Zampini Collective 43453cdbc3dSStefano Zampini 43553cdbc3dSStefano Zampini Input Parameters: 43653cdbc3dSStefano Zampini + pc - the preconditioning context 43782ba6b80SStefano Zampini - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering) 43853cdbc3dSStefano Zampini 43953cdbc3dSStefano Zampini Level: intermediate 44053cdbc3dSStefano Zampini 4419c0446d6SStefano Zampini Notes: 44253cdbc3dSStefano Zampini 44353cdbc3dSStefano Zampini .seealso: PCBDDC 44453cdbc3dSStefano Zampini @*/ 44582ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries) 4460c7d97c5SJed Brown { 4470c7d97c5SJed Brown PetscErrorCode ierr; 4480c7d97c5SJed Brown 4490c7d97c5SJed Brown PetscFunctionBegin; 4500c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4510c7d97c5SJed Brown PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 45282ba6b80SStefano Zampini PetscCheckSameComm(pc,1,DirichletBoundaries,2); 45382ba6b80SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 4540c7d97c5SJed Brown PetscFunctionReturn(0); 4550c7d97c5SJed Brown } 4560c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 4570c7d97c5SJed Brown 4580c7d97c5SJed Brown #undef __FUNCT__ 4590c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 46053cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 4610c7d97c5SJed Brown { 4620c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 46353cdbc3dSStefano Zampini PetscErrorCode ierr; 4640c7d97c5SJed Brown 4650c7d97c5SJed Brown PetscFunctionBegin; 466785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 467785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 46853cdbc3dSStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 46936e030ebSStefano Zampini ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 47036e030ebSStefano Zampini pcbddc->NeumannBoundaries = NeumannBoundaries; 471fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 4720c7d97c5SJed Brown PetscFunctionReturn(0); 4730c7d97c5SJed Brown } 4741e6b0712SBarry Smith 4750c7d97c5SJed Brown #undef __FUNCT__ 4760c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 47757527edcSJed Brown /*@ 47828509bceSStefano Zampini PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 47957527edcSJed Brown 480785d1243SStefano Zampini Collective 48157527edcSJed Brown 48257527edcSJed Brown Input Parameters: 48357527edcSJed Brown + pc - the preconditioning context 484785d1243SStefano Zampini - NeumannBoundaries - parallel IS defining the Neumann boundaries 48557527edcSJed Brown 48657527edcSJed Brown Level: intermediate 48757527edcSJed Brown 488785d1243SStefano Zampini Notes: Any process can list any global node 48957527edcSJed Brown 49057527edcSJed Brown .seealso: PCBDDC 49157527edcSJed Brown @*/ 49253cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 4930c7d97c5SJed Brown { 4940c7d97c5SJed Brown PetscErrorCode ierr; 4950c7d97c5SJed Brown 4960c7d97c5SJed Brown PetscFunctionBegin; 4970c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 498674ae819SStefano Zampini PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 499785d1243SStefano Zampini PetscCheckSameComm(pc,1,NeumannBoundaries,2); 50053cdbc3dSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 50153cdbc3dSStefano Zampini PetscFunctionReturn(0); 50253cdbc3dSStefano Zampini } 50353cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 5041e6b0712SBarry Smith 50553cdbc3dSStefano Zampini #undef __FUNCT__ 50682ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC" 50782ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries) 5080c7d97c5SJed Brown { 5090c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 5100c7d97c5SJed Brown PetscErrorCode ierr; 5110c7d97c5SJed Brown 5120c7d97c5SJed Brown PetscFunctionBegin; 513785d1243SStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 514785d1243SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 5150c7d97c5SJed Brown ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 5160c7d97c5SJed Brown ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 517785d1243SStefano Zampini pcbddc->NeumannBoundariesLocal = NeumannBoundaries; 5187d24d155SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 5190c7d97c5SJed Brown PetscFunctionReturn(0); 5200c7d97c5SJed Brown } 5210c7d97c5SJed Brown 5220c7d97c5SJed Brown #undef __FUNCT__ 52382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal" 5240c7d97c5SJed Brown /*@ 52582ba6b80SStefano Zampini PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering. 5260c7d97c5SJed Brown 527785d1243SStefano Zampini Collective 5280c7d97c5SJed Brown 5290c7d97c5SJed Brown Input Parameters: 5300c7d97c5SJed Brown + pc - the preconditioning context 53182ba6b80SStefano Zampini - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering) 5320c7d97c5SJed Brown 5330c7d97c5SJed Brown Level: intermediate 5340c7d97c5SJed Brown 5350c7d97c5SJed Brown Notes: 5360c7d97c5SJed Brown 5370c7d97c5SJed Brown .seealso: PCBDDC 5380c7d97c5SJed Brown @*/ 53982ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries) 5400c7d97c5SJed Brown { 5410c7d97c5SJed Brown PetscErrorCode ierr; 5420c7d97c5SJed Brown 5430c7d97c5SJed Brown PetscFunctionBegin; 5440c7d97c5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5450c7d97c5SJed Brown PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 54682ba6b80SStefano Zampini PetscCheckSameComm(pc,1,NeumannBoundaries,2); 54782ba6b80SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 54853cdbc3dSStefano Zampini PetscFunctionReturn(0); 54953cdbc3dSStefano Zampini } 55053cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */ 55153cdbc3dSStefano Zampini 55253cdbc3dSStefano Zampini #undef __FUNCT__ 553da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 554da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 555da1bb401SStefano Zampini { 556da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 557da1bb401SStefano Zampini 558da1bb401SStefano Zampini PetscFunctionBegin; 559da1bb401SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundaries; 560da1bb401SStefano Zampini PetscFunctionReturn(0); 561da1bb401SStefano Zampini } 5621e6b0712SBarry Smith 563da1bb401SStefano Zampini #undef __FUNCT__ 564da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 565da1bb401SStefano Zampini /*@ 566785d1243SStefano Zampini PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries 567da1bb401SStefano Zampini 568785d1243SStefano Zampini Collective 569785d1243SStefano Zampini 570785d1243SStefano Zampini Input Parameters: 571785d1243SStefano Zampini . pc - the preconditioning context 572785d1243SStefano Zampini 573785d1243SStefano Zampini Output Parameters: 574785d1243SStefano Zampini . DirichletBoundaries - index set defining the Dirichlet boundaries 575785d1243SStefano Zampini 576785d1243SStefano Zampini Level: intermediate 577785d1243SStefano Zampini 578785d1243SStefano Zampini Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries 579785d1243SStefano Zampini 580785d1243SStefano Zampini .seealso: PCBDDC 581785d1243SStefano Zampini @*/ 582785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 583785d1243SStefano Zampini { 584785d1243SStefano Zampini PetscErrorCode ierr; 585785d1243SStefano Zampini 586785d1243SStefano Zampini PetscFunctionBegin; 587785d1243SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 588785d1243SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 589785d1243SStefano Zampini PetscFunctionReturn(0); 590785d1243SStefano Zampini } 591785d1243SStefano Zampini /* -------------------------------------------------------------------------- */ 592785d1243SStefano Zampini 593785d1243SStefano Zampini #undef __FUNCT__ 594785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC" 595785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries) 596785d1243SStefano Zampini { 597785d1243SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 598785d1243SStefano Zampini 599785d1243SStefano Zampini PetscFunctionBegin; 600785d1243SStefano Zampini *DirichletBoundaries = pcbddc->DirichletBoundariesLocal; 601785d1243SStefano Zampini PetscFunctionReturn(0); 602785d1243SStefano Zampini } 603785d1243SStefano Zampini 604785d1243SStefano Zampini #undef __FUNCT__ 60582ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal" 606da1bb401SStefano Zampini /*@ 60782ba6b80SStefano Zampini PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering) 608da1bb401SStefano Zampini 609785d1243SStefano Zampini Collective 610da1bb401SStefano Zampini 611da1bb401SStefano Zampini Input Parameters: 61228509bceSStefano Zampini . pc - the preconditioning context 613da1bb401SStefano Zampini 614da1bb401SStefano Zampini Output Parameters: 61528509bceSStefano Zampini . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 616da1bb401SStefano Zampini 617da1bb401SStefano Zampini Level: intermediate 618da1bb401SStefano Zampini 619da1bb401SStefano Zampini Notes: 620da1bb401SStefano Zampini 621da1bb401SStefano Zampini .seealso: PCBDDC 622da1bb401SStefano Zampini @*/ 62382ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries) 624da1bb401SStefano Zampini { 625da1bb401SStefano Zampini PetscErrorCode ierr; 626da1bb401SStefano Zampini 627da1bb401SStefano Zampini PetscFunctionBegin; 628da1bb401SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 62982ba6b80SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 630da1bb401SStefano Zampini PetscFunctionReturn(0); 631da1bb401SStefano Zampini } 632da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 6331e6b0712SBarry Smith 634da1bb401SStefano Zampini #undef __FUNCT__ 63553cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 63653cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 63753cdbc3dSStefano Zampini { 63853cdbc3dSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 63953cdbc3dSStefano Zampini 64053cdbc3dSStefano Zampini PetscFunctionBegin; 64153cdbc3dSStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundaries; 64253cdbc3dSStefano Zampini PetscFunctionReturn(0); 64353cdbc3dSStefano Zampini } 6441e6b0712SBarry Smith 64553cdbc3dSStefano Zampini #undef __FUNCT__ 64653cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 64753cdbc3dSStefano Zampini /*@ 648785d1243SStefano Zampini PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries 64953cdbc3dSStefano Zampini 650785d1243SStefano Zampini Collective 651785d1243SStefano Zampini 652785d1243SStefano Zampini Input Parameters: 653785d1243SStefano Zampini . pc - the preconditioning context 654785d1243SStefano Zampini 655785d1243SStefano Zampini Output Parameters: 656785d1243SStefano Zampini . NeumannBoundaries - index set defining the Neumann boundaries 657785d1243SStefano Zampini 658785d1243SStefano Zampini Level: intermediate 659785d1243SStefano Zampini 660785d1243SStefano Zampini Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries 661785d1243SStefano Zampini 662785d1243SStefano Zampini .seealso: PCBDDC 663785d1243SStefano Zampini @*/ 664785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 665785d1243SStefano Zampini { 666785d1243SStefano Zampini PetscErrorCode ierr; 667785d1243SStefano Zampini 668785d1243SStefano Zampini PetscFunctionBegin; 669785d1243SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 670785d1243SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 671785d1243SStefano Zampini PetscFunctionReturn(0); 672785d1243SStefano Zampini } 673785d1243SStefano Zampini /* -------------------------------------------------------------------------- */ 674785d1243SStefano Zampini 675785d1243SStefano Zampini #undef __FUNCT__ 676785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC" 677785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries) 678785d1243SStefano Zampini { 679785d1243SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 680785d1243SStefano Zampini 681785d1243SStefano Zampini PetscFunctionBegin; 682785d1243SStefano Zampini *NeumannBoundaries = pcbddc->NeumannBoundariesLocal; 683785d1243SStefano Zampini PetscFunctionReturn(0); 684785d1243SStefano Zampini } 685785d1243SStefano Zampini 686785d1243SStefano Zampini #undef __FUNCT__ 68782ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal" 68853cdbc3dSStefano Zampini /*@ 68982ba6b80SStefano Zampini PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering) 69053cdbc3dSStefano Zampini 691785d1243SStefano Zampini Collective 69253cdbc3dSStefano Zampini 69353cdbc3dSStefano Zampini Input Parameters: 69428509bceSStefano Zampini . pc - the preconditioning context 69553cdbc3dSStefano Zampini 69653cdbc3dSStefano Zampini Output Parameters: 69728509bceSStefano Zampini . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 69853cdbc3dSStefano Zampini 69953cdbc3dSStefano Zampini Level: intermediate 70053cdbc3dSStefano Zampini 70153cdbc3dSStefano Zampini Notes: 70253cdbc3dSStefano Zampini 70353cdbc3dSStefano Zampini .seealso: PCBDDC 70453cdbc3dSStefano Zampini @*/ 70582ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries) 70653cdbc3dSStefano Zampini { 70753cdbc3dSStefano Zampini PetscErrorCode ierr; 70853cdbc3dSStefano Zampini 70953cdbc3dSStefano Zampini PetscFunctionBegin; 71053cdbc3dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 71182ba6b80SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 7120c7d97c5SJed Brown PetscFunctionReturn(0); 7130c7d97c5SJed Brown } 71436e030ebSStefano Zampini /* -------------------------------------------------------------------------- */ 7151e6b0712SBarry Smith 71636e030ebSStefano Zampini #undef __FUNCT__ 717da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 7181a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 71936e030ebSStefano Zampini { 72036e030ebSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 721da1bb401SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 722da1bb401SStefano Zampini PetscErrorCode ierr; 72336e030ebSStefano Zampini 72436e030ebSStefano Zampini PetscFunctionBegin; 725674ae819SStefano Zampini /* free old CSR */ 726674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 727fb223d50SStefano Zampini /* TODO: PCBDDCGraphSetAdjacency */ 728674ae819SStefano Zampini /* get CSR into graph structure */ 729da1bb401SStefano Zampini if (copymode == PETSC_COPY_VALUES) { 730785e854fSJed Brown ierr = PetscMalloc1((nvtxs+1),&mat_graph->xadj);CHKERRQ(ierr); 731785e854fSJed Brown ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr); 732674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 733674ae819SStefano Zampini ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 734da1bb401SStefano Zampini } else if (copymode == PETSC_OWN_POINTER) { 7351a83f524SJed Brown mat_graph->xadj = (PetscInt*)xadj; 7361a83f524SJed Brown mat_graph->adjncy = (PetscInt*)adjncy; 737674ae819SStefano Zampini } 738575ad6abSStefano Zampini mat_graph->nvtxs_csr = nvtxs; 73936e030ebSStefano Zampini PetscFunctionReturn(0); 74036e030ebSStefano Zampini } 7411e6b0712SBarry Smith 74236e030ebSStefano Zampini #undef __FUNCT__ 743da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 74436e030ebSStefano Zampini /*@ 74528509bceSStefano Zampini PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix 74636e030ebSStefano Zampini 74736e030ebSStefano Zampini Not collective 74836e030ebSStefano Zampini 74936e030ebSStefano Zampini Input Parameters: 75036e030ebSStefano Zampini + pc - the preconditioning context 75128509bceSStefano Zampini . nvtxs - number of local vertices of the graph (i.e., the local size of your problem) 75228509bceSStefano Zampini . xadj, adjncy - the CSR graph 75328509bceSStefano Zampini - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 75436e030ebSStefano Zampini 75536e030ebSStefano Zampini Level: intermediate 75636e030ebSStefano Zampini 75736e030ebSStefano Zampini Notes: 75836e030ebSStefano Zampini 75928509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode 76036e030ebSStefano Zampini @*/ 7611a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 76236e030ebSStefano Zampini { 763575ad6abSStefano Zampini void (*f)(void) = 0; 76436e030ebSStefano Zampini PetscErrorCode ierr; 76536e030ebSStefano Zampini 76636e030ebSStefano Zampini PetscFunctionBegin; 76736e030ebSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 768674ae819SStefano Zampini PetscValidIntPointer(xadj,3); 769d7de1dd9SStefano Zampini PetscValidIntPointer(adjncy,4); 770674ae819SStefano Zampini if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 771674ae819SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 772da1bb401SStefano Zampini } 77336e030ebSStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 774575ad6abSStefano Zampini /* free arrays if PCBDDC is not the PC type */ 775575ad6abSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 776575ad6abSStefano Zampini if (!f && copymode == PETSC_OWN_POINTER) { 777575ad6abSStefano Zampini ierr = PetscFree(xadj);CHKERRQ(ierr); 778575ad6abSStefano Zampini ierr = PetscFree(adjncy);CHKERRQ(ierr); 77936e030ebSStefano Zampini } 78036e030ebSStefano Zampini PetscFunctionReturn(0); 78136e030ebSStefano Zampini } 7829c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */ 7831e6b0712SBarry Smith 7849c0446d6SStefano Zampini #undef __FUNCT__ 78563602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC" 78663602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 78763602bcaSStefano Zampini { 78863602bcaSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 78963602bcaSStefano Zampini PetscInt i; 79063602bcaSStefano Zampini PetscErrorCode ierr; 79163602bcaSStefano Zampini 79263602bcaSStefano Zampini PetscFunctionBegin; 79363602bcaSStefano Zampini /* Destroy ISes if they were already set */ 79463602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 79563602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 79663602bcaSStefano Zampini } 79763602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 79863602bcaSStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 79963602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 80063602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 80163602bcaSStefano Zampini } 80263602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 80363602bcaSStefano Zampini pcbddc->n_ISForDofs = 0; 80463602bcaSStefano Zampini /* allocate space then set */ 805d02579f5SStefano Zampini if (n_is) { 806d02579f5SStefano Zampini ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr); 807d02579f5SStefano Zampini } 80863602bcaSStefano Zampini for (i=0;i<n_is;i++) { 80963602bcaSStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 81063602bcaSStefano Zampini pcbddc->ISForDofsLocal[i]=ISForDofs[i]; 81163602bcaSStefano Zampini } 81263602bcaSStefano Zampini pcbddc->n_ISForDofsLocal=n_is; 81363602bcaSStefano Zampini if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 8141a8a16d1SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 81563602bcaSStefano Zampini PetscFunctionReturn(0); 81663602bcaSStefano Zampini } 81763602bcaSStefano Zampini 81863602bcaSStefano Zampini #undef __FUNCT__ 81963602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal" 82063602bcaSStefano Zampini /*@ 82163602bcaSStefano Zampini PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix 82263602bcaSStefano Zampini 82363602bcaSStefano Zampini Collective 82463602bcaSStefano Zampini 82563602bcaSStefano Zampini Input Parameters: 82663602bcaSStefano Zampini + pc - the preconditioning context 82763602bcaSStefano Zampini - n_is - number of index sets defining the fields 82863602bcaSStefano Zampini . ISForDofs - array of IS describing the fields in local ordering 82963602bcaSStefano Zampini 83063602bcaSStefano Zampini Level: intermediate 83163602bcaSStefano Zampini 83263602bcaSStefano 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. 83363602bcaSStefano Zampini 83463602bcaSStefano Zampini .seealso: PCBDDC 83563602bcaSStefano Zampini @*/ 83663602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[]) 83763602bcaSStefano Zampini { 83863602bcaSStefano Zampini PetscInt i; 83963602bcaSStefano Zampini PetscErrorCode ierr; 84063602bcaSStefano Zampini 84163602bcaSStefano Zampini PetscFunctionBegin; 84263602bcaSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 84363602bcaSStefano Zampini PetscValidLogicalCollectiveInt(pc,n_is,2); 84463602bcaSStefano Zampini for (i=0;i<n_is;i++) { 84563602bcaSStefano Zampini PetscCheckSameComm(pc,1,ISForDofs[i],3); 84663602bcaSStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 84763602bcaSStefano Zampini } 848e71e7a71SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 84963602bcaSStefano Zampini PetscFunctionReturn(0); 85063602bcaSStefano Zampini } 85163602bcaSStefano Zampini /* -------------------------------------------------------------------------- */ 85263602bcaSStefano Zampini 85363602bcaSStefano Zampini #undef __FUNCT__ 8549c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 8559c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 8569c0446d6SStefano Zampini { 8579c0446d6SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 8589c0446d6SStefano Zampini PetscInt i; 8599c0446d6SStefano Zampini PetscErrorCode ierr; 8609c0446d6SStefano Zampini 8619c0446d6SStefano Zampini PetscFunctionBegin; 862da1bb401SStefano Zampini /* Destroy ISes if they were already set */ 8639c0446d6SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 8649c0446d6SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 8659c0446d6SStefano Zampini } 866d11ae9bbSstefano_zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 86763602bcaSStefano Zampini /* last user setting takes precendence -> destroy any other customization */ 86863602bcaSStefano Zampini for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 86963602bcaSStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 87063602bcaSStefano Zampini } 87163602bcaSStefano Zampini ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 87263602bcaSStefano Zampini pcbddc->n_ISForDofsLocal = 0; 873da1bb401SStefano Zampini /* allocate space then set */ 874d02579f5SStefano Zampini if (n_is) { 875785e854fSJed Brown ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr); 876d02579f5SStefano Zampini } 8779c0446d6SStefano Zampini for (i=0;i<n_is;i++) { 878da1bb401SStefano Zampini ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 879da1bb401SStefano Zampini pcbddc->ISForDofs[i]=ISForDofs[i]; 8809c0446d6SStefano Zampini } 8819c0446d6SStefano Zampini pcbddc->n_ISForDofs=n_is; 88263602bcaSStefano Zampini if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 8831a8a16d1SStefano Zampini pcbddc->recompute_topography = PETSC_TRUE; 8849c0446d6SStefano Zampini PetscFunctionReturn(0); 8859c0446d6SStefano Zampini } 8861e6b0712SBarry Smith 8879c0446d6SStefano Zampini #undef __FUNCT__ 8889c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting" 8899c0446d6SStefano Zampini /*@ 89063602bcaSStefano Zampini PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix 8919c0446d6SStefano Zampini 89263602bcaSStefano Zampini Collective 8939c0446d6SStefano Zampini 8949c0446d6SStefano Zampini Input Parameters: 8959c0446d6SStefano Zampini + pc - the preconditioning context 89628509bceSStefano Zampini - n_is - number of index sets defining the fields 89763602bcaSStefano Zampini . ISForDofs - array of IS describing the fields in global ordering 8989c0446d6SStefano Zampini 8999c0446d6SStefano Zampini Level: intermediate 9009c0446d6SStefano Zampini 90163602bcaSStefano Zampini Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field. 9029c0446d6SStefano Zampini 9039c0446d6SStefano Zampini .seealso: PCBDDC 9049c0446d6SStefano Zampini @*/ 9059c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 9069c0446d6SStefano Zampini { 9072b510759SStefano Zampini PetscInt i; 9089c0446d6SStefano Zampini PetscErrorCode ierr; 9099c0446d6SStefano Zampini 9109c0446d6SStefano Zampini PetscFunctionBegin; 9119c0446d6SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 91263602bcaSStefano Zampini PetscValidLogicalCollectiveInt(pc,n_is,2); 9132b510759SStefano Zampini for (i=0;i<n_is;i++) { 91463602bcaSStefano Zampini PetscCheckSameComm(pc,1,ISForDofs[i],3); 91563602bcaSStefano Zampini PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 9162b510759SStefano Zampini } 9179c0446d6SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 9189c0446d6SStefano Zampini PetscFunctionReturn(0); 9199c0446d6SStefano Zampini } 920906d46d4SStefano Zampini 921da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 922534831adSStefano Zampini #undef __FUNCT__ 923534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC" 924534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 925534831adSStefano Zampini /* 926534831adSStefano Zampini PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 927534831adSStefano Zampini guess if a transformation of basis approach has been selected. 9289c0446d6SStefano Zampini 929534831adSStefano Zampini Input Parameter: 930534831adSStefano Zampini + pc - the preconditioner contex 931534831adSStefano Zampini 932534831adSStefano Zampini Application Interface Routine: PCPreSolve() 933534831adSStefano Zampini 934534831adSStefano Zampini Notes: 935534831adSStefano Zampini The interface routine PCPreSolve() is not usually called directly by 936534831adSStefano Zampini the user, but instead is called by KSPSolve(). 937534831adSStefano Zampini */ 938534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 939534831adSStefano Zampini { 940534831adSStefano Zampini PetscErrorCode ierr; 941534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 942534831adSStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 9433972b0daSStefano Zampini IS dirIS; 9443972b0daSStefano Zampini Vec used_vec; 9458d00608fSStefano Zampini PetscBool copy_rhs = PETSC_TRUE; 946534831adSStefano Zampini 947534831adSStefano Zampini PetscFunctionBegin; 94885c4d303SStefano Zampini /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 94985c4d303SStefano Zampini if (ksp) { 95085c4d303SStefano Zampini PetscBool iscg; 95185c4d303SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 95285c4d303SStefano Zampini if (!iscg) { 95385c4d303SStefano Zampini ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 95485c4d303SStefano Zampini } 95585c4d303SStefano Zampini } 95685c4d303SStefano Zampini /* Creates parallel work vectors used in presolve */ 95762a6ff1dSStefano Zampini if (!pcbddc->original_rhs) { 95862a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 95962a6ff1dSStefano Zampini } 96062a6ff1dSStefano Zampini if (!pcbddc->temp_solution) { 96162a6ff1dSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 96262a6ff1dSStefano Zampini } 9638d00608fSStefano Zampini 9643972b0daSStefano Zampini if (x) { 9653972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 9663972b0daSStefano Zampini used_vec = x; 9678d00608fSStefano Zampini } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */ 9683972b0daSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 9693972b0daSStefano Zampini used_vec = pcbddc->temp_solution; 9703972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 9713972b0daSStefano Zampini } 9728efcfb23SStefano Zampini 9738efcfb23SStefano Zampini /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */ 9743972b0daSStefano Zampini if (ksp) { 9758efcfb23SStefano Zampini ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr); 9768efcfb23SStefano Zampini if (!pcbddc->ksp_guess_nonzero) { 9773972b0daSStefano Zampini ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 9783972b0daSStefano Zampini } 9793972b0daSStefano Zampini } 9803308cffdSStefano Zampini 9818d00608fSStefano Zampini pcbddc->rhs_change = PETSC_FALSE; 9823972b0daSStefano Zampini 9833972b0daSStefano Zampini /* Take into account zeroed rows -> change rhs and store solution removed */ 984785d1243SStefano Zampini /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */ 98582ba6b80SStefano Zampini ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr); 98682ba6b80SStefano Zampini if (rhs && dirIS) { 987906d46d4SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 988785d1243SStefano Zampini PetscInt dirsize,i,*is_indices; 989785d1243SStefano Zampini PetscScalar *array_x,*array_diagonal; 990785d1243SStefano Zampini 9913972b0daSStefano Zampini ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 9923972b0daSStefano Zampini ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 9933972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9943972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9953972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9963972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 99782ba6b80SStefano Zampini ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr); 9983972b0daSStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 9993972b0daSStefano Zampini ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 10003972b0daSStefano Zampini ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 10012fa5cd67SKarl Rupp for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 10023972b0daSStefano Zampini ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 10033972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 10043972b0daSStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 10053972b0daSStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10063972b0daSStefano Zampini ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10078d00608fSStefano Zampini pcbddc->rhs_change = PETSC_TRUE; 10088efcfb23SStefano Zampini } 1009b76ba322SStefano Zampini 10108efcfb23SStefano Zampini /* remove the computed solution or the initial guess from the rhs */ 10118d00608fSStefano Zampini if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) { 10128d00608fSStefano Zampini /* store the original rhs */ 10138d00608fSStefano Zampini if (copy_rhs) { 10148d00608fSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 10158d00608fSStefano Zampini copy_rhs = PETSC_FALSE; 10168d00608fSStefano Zampini } 10178d00608fSStefano Zampini pcbddc->rhs_change = PETSC_TRUE; 10183972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 10193972b0daSStefano Zampini ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 10203972b0daSStefano Zampini ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 10218efcfb23SStefano Zampini ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 10223308cffdSStefano Zampini } 10238efcfb23SStefano Zampini ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 1024b76ba322SStefano Zampini 1025b76ba322SStefano Zampini /* store partially computed solution and set initial guess */ 10263972b0daSStefano Zampini if (x) { 10278efcfb23SStefano Zampini ierr = VecSet(x,0.0);CHKERRQ(ierr); 102885c4d303SStefano Zampini if (pcbddc->use_exact_dirichlet_trick) { 1029b76ba322SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1030b76ba322SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1031b76ba322SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 10328efcfb23SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10338efcfb23SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10348efcfb23SStefano Zampini if (ksp && !pcbddc->ksp_guess_nonzero) { 1035b76ba322SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 1036b76ba322SStefano Zampini } 1037b76ba322SStefano Zampini } 10383972b0daSStefano Zampini } 1039b76ba322SStefano Zampini 1040b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix) { 1041906d46d4SStefano Zampini PCBDDCChange_ctx change_ctx; 1042906d46d4SStefano Zampini 1043906d46d4SStefano Zampini /* get change ctx */ 1044906d46d4SStefano Zampini ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr); 1045906d46d4SStefano Zampini 1046906d46d4SStefano Zampini /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */ 1047906d46d4SStefano Zampini ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr); 1048906d46d4SStefano Zampini ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr); 1049906d46d4SStefano Zampini change_ctx->original_mat = pc->mat; 1050906d46d4SStefano Zampini 1051906d46d4SStefano Zampini /* change iteration matrix */ 1052906d46d4SStefano Zampini ierr = MatDestroy(&pc->mat);CHKERRQ(ierr); 1053906d46d4SStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr); 1054906d46d4SStefano Zampini pc->mat = pcbddc->new_global_mat; 1055906d46d4SStefano Zampini 10568d00608fSStefano Zampini /* store the original rhs */ 10578d00608fSStefano Zampini if (copy_rhs) { 10588d00608fSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 10598d00608fSStefano Zampini copy_rhs = PETSC_FALSE; 10608d00608fSStefano Zampini } 10618d00608fSStefano Zampini 1062906d46d4SStefano Zampini /* change rhs */ 1063906d46d4SStefano Zampini ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr); 1064906d46d4SStefano Zampini ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr); 10658d00608fSStefano Zampini pcbddc->rhs_change = PETSC_TRUE; 106692e3dcfbSStefano Zampini } 106792e3dcfbSStefano Zampini 106892e3dcfbSStefano Zampini /* remove nullspace if present */ 10698efcfb23SStefano Zampini if (ksp && x && pcbddc->NullSpace) { 10708efcfb23SStefano Zampini ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr); 10718d00608fSStefano Zampini /* store the original rhs */ 10728d00608fSStefano Zampini if (copy_rhs) { 10738d00608fSStefano Zampini ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 10748d00608fSStefano Zampini copy_rhs = PETSC_FALSE; 10758d00608fSStefano Zampini } 10768d00608fSStefano Zampini pcbddc->rhs_change = PETSC_TRUE; 1077d0195637SJed Brown ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 1078b76ba322SStefano Zampini } 1079534831adSStefano Zampini PetscFunctionReturn(0); 1080534831adSStefano Zampini } 1081906d46d4SStefano Zampini 1082534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 1083534831adSStefano Zampini #undef __FUNCT__ 1084534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC" 1085534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 1086534831adSStefano Zampini /* 1087534831adSStefano Zampini PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 1088534831adSStefano Zampini approach has been selected. Also, restores rhs to its original state. 1089534831adSStefano Zampini 1090534831adSStefano Zampini Input Parameter: 1091534831adSStefano Zampini + pc - the preconditioner contex 1092534831adSStefano Zampini 1093534831adSStefano Zampini Application Interface Routine: PCPostSolve() 1094534831adSStefano Zampini 1095534831adSStefano Zampini Notes: 1096534831adSStefano Zampini The interface routine PCPostSolve() is not usually called directly by 1097534831adSStefano Zampini the user, but instead is called by KSPSolve(). 1098534831adSStefano Zampini */ 1099534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1100534831adSStefano Zampini { 1101534831adSStefano Zampini PetscErrorCode ierr; 1102534831adSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1103534831adSStefano Zampini 1104534831adSStefano Zampini PetscFunctionBegin; 1105b9b85e73SStefano Zampini if (pcbddc->ChangeOfBasisMatrix) { 1106906d46d4SStefano Zampini PCBDDCChange_ctx change_ctx; 1107906d46d4SStefano Zampini 1108906d46d4SStefano Zampini /* get change ctx */ 1109906d46d4SStefano Zampini ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr); 1110906d46d4SStefano Zampini 1111906d46d4SStefano Zampini /* restore iteration matrix */ 1112906d46d4SStefano Zampini ierr = MatDestroy(&pc->mat);CHKERRQ(ierr); 1113906d46d4SStefano Zampini ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr); 1114906d46d4SStefano Zampini pc->mat = change_ctx->original_mat; 1115906d46d4SStefano Zampini 1116906d46d4SStefano Zampini /* get solution in original basis */ 1117906d46d4SStefano Zampini if (x) { 1118906d46d4SStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 1119906d46d4SStefano Zampini ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr); 1120906d46d4SStefano Zampini ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr); 11213425bc38SStefano Zampini } 1122534831adSStefano Zampini } 1123906d46d4SStefano Zampini 11243972b0daSStefano Zampini /* add solution removed in presolve */ 11256bcfc461SStefano Zampini if (x && pcbddc->rhs_change) { 11263425bc38SStefano Zampini ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 11273425bc38SStefano Zampini } 1128906d46d4SStefano Zampini 1129fb223d50SStefano Zampini /* restore rhs to its original state */ 11308d00608fSStefano Zampini if (rhs && pcbddc->rhs_change) { 1131fb223d50SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 1132fb223d50SStefano Zampini } 11338d00608fSStefano Zampini pcbddc->rhs_change = PETSC_FALSE; 11348efcfb23SStefano Zampini 11358efcfb23SStefano Zampini /* restore ksp guess state */ 11368efcfb23SStefano Zampini if (ksp) { 11378efcfb23SStefano Zampini ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr); 11388efcfb23SStefano Zampini } 1139534831adSStefano Zampini PetscFunctionReturn(0); 1140534831adSStefano Zampini } 1141534831adSStefano Zampini /* -------------------------------------------------------------------------- */ 114253cdbc3dSStefano Zampini #undef __FUNCT__ 114353cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC" 11440c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 11450c7d97c5SJed Brown /* 11460c7d97c5SJed Brown PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 11470c7d97c5SJed Brown by setting data structures and options. 11480c7d97c5SJed Brown 11490c7d97c5SJed Brown Input Parameter: 115053cdbc3dSStefano Zampini + pc - the preconditioner context 11510c7d97c5SJed Brown 11520c7d97c5SJed Brown Application Interface Routine: PCSetUp() 11530c7d97c5SJed Brown 11540c7d97c5SJed Brown Notes: 11550c7d97c5SJed Brown The interface routine PCSetUp() is not usually called directly by 11560c7d97c5SJed Brown the user, but instead is called by PCApply() if necessary. 11570c7d97c5SJed Brown */ 115853cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc) 11590c7d97c5SJed Brown { 11600c7d97c5SJed Brown PetscErrorCode ierr; 11610c7d97c5SJed Brown PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 1162f4ddd8eeSStefano Zampini MatNullSpace nearnullspace; 1163674ae819SStefano Zampini PetscBool computeis,computetopography,computesolvers; 1164165b64e2SStefano Zampini PetscBool new_nearnullspace_provided; 11650c7d97c5SJed Brown 11660c7d97c5SJed Brown PetscFunctionBegin; 1167f4ddd8eeSStefano Zampini /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 11683b03a366Sstefano_zampini /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1169fcd409f5SStefano Zampini Also, BDDC directly build the Dirichlet problem */ 1170f4ddd8eeSStefano Zampini /* split work */ 1171674ae819SStefano Zampini if (pc->setupcalled) { 1172674ae819SStefano Zampini computeis = PETSC_FALSE; 1173d1e9a80fSBarry Smith if (pc->flag == SAME_NONZERO_PATTERN) { 1174674ae819SStefano Zampini computetopography = PETSC_FALSE; 1175674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1176674ae819SStefano Zampini } else { /* DIFFERENT_NONZERO_PATTERN */ 1177674ae819SStefano Zampini computetopography = PETSC_TRUE; 1178674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1179674ae819SStefano Zampini } 1180674ae819SStefano Zampini } else { 1181674ae819SStefano Zampini computeis = PETSC_TRUE; 1182674ae819SStefano Zampini computetopography = PETSC_TRUE; 1183674ae819SStefano Zampini computesolvers = PETSC_TRUE; 1184674ae819SStefano Zampini } 1185fb180af4SStefano Zampini if (pcbddc->recompute_topography) { 1186fb180af4SStefano Zampini computetopography = PETSC_TRUE; 1187fb180af4SStefano Zampini } 1188f4ddd8eeSStefano Zampini 1189f4ddd8eeSStefano Zampini /* Get stdout for dbg */ 119070cf5478SStefano Zampini if (pcbddc->dbg_flag) { 119170cf5478SStefano Zampini if (!pcbddc->dbg_viewer) { 119258a03d70SStefano Zampini pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)); 1193f4ddd8eeSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 1194f4ddd8eeSStefano Zampini } 119558a03d70SStefano Zampini ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1196f4ddd8eeSStefano Zampini } 1197f4ddd8eeSStefano Zampini 1198fcd409f5SStefano Zampini /* Set up all the "iterative substructuring" common block without computing solvers */ 1199674ae819SStefano Zampini if (computeis) { 1200fcd409f5SStefano Zampini /* HACK INTO PCIS */ 1201fcd409f5SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 1202fcd409f5SStefano Zampini pcis->computesolvers = PETSC_FALSE; 1203674ae819SStefano Zampini ierr = PCISSetUp(pc);CHKERRQ(ierr); 120439e2fb2aSStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr); 1205674ae819SStefano Zampini } 1206f4ddd8eeSStefano Zampini 1207b9b85e73SStefano Zampini /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */ 1208906d46d4SStefano Zampini if (pcbddc->user_ChangeOfBasisMatrix) { 1209b9b85e73SStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 1210b9b85e73SStefano Zampini } 1211906d46d4SStefano Zampini 1212f4ddd8eeSStefano Zampini /* Analyze interface */ 1213674ae819SStefano Zampini if (computetopography) { 1214674ae819SStefano Zampini ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 121534a97f8cSStefano Zampini /* Schurs on subsets should be reset */ 121634a97f8cSStefano Zampini if (pcbddc->deluxe_ctx) { 121734a97f8cSStefano Zampini ierr = PCBDDCSubSchursReset(pcbddc->deluxe_ctx->sub_schurs);CHKERRQ(ierr); 121834a97f8cSStefano Zampini } 1219674ae819SStefano Zampini } 1220fb8d54d4SStefano Zampini 1221f4ddd8eeSStefano Zampini /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1222fb8d54d4SStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 1223f4ddd8eeSStefano Zampini ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr); 1224f4ddd8eeSStefano Zampini if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1225f4ddd8eeSStefano Zampini if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1226f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1227f4ddd8eeSStefano Zampini } else { 1228f4ddd8eeSStefano Zampini /* determine if the two nullspaces are different (should be lightweight) */ 1229f4ddd8eeSStefano Zampini if (nearnullspace != pcbddc->onearnullspace) { 1230f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1231165b64e2SStefano Zampini } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1232f4ddd8eeSStefano Zampini PetscInt i; 1233165b64e2SStefano Zampini const Vec *nearnullvecs; 1234165b64e2SStefano Zampini PetscObjectState state; 1235165b64e2SStefano Zampini PetscInt nnsp_size; 1236165b64e2SStefano Zampini ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 1237f4ddd8eeSStefano Zampini for (i=0;i<nnsp_size;i++) { 1238f4ddd8eeSStefano Zampini ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr); 1239165b64e2SStefano Zampini if (pcbddc->onearnullvecs_state[i] != state) { 1240f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1241f4ddd8eeSStefano Zampini break; 1242f4ddd8eeSStefano Zampini } 1243f4ddd8eeSStefano Zampini } 1244f4ddd8eeSStefano Zampini } 1245f4ddd8eeSStefano Zampini } 1246f4ddd8eeSStefano Zampini } else { 1247f4ddd8eeSStefano Zampini if (!nearnullspace) { /* both nearnullspaces are null */ 1248f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_FALSE; 1249f4ddd8eeSStefano Zampini } else { /* nearnullspace attached later */ 1250f4ddd8eeSStefano Zampini new_nearnullspace_provided = PETSC_TRUE; 1251f4ddd8eeSStefano Zampini } 1252f4ddd8eeSStefano Zampini } 1253f4ddd8eeSStefano Zampini 1254f4ddd8eeSStefano Zampini /* Setup constraints and related work vectors */ 1255727cdba6SStefano Zampini /* reset primal space flags */ 1256f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1257727cdba6SStefano Zampini pcbddc->new_primal_space_local = PETSC_FALSE; 1258fb8d54d4SStefano Zampini if (computetopography || new_nearnullspace_provided) { 1259727cdba6SStefano Zampini /* It also sets the primal space flags */ 1260674ae819SStefano Zampini ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 1261e7b262bdSStefano Zampini /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1262f4ddd8eeSStefano Zampini ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr); 126334a97f8cSStefano Zampini /* Schurs on subsets should be reset */ 126434a97f8cSStefano Zampini if (pcbddc->deluxe_ctx) { 126534a97f8cSStefano Zampini ierr = PCBDDCSubSchursReset(pcbddc->deluxe_ctx->sub_schurs);CHKERRQ(ierr); 126634a97f8cSStefano Zampini } 1267674ae819SStefano Zampini } 1268fb8d54d4SStefano Zampini 1269f4ddd8eeSStefano Zampini if (computesolvers || pcbddc->new_primal_space) { 127099cc7994SStefano Zampini /* Create coarse and local stuffs */ 127199cc7994SStefano Zampini ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 127234a97f8cSStefano Zampini /* Create scaling operators */ 1273674ae819SStefano Zampini ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 12740c7d97c5SJed Brown } 127570cf5478SStefano Zampini 127658a03d70SStefano Zampini if (pcbddc->dbg_flag) { 127758a03d70SStefano Zampini ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 12782b510759SStefano Zampini } 12790c7d97c5SJed Brown PetscFunctionReturn(0); 12800c7d97c5SJed Brown } 12810c7d97c5SJed Brown 12820c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 12830c7d97c5SJed Brown /* 128450efa1b5SStefano Zampini PCApply_BDDC - Applies the BDDC operator to a vector. 12850c7d97c5SJed Brown 12860c7d97c5SJed Brown Input Parameters: 12870c7d97c5SJed Brown . pc - the preconditioner context 12880c7d97c5SJed Brown . r - input vector (global) 12890c7d97c5SJed Brown 12900c7d97c5SJed Brown Output Parameter: 12910c7d97c5SJed Brown . z - output vector (global) 12920c7d97c5SJed Brown 12930c7d97c5SJed Brown Application Interface Routine: PCApply() 12940c7d97c5SJed Brown */ 12950c7d97c5SJed Brown #undef __FUNCT__ 12960c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC" 129753cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 12980c7d97c5SJed Brown { 12990c7d97c5SJed Brown PC_IS *pcis = (PC_IS*)(pc->data); 13000c7d97c5SJed Brown PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 13010c7d97c5SJed Brown PetscErrorCode ierr; 13023b03a366Sstefano_zampini const PetscScalar one = 1.0; 13033b03a366Sstefano_zampini const PetscScalar m_one = -1.0; 13042617d88aSStefano Zampini const PetscScalar zero = 0.0; 13050c7d97c5SJed Brown 13060c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN 13070c7d97c5SJed Brown NN interface preconditioner changed to BDDC 13088eeda7d8SStefano Zampini Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 13090c7d97c5SJed Brown 13100c7d97c5SJed Brown PetscFunctionBegin; 131185c4d303SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 13120c7d97c5SJed Brown /* First Dirichlet solve */ 13130c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 13140c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 131553cdbc3dSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 13160c7d97c5SJed Brown /* 13170c7d97c5SJed Brown Assembling right hand side for BDDC operator 1318674ae819SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 1319674ae819SStefano Zampini - pcis->vec1_B the interface part of the global vector z 13200c7d97c5SJed Brown */ 13210c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 13220c7d97c5SJed Brown ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 13238eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 13240c7d97c5SJed Brown ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 13250c7d97c5SJed Brown ierr = VecCopy(r,z);CHKERRQ(ierr); 13260c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 13270c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1328674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1329b76ba322SStefano Zampini } else { 13300bdf917eSStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 1331b76ba322SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 1332674ae819SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1333b76ba322SStefano Zampini } 1334b76ba322SStefano Zampini 13352617d88aSStefano Zampini /* Apply interface preconditioner 13362617d88aSStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1337dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr); 13382617d88aSStefano Zampini 1339674ae819SStefano Zampini /* Apply transpose of partition of unity operator */ 1340674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 13410c7d97c5SJed Brown 13423b03a366Sstefano_zampini /* Second Dirichlet solve and assembling of output */ 13430c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 13440c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 13450c7d97c5SJed Brown ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 13468eeda7d8SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1347df187020SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1348df187020SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1349df187020SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1350df187020SStefano Zampini ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 13510c7d97c5SJed Brown ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 13520c7d97c5SJed Brown ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 13530c7d97c5SJed Brown PetscFunctionReturn(0); 13540c7d97c5SJed Brown } 135550efa1b5SStefano Zampini 135650efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */ 135750efa1b5SStefano Zampini /* 135850efa1b5SStefano Zampini PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector. 135950efa1b5SStefano Zampini 136050efa1b5SStefano Zampini Input Parameters: 136150efa1b5SStefano Zampini . pc - the preconditioner context 136250efa1b5SStefano Zampini . r - input vector (global) 136350efa1b5SStefano Zampini 136450efa1b5SStefano Zampini Output Parameter: 136550efa1b5SStefano Zampini . z - output vector (global) 136650efa1b5SStefano Zampini 136750efa1b5SStefano Zampini Application Interface Routine: PCApplyTranspose() 136850efa1b5SStefano Zampini */ 136950efa1b5SStefano Zampini #undef __FUNCT__ 137050efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC" 137150efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z) 137250efa1b5SStefano Zampini { 137350efa1b5SStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 137450efa1b5SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 137550efa1b5SStefano Zampini PetscErrorCode ierr; 137650efa1b5SStefano Zampini const PetscScalar one = 1.0; 137750efa1b5SStefano Zampini const PetscScalar m_one = -1.0; 137850efa1b5SStefano Zampini const PetscScalar zero = 0.0; 137950efa1b5SStefano Zampini 138050efa1b5SStefano Zampini PetscFunctionBegin; 138150efa1b5SStefano Zampini if (!pcbddc->use_exact_dirichlet_trick) { 138250efa1b5SStefano Zampini /* First Dirichlet solve */ 138350efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 138450efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 138550efa1b5SStefano Zampini ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 138650efa1b5SStefano Zampini /* 138750efa1b5SStefano Zampini Assembling right hand side for BDDC operator 138850efa1b5SStefano Zampini - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 138950efa1b5SStefano Zampini - pcis->vec1_B the interface part of the global vector z 139050efa1b5SStefano Zampini */ 139150efa1b5SStefano Zampini ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 139250efa1b5SStefano Zampini ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 139350efa1b5SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 139450efa1b5SStefano Zampini ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 139550efa1b5SStefano Zampini ierr = VecCopy(r,z);CHKERRQ(ierr); 139650efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 139750efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 139850efa1b5SStefano Zampini ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 139950efa1b5SStefano Zampini } else { 140050efa1b5SStefano Zampini ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 140150efa1b5SStefano Zampini ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 140250efa1b5SStefano Zampini ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 140350efa1b5SStefano Zampini } 140450efa1b5SStefano Zampini 140550efa1b5SStefano Zampini /* Apply interface preconditioner 140650efa1b5SStefano Zampini input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1407dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr); 140850efa1b5SStefano Zampini 140950efa1b5SStefano Zampini /* Apply transpose of partition of unity operator */ 141050efa1b5SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 141150efa1b5SStefano Zampini 141250efa1b5SStefano Zampini /* Second Dirichlet solve and assembling of output */ 141350efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 141450efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 141550efa1b5SStefano Zampini ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 141650efa1b5SStefano Zampini if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1417b0147a47SStefano Zampini ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1418b0147a47SStefano Zampini ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1419b0147a47SStefano Zampini if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1420b0147a47SStefano Zampini ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 142150efa1b5SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 142250efa1b5SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 142350efa1b5SStefano Zampini PetscFunctionReturn(0); 142450efa1b5SStefano Zampini } 1425da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */ 1426674ae819SStefano Zampini 1427da1bb401SStefano Zampini #undef __FUNCT__ 1428da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC" 1429da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc) 1430da1bb401SStefano Zampini { 1431da1bb401SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1432da1bb401SStefano Zampini PetscErrorCode ierr; 1433da1bb401SStefano Zampini 1434da1bb401SStefano Zampini PetscFunctionBegin; 1435da1bb401SStefano Zampini /* free data created by PCIS */ 1436da1bb401SStefano Zampini ierr = PCISDestroy(pc);CHKERRQ(ierr); 1437674ae819SStefano Zampini /* free BDDC custom data */ 1438674ae819SStefano Zampini ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 1439674ae819SStefano Zampini /* destroy objects related to topography */ 1440674ae819SStefano Zampini ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 1441674ae819SStefano Zampini /* free allocated graph structure */ 1442da1bb401SStefano Zampini ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 144334a97f8cSStefano Zampini /* destroy objects for scaling operator */ 1444674ae819SStefano Zampini ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 144534a97f8cSStefano Zampini ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr); 1446674ae819SStefano Zampini /* free solvers stuff */ 1447674ae819SStefano Zampini ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 144862a6ff1dSStefano Zampini /* free global vectors needed in presolve */ 144962a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 145062a6ff1dSStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 1451906d46d4SStefano Zampini /* free stuff for change of basis hooks */ 1452906d46d4SStefano Zampini if (pcbddc->new_global_mat) { 1453906d46d4SStefano Zampini PCBDDCChange_ctx change_ctx; 1454906d46d4SStefano Zampini ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr); 1455906d46d4SStefano Zampini ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr); 1456906d46d4SStefano Zampini ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr); 1457906d46d4SStefano Zampini ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr); 1458906d46d4SStefano Zampini ierr = PetscFree(change_ctx);CHKERRQ(ierr); 1459906d46d4SStefano Zampini } 1460906d46d4SStefano Zampini ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr); 1461906d46d4SStefano Zampini /* remove map from local boundary to local numbering */ 146239e2fb2aSStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr); 14633425bc38SStefano Zampini /* remove functions */ 1464906d46d4SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr); 1465674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 1466bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 14672b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 1468b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 14692b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 1470bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 1471bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 147282ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1473bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 147482ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1475bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 147682ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1477bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1478785d1243SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1479bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 148063602bcaSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr); 1481bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 1482bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 1483bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 1484bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 1485674ae819SStefano Zampini /* Free the private data structure */ 1486674ae819SStefano Zampini ierr = PetscFree(pc->data);CHKERRQ(ierr); 1487da1bb401SStefano Zampini PetscFunctionReturn(0); 1488da1bb401SStefano Zampini } 14893425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 14901e6b0712SBarry Smith 14913425bc38SStefano Zampini #undef __FUNCT__ 14923425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 14933425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 14943425bc38SStefano Zampini { 1495674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 14963425bc38SStefano Zampini PC_IS* pcis; 14973425bc38SStefano Zampini PC_BDDC* pcbddc; 14983425bc38SStefano Zampini PetscErrorCode ierr; 14990c7d97c5SJed Brown 15003425bc38SStefano Zampini PetscFunctionBegin; 15013425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15023425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 15033425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 15043425bc38SStefano Zampini 15053425bc38SStefano Zampini /* change of basis for physical rhs if needed 15063425bc38SStefano Zampini It also changes the rhs in case of dirichlet boundaries */ 15073308cffdSStefano Zampini ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 15083425bc38SStefano Zampini /* store vectors for computation of fetidp final solution */ 15093425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 15103425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1511fb223d50SStefano Zampini /* scale rhs since it should be unassembled */ 1512fb223d50SStefano Zampini /* TODO use counter scaling? (also below) */ 15133425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 15143425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1515674ae819SStefano Zampini /* Apply partition of unity */ 15163425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1517674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 15188eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 15193425bc38SStefano Zampini /* compute partially subassembled Schur complement right-hand side */ 15203425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 15213425bc38SStefano Zampini ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 15223425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 15233425bc38SStefano Zampini ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 15243425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15253425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1526674ae819SStefano Zampini /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 15273425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 15283425bc38SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 15293425bc38SStefano Zampini ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 15303425bc38SStefano Zampini } 15313425bc38SStefano Zampini /* BDDC rhs */ 15323425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 15338eeda7d8SStefano Zampini if (pcbddc->switch_static) { 15343425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 15353425bc38SStefano Zampini } 15363425bc38SStefano Zampini /* apply BDDC */ 1537dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 15383425bc38SStefano Zampini /* Application of B_delta and assembling of rhs for fetidp fluxes */ 15393425bc38SStefano Zampini ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 15403425bc38SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 15413425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 15423425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 15433425bc38SStefano Zampini /* restore original rhs */ 15443425bc38SStefano Zampini ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 15453425bc38SStefano Zampini PetscFunctionReturn(0); 15463425bc38SStefano Zampini } 15471e6b0712SBarry Smith 15483425bc38SStefano Zampini #undef __FUNCT__ 15493425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 15503425bc38SStefano Zampini /*@ 155128509bceSStefano Zampini PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system 15523425bc38SStefano Zampini 15533425bc38SStefano Zampini Collective 15543425bc38SStefano Zampini 15553425bc38SStefano Zampini Input Parameters: 155628509bceSStefano Zampini + fetidp_mat - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators 155728509bceSStefano Zampini . standard_rhs - the right-hand side for your linear system 15583425bc38SStefano Zampini 15593425bc38SStefano Zampini Output Parameters: 156028509bceSStefano Zampini - fetidp_flux_rhs - the right-hand side for the FETIDP linear system 15613425bc38SStefano Zampini 15623425bc38SStefano Zampini Level: developer 15633425bc38SStefano Zampini 15643425bc38SStefano Zampini Notes: 15653425bc38SStefano Zampini 156628509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 15673425bc38SStefano Zampini @*/ 15683425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 15693425bc38SStefano Zampini { 1570674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 15713425bc38SStefano Zampini PetscErrorCode ierr; 15723425bc38SStefano Zampini 15733425bc38SStefano Zampini PetscFunctionBegin; 15743425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15753425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 15763425bc38SStefano Zampini PetscFunctionReturn(0); 15773425bc38SStefano Zampini } 15783425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 15791e6b0712SBarry Smith 15803425bc38SStefano Zampini #undef __FUNCT__ 15813425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 15823425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 15833425bc38SStefano Zampini { 1584674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 15853425bc38SStefano Zampini PC_IS* pcis; 15863425bc38SStefano Zampini PC_BDDC* pcbddc; 15873425bc38SStefano Zampini PetscErrorCode ierr; 15883425bc38SStefano Zampini 15893425bc38SStefano Zampini PetscFunctionBegin; 15903425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 15913425bc38SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 15923425bc38SStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 15933425bc38SStefano Zampini 15943425bc38SStefano Zampini /* apply B_delta^T */ 15953425bc38SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15963425bc38SStefano Zampini ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15973425bc38SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 15983425bc38SStefano Zampini /* compute rhs for BDDC application */ 15993425bc38SStefano Zampini ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 16008eeda7d8SStefano Zampini if (pcbddc->switch_static) { 16013425bc38SStefano Zampini ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 16023425bc38SStefano Zampini } 16033425bc38SStefano Zampini /* apply BDDC */ 1604dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 16053425bc38SStefano Zampini /* put values into standard global vector */ 16063425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 16073425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 16088eeda7d8SStefano Zampini if (!pcbddc->switch_static) { 16093425bc38SStefano Zampini /* compute values into the interior if solved for the partially subassembled Schur complement */ 16103425bc38SStefano Zampini ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 16113425bc38SStefano Zampini ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 16123425bc38SStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 16133425bc38SStefano Zampini } 16143425bc38SStefano Zampini ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 16153425bc38SStefano Zampini ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 16163425bc38SStefano Zampini /* final change of basis if needed 16173425bc38SStefano Zampini Is also sums the dirichlet part removed during RHS assembling */ 16183308cffdSStefano Zampini ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 16193425bc38SStefano Zampini PetscFunctionReturn(0); 16203425bc38SStefano Zampini } 16211e6b0712SBarry Smith 16223425bc38SStefano Zampini #undef __FUNCT__ 16233425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 16243425bc38SStefano Zampini /*@ 162528509bceSStefano Zampini PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system 16263425bc38SStefano Zampini 16273425bc38SStefano Zampini Collective 16283425bc38SStefano Zampini 16293425bc38SStefano Zampini Input Parameters: 163028509bceSStefano Zampini + fetidp_mat - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators 163128509bceSStefano Zampini . fetidp_flux_sol - the solution of the FETIDP linear system 16323425bc38SStefano Zampini 16333425bc38SStefano Zampini Output Parameters: 163428509bceSStefano Zampini - standard_sol - the solution defined on the physical domain 16353425bc38SStefano Zampini 16363425bc38SStefano Zampini Level: developer 16373425bc38SStefano Zampini 16383425bc38SStefano Zampini Notes: 16393425bc38SStefano Zampini 164028509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 16413425bc38SStefano Zampini @*/ 16423425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 16433425bc38SStefano Zampini { 1644674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 16453425bc38SStefano Zampini PetscErrorCode ierr; 16463425bc38SStefano Zampini 16473425bc38SStefano Zampini PetscFunctionBegin; 16483425bc38SStefano Zampini ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 16493425bc38SStefano Zampini ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 16503425bc38SStefano Zampini PetscFunctionReturn(0); 16513425bc38SStefano Zampini } 16523425bc38SStefano Zampini /* -------------------------------------------------------------------------- */ 16531e6b0712SBarry Smith 1654f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1655edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec); 1656f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1657f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1658edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec); 1659f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1660674ae819SStefano Zampini 16613425bc38SStefano Zampini #undef __FUNCT__ 16623425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 16633425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 16643425bc38SStefano Zampini { 1665674ae819SStefano Zampini 1666674ae819SStefano Zampini FETIDPMat_ctx fetidpmat_ctx; 16673425bc38SStefano Zampini Mat newmat; 1668674ae819SStefano Zampini FETIDPPC_ctx fetidppc_ctx; 16693425bc38SStefano Zampini PC newpc; 1670ce94432eSBarry Smith MPI_Comm comm; 16713425bc38SStefano Zampini PetscErrorCode ierr; 16723425bc38SStefano Zampini 16733425bc38SStefano Zampini PetscFunctionBegin; 1674ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 16753425bc38SStefano Zampini /* FETIDP linear matrix */ 16763425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 16773425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 16783425bc38SStefano Zampini ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 16793425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 1680edf7251bSStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr); 16813425bc38SStefano Zampini ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 16823425bc38SStefano Zampini ierr = MatSetUp(newmat);CHKERRQ(ierr); 16833425bc38SStefano Zampini /* FETIDP preconditioner */ 16843425bc38SStefano Zampini ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 16853425bc38SStefano Zampini ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 16863425bc38SStefano Zampini ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 16873425bc38SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 16883425bc38SStefano Zampini ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 16893425bc38SStefano Zampini ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 1690edf7251bSStefano Zampini ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr); 16913425bc38SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 169223ee1639SBarry Smith ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr); 16933425bc38SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 16943425bc38SStefano Zampini /* return pointers for objects created */ 16953425bc38SStefano Zampini *fetidp_mat=newmat; 16963425bc38SStefano Zampini *fetidp_pc=newpc; 16973425bc38SStefano Zampini PetscFunctionReturn(0); 16983425bc38SStefano Zampini } 16991e6b0712SBarry Smith 17003425bc38SStefano Zampini #undef __FUNCT__ 17013425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 17023425bc38SStefano Zampini /*@ 170328509bceSStefano Zampini PCBDDCCreateFETIDPOperators - Create operators for FETIDP 17043425bc38SStefano Zampini 17053425bc38SStefano Zampini Collective 17063425bc38SStefano Zampini 17073425bc38SStefano Zampini Input Parameters: 170828509bceSStefano Zampini + pc - the BDDC preconditioning context already setup 170928509bceSStefano Zampini 171028509bceSStefano Zampini Output Parameters: 171128509bceSStefano Zampini . fetidp_mat - shell FETIDP matrix object 171228509bceSStefano Zampini . fetidp_pc - shell Dirichlet preconditioner for FETIDP matrix 171328509bceSStefano Zampini 171428509bceSStefano Zampini Options Database Keys: 171528509bceSStefano Zampini - -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers 17163425bc38SStefano Zampini 17173425bc38SStefano Zampini Level: developer 17183425bc38SStefano Zampini 17193425bc38SStefano Zampini Notes: 172028509bceSStefano Zampini Currently the only operation provided for FETIDP matrix is MatMult 17213425bc38SStefano Zampini 17223425bc38SStefano Zampini .seealso: PCBDDC 17233425bc38SStefano Zampini @*/ 17243425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 17253425bc38SStefano Zampini { 17263425bc38SStefano Zampini PetscErrorCode ierr; 17273425bc38SStefano Zampini 17283425bc38SStefano Zampini PetscFunctionBegin; 17293425bc38SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 17303425bc38SStefano Zampini if (pc->setupcalled) { 1731516d51a7SStefano Zampini ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1732f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 17333425bc38SStefano Zampini PetscFunctionReturn(0); 17343425bc38SStefano Zampini } 17350c7d97c5SJed Brown /* -------------------------------------------------------------------------- */ 1736da1bb401SStefano Zampini /*MC 1737da1bb401SStefano Zampini PCBDDC - Balancing Domain Decomposition by Constraints. 17380c7d97c5SJed Brown 173928509bceSStefano Zampini An implementation of the BDDC preconditioner based on 174028509bceSStefano Zampini 174128509bceSStefano Zampini .vb 174228509bceSStefano Zampini [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 174328509bceSStefano 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 174428509bceSStefano Zampini [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 174528509bceSStefano Zampini .ve 174628509bceSStefano Zampini 174728509bceSStefano Zampini The matrix to be preconditioned (Pmat) must be of type MATIS. 174828509bceSStefano Zampini 1749b6fdb6dfSStefano Zampini Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 175028509bceSStefano Zampini 175128509bceSStefano Zampini It also works with unsymmetric and indefinite problems. 175228509bceSStefano Zampini 1753b6fdb6dfSStefano 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. 1754b6fdb6dfSStefano Zampini 175528509bceSStefano Zampini Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace 175628509bceSStefano Zampini 17576a818285SBarry Smith Boundary nodes are split in vertices, edges and faces using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph() 175828509bceSStefano Zampini 17596a818285SBarry Smith Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). 176028509bceSStefano Zampini 1761b6fdb6dfSStefano 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. 176228509bceSStefano Zampini 176328509bceSStefano Zampini The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object. 176428509bceSStefano Zampini 1765da1bb401SStefano Zampini Options Database Keys: 176628509bceSStefano Zampini 176728509bceSStefano Zampini . -pc_bddc_use_vertices <1> - use or not vertices in primal space 176828509bceSStefano Zampini . -pc_bddc_use_edges <1> - use or not edges in primal space 1769b6fdb6dfSStefano Zampini . -pc_bddc_use_faces <0> - use or not faces in primal space 177028509bceSStefano Zampini . -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only) 177128509bceSStefano Zampini . -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested 177228509bceSStefano Zampini . -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1]) 177328509bceSStefano Zampini . -pc_bddc_levels <0> - maximum number of levels for multilevel 177428509bceSStefano Zampini . -pc_bddc_coarsening_ratio - H/h ratio at the coarser level 177528509bceSStefano Zampini - -pc_bddc_check_level <0> - set verbosity level of debugging output 177628509bceSStefano Zampini 177728509bceSStefano Zampini Options for Dirichlet, Neumann or coarse solver can be set with 177828509bceSStefano Zampini .vb 177928509bceSStefano Zampini -pc_bddc_dirichlet_ 178028509bceSStefano Zampini -pc_bddc_neumann_ 178128509bceSStefano Zampini -pc_bddc_coarse_ 178228509bceSStefano Zampini .ve 178328509bceSStefano Zampini e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg 178428509bceSStefano Zampini 178528509bceSStefano Zampini When using a multilevel approach, solvers' options at the N-th level can be specified as 178628509bceSStefano Zampini .vb 1787312be037SStefano Zampini -pc_bddc_dirichlet_lN_ 1788312be037SStefano Zampini -pc_bddc_neumann_lN_ 1789312be037SStefano Zampini -pc_bddc_coarse_lN_ 179028509bceSStefano Zampini .ve 1791312be037SStefano Zampini Note that level number ranges from the finest 0 to the coarsest N. 1792da1bb401SStefano Zampini 1793da1bb401SStefano Zampini Level: intermediate 1794da1bb401SStefano Zampini 1795b6fdb6dfSStefano Zampini Developer notes: 1796da1bb401SStefano Zampini 179728509bceSStefano Zampini New deluxe scaling operator will be available soon. 1798da1bb401SStefano Zampini 1799da1bb401SStefano Zampini Contributed by Stefano Zampini 1800da1bb401SStefano Zampini 1801da1bb401SStefano Zampini .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1802da1bb401SStefano Zampini M*/ 1803b2573a8aSBarry Smith 1804da1bb401SStefano Zampini #undef __FUNCT__ 1805da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC" 18068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1807da1bb401SStefano Zampini { 1808da1bb401SStefano Zampini PetscErrorCode ierr; 1809da1bb401SStefano Zampini PC_BDDC *pcbddc; 1810da1bb401SStefano Zampini 1811da1bb401SStefano Zampini PetscFunctionBegin; 1812da1bb401SStefano Zampini /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1813b00a9115SJed Brown ierr = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr); 1814da1bb401SStefano Zampini pc->data = (void*)pcbddc; 1815da1bb401SStefano Zampini 1816da1bb401SStefano Zampini /* create PCIS data structure */ 1817da1bb401SStefano Zampini ierr = PCISCreate(pc);CHKERRQ(ierr); 1818da1bb401SStefano Zampini 181947d04d0dSStefano Zampini /* BDDC customization */ 182008a5cf49SStefano Zampini pcbddc->use_local_adj = PETSC_TRUE; 182147d04d0dSStefano Zampini pcbddc->use_vertices = PETSC_TRUE; 182247d04d0dSStefano Zampini pcbddc->use_edges = PETSC_TRUE; 182347d04d0dSStefano Zampini pcbddc->use_faces = PETSC_FALSE; 182447d04d0dSStefano Zampini pcbddc->use_change_of_basis = PETSC_FALSE; 182547d04d0dSStefano Zampini pcbddc->use_change_on_faces = PETSC_FALSE; 182647d04d0dSStefano Zampini pcbddc->switch_static = PETSC_FALSE; 182747d04d0dSStefano Zampini pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 182847d04d0dSStefano Zampini pcbddc->dbg_flag = 0; 1829b9d89cd5SStefano Zampini /* private */ 1830b9d89cd5SStefano Zampini pcbddc->issym = PETSC_FALSE; 183139e2fb2aSStefano Zampini pcbddc->BtoNmap = 0; 1832727cdba6SStefano Zampini pcbddc->local_primal_size = 0; 1833e9189074SStefano Zampini pcbddc->n_vertices = 0; 1834e9189074SStefano Zampini pcbddc->n_actual_vertices = 0; 1835e9189074SStefano Zampini pcbddc->n_constraints = 0; 1836727cdba6SStefano Zampini pcbddc->primal_indices_local_idxs = 0; 1837fb180af4SStefano Zampini pcbddc->recompute_topography = PETSC_FALSE; 183868457ee5SStefano Zampini pcbddc->coarse_size = -1; 1839f4ddd8eeSStefano Zampini pcbddc->new_primal_space = PETSC_FALSE; 1840727cdba6SStefano Zampini pcbddc->new_primal_space_local = PETSC_FALSE; 1841f4ddd8eeSStefano Zampini pcbddc->global_primal_indices = 0; 1842f4ddd8eeSStefano Zampini pcbddc->onearnullspace = 0; 1843f4ddd8eeSStefano Zampini pcbddc->onearnullvecs_state = 0; 1844674ae819SStefano Zampini pcbddc->user_primal_vertices = 0; 18450bdf917eSStefano Zampini pcbddc->NullSpace = 0; 18463972b0daSStefano Zampini pcbddc->temp_solution = 0; 1847534831adSStefano Zampini pcbddc->original_rhs = 0; 1848534831adSStefano Zampini pcbddc->local_mat = 0; 1849534831adSStefano Zampini pcbddc->ChangeOfBasisMatrix = 0; 1850b9b85e73SStefano Zampini pcbddc->user_ChangeOfBasisMatrix = 0; 1851906d46d4SStefano Zampini pcbddc->new_global_mat = 0; 1852da1bb401SStefano Zampini pcbddc->coarse_vec = 0; 1853da1bb401SStefano Zampini pcbddc->coarse_rhs = 0; 1854da1bb401SStefano Zampini pcbddc->coarse_ksp = 0; 1855da1bb401SStefano Zampini pcbddc->coarse_phi_B = 0; 1856da1bb401SStefano Zampini pcbddc->coarse_phi_D = 0; 185715aaf578SStefano Zampini pcbddc->coarse_psi_B = 0; 185815aaf578SStefano Zampini pcbddc->coarse_psi_D = 0; 1859da1bb401SStefano Zampini pcbddc->vec1_P = 0; 1860da1bb401SStefano Zampini pcbddc->vec1_R = 0; 1861da1bb401SStefano Zampini pcbddc->vec2_R = 0; 1862da1bb401SStefano Zampini pcbddc->local_auxmat1 = 0; 1863da1bb401SStefano Zampini pcbddc->local_auxmat2 = 0; 1864da1bb401SStefano Zampini pcbddc->R_to_B = 0; 1865da1bb401SStefano Zampini pcbddc->R_to_D = 0; 1866da1bb401SStefano Zampini pcbddc->ksp_D = 0; 1867da1bb401SStefano Zampini pcbddc->ksp_R = 0; 1868da1bb401SStefano Zampini pcbddc->NeumannBoundaries = 0; 1869785d1243SStefano Zampini pcbddc->NeumannBoundariesLocal = 0; 1870785d1243SStefano Zampini pcbddc->DirichletBoundaries = 0; 1871785d1243SStefano Zampini pcbddc->DirichletBoundariesLocal = 0; 187260d44989SStefano Zampini pcbddc->user_provided_isfordofs = PETSC_FALSE; 187360d44989SStefano Zampini pcbddc->n_ISForDofs = 0; 187463602bcaSStefano Zampini pcbddc->n_ISForDofsLocal = 0; 1875da1bb401SStefano Zampini pcbddc->ISForDofs = 0; 187663602bcaSStefano Zampini pcbddc->ISForDofsLocal = 0; 1877da1bb401SStefano Zampini pcbddc->ConstraintMatrix = 0; 187885c4d303SStefano Zampini pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 187947d04d0dSStefano Zampini pcbddc->coarse_loc_to_glob = 0; 188047d04d0dSStefano Zampini pcbddc->coarsening_ratio = 8; 18814fad6a16SStefano Zampini pcbddc->current_level = 0; 18822b510759SStefano Zampini pcbddc->max_levels = 0; 1883323d395dSStefano Zampini pcbddc->use_coarse_estimates = PETSC_FALSE; 188474e2c79eSStefano Zampini pcbddc->redistribute_coarse = 0; 1885f3bde8b3SStefano Zampini pcbddc->coarse_subassembling = 0; 1886323d395dSStefano Zampini pcbddc->coarse_subassembling_init = 0; 1887da1bb401SStefano Zampini 1888674ae819SStefano Zampini /* create local graph structure */ 1889674ae819SStefano Zampini ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1890674ae819SStefano Zampini 1891674ae819SStefano Zampini /* scaling */ 1892674ae819SStefano Zampini pcbddc->work_scaling = 0; 189334a97f8cSStefano Zampini pcbddc->use_deluxe_scaling = PETSC_FALSE; 189434a97f8cSStefano Zampini pcbddc->deluxe_threshold = 1; 189534a97f8cSStefano Zampini pcbddc->deluxe_rebuild = PETSC_FALSE; 189634a97f8cSStefano Zampini pcbddc->deluxe_layers = -1; 189734a97f8cSStefano Zampini pcbddc->deluxe_use_useradj = PETSC_FALSE; 189834a97f8cSStefano Zampini pcbddc->deluxe_compute_rowadj = PETSC_TRUE; 1899da1bb401SStefano Zampini 1900da1bb401SStefano Zampini /* function pointers */ 1901da1bb401SStefano Zampini pc->ops->apply = PCApply_BDDC; 190293bd9ae7SStefano Zampini pc->ops->applytranspose = PCApplyTranspose_BDDC; 1903da1bb401SStefano Zampini pc->ops->setup = PCSetUp_BDDC; 1904da1bb401SStefano Zampini pc->ops->destroy = PCDestroy_BDDC; 1905da1bb401SStefano Zampini pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1906da1bb401SStefano Zampini pc->ops->view = 0; 1907da1bb401SStefano Zampini pc->ops->applyrichardson = 0; 1908da1bb401SStefano Zampini pc->ops->applysymmetricleft = 0; 1909da1bb401SStefano Zampini pc->ops->applysymmetricright = 0; 1910534831adSStefano Zampini pc->ops->presolve = PCPreSolve_BDDC; 1911534831adSStefano Zampini pc->ops->postsolve = PCPostSolve_BDDC; 1912da1bb401SStefano Zampini 1913da1bb401SStefano Zampini /* composing function */ 1914906d46d4SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr); 1915674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1916bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 19172b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1918b8ffe317SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 19192b510759SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1920bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1921bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 192282ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1923bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 192482ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1925bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 192682ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1927bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 192882ba6b80SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1929bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 193063602bcaSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr); 1931bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1932bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1933bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1934bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1935da1bb401SStefano Zampini PetscFunctionReturn(0); 1936da1bb401SStefano Zampini } 19373425bc38SStefano Zampini 1938