1 /* TODOLIST 2 3 Solvers 4 - Add support for cholesky for coarse solver (similar to local solvers) 5 - Propagate ksp prefixes for solvers to mat objects? 6 7 User interface 8 - ** DM attached to pc? 9 10 Debugging output 11 - * Better management of verbosity levels of debugging output 12 13 Extra 14 - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)? 15 - BDDC with MG framework? 16 17 FETIDP 18 - Move FETIDP code to its own classes 19 20 MATIS related operations contained in BDDC code 21 - Provide general case for subassembling 22 23 */ 24 25 #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 26 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 27 #include <petscblaslapack.h> 28 29 /* temporarily declare it */ 30 PetscErrorCode PCApply_BDDC(PC,Vec,Vec); 31 32 /* -------------------------------------------------------------------------- */ 33 #undef __FUNCT__ 34 #define __FUNCT__ "PCSetFromOptions_BDDC" 35 PetscErrorCode PCSetFromOptions_BDDC(PetscOptionItems *PetscOptionsObject,PC pc) 36 { 37 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 38 PetscErrorCode ierr; 39 40 PetscFunctionBegin; 41 ierr = PetscOptionsHead(PetscOptionsObject,"BDDC options");CHKERRQ(ierr); 42 /* Verbose debugging */ 43 ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 44 /* Primal space customization */ 45 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); 46 ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr); 47 ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr); 48 ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr); 49 ierr = PetscOptionsInt("-pc_bddc_vertex_size","Connected components smaller or equal to vertex size will be considered as primal vertices","none",pcbddc->vertex_size,&pcbddc->vertex_size,NULL);CHKERRQ(ierr); 50 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); 51 ierr = PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is always used when multiple constraints are present)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);CHKERRQ(ierr); 52 /* Change of basis */ 53 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); 54 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); 55 if (!pcbddc->use_change_of_basis) { 56 pcbddc->use_change_on_faces = PETSC_FALSE; 57 } 58 /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 59 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); 60 ierr = PetscOptionsInt("-pc_bddc_coarse_eqs_per_proc","Number of equations per process for coarse problem redistribution (significant only at the coarsest level)","none",pcbddc->coarse_eqs_per_proc,&pcbddc->coarse_eqs_per_proc,NULL);CHKERRQ(ierr); 61 ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 62 ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 63 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); 64 ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr); 65 ierr = PetscOptionsBool("-pc_bddc_schur_rebuild","Whether or not the interface graph for Schur principal minors has to be rebuilt (i.e. define the interface without any adjacency)","none",pcbddc->sub_schurs_rebuild,&pcbddc->sub_schurs_rebuild,NULL);CHKERRQ(ierr); 66 ierr = PetscOptionsInt("-pc_bddc_schur_layers","Number of dofs' layers for the computation of principal minors (i.e. -1 uses all dofs)","none",pcbddc->sub_schurs_layers,&pcbddc->sub_schurs_layers,NULL);CHKERRQ(ierr); 67 ierr = PetscOptionsBool("-pc_bddc_schur_use_useradj","Whether or not the CSR graph specified by the user should be used for computing successive layers (default is to use adj of local mat)","none",pcbddc->sub_schurs_use_useradj,&pcbddc->sub_schurs_use_useradj,NULL);CHKERRQ(ierr); 68 ierr = PetscOptionsBool("-pc_bddc_schur_exact","Whether or not to use the exact Schur complement instead of the reduced one (which excludes size 1 cc)","none",pcbddc->sub_schurs_exact_schur,&pcbddc->sub_schurs_exact_schur,NULL);CHKERRQ(ierr); 69 ierr = PetscOptionsBool("-pc_bddc_deluxe_zerorows","Zero rows and columns of deluxe operators associated with primal dofs","none",pcbddc->deluxe_zerorows,&pcbddc->deluxe_zerorows,NULL);CHKERRQ(ierr); 70 ierr = PetscOptionsBool("-pc_bddc_adaptive_userdefined","Use user-defined constraints (should be attached via MatSetNearNullSpace to pmat) in addition to those adaptively generated","none",pcbddc->adaptive_userdefined,&pcbddc->adaptive_userdefined,NULL);CHKERRQ(ierr); 71 ierr = PetscOptionsReal("-pc_bddc_adaptive_threshold","Threshold to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&pcbddc->adaptive_threshold,NULL);CHKERRQ(ierr); 72 ierr = PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);CHKERRQ(ierr); 73 ierr = PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);CHKERRQ(ierr); 74 ierr = PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);CHKERRQ(ierr); 75 ierr = PetscOptionsInt("-pc_bddc_coarse_adj","Number of processors where to map the coarse adjacency list","none",pcbddc->coarse_adj_red,&pcbddc->coarse_adj_red,NULL);CHKERRQ(ierr); 76 ierr = PetscOptionsBool("-pc_bddc_benign_trick","Apply the benign subspace trick to saddle point problems with discontinuous pressures","none",pcbddc->benign_saddle_point,&pcbddc->benign_saddle_point,NULL);CHKERRQ(ierr); 77 ierr = PetscOptionsBool("-pc_bddc_benign_nonetflux","Automatic computation of no-net-flux quadrature weights","none",pcbddc->benign_compute_nonetflux,&pcbddc->benign_compute_nonetflux,NULL);CHKERRQ(ierr); 78 ierr = PetscOptionsBool("-pc_bddc_benign_change","Compute the pressure change of basis explicitly","none",pcbddc->benign_change_explicit,&pcbddc->benign_change_explicit,NULL);CHKERRQ(ierr); 79 ierr = PetscOptionsBool("-pc_bddc_detect_disconnected","Detects disconnected subdomains","none",pcbddc->detect_disconnected,&pcbddc->detect_disconnected,NULL);CHKERRQ(ierr); 80 ierr = PetscOptionsBool("-pc_bddc_eliminate_dirichlet","Whether or not we want to eliminate dirichlet dofs during presolve","none",pcbddc->eliminate_dirdofs,&pcbddc->eliminate_dirdofs,NULL);CHKERRQ(ierr); 81 ierr = PetscOptionsTail();CHKERRQ(ierr); 82 PetscFunctionReturn(0); 83 } 84 85 /* -------------------------------------------------------------------------- */ 86 #undef __FUNCT__ 87 #define __FUNCT__ "PCView_BDDC" 88 static PetscErrorCode PCView_BDDC(PC pc,PetscViewer viewer) 89 { 90 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 91 PC_IS *pcis = (PC_IS*)pc->data; 92 PetscErrorCode ierr; 93 PetscBool isascii,isstring; 94 PetscSubcomm subcomm; 95 PetscViewer subviewer; 96 97 PetscFunctionBegin; 98 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 99 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 100 /* Nothing printed for the String viewer */ 101 /* ASCII viewer */ 102 if (isascii) { 103 PetscMPIInt color,rank; 104 Petsc64bitInt loc[4],gsum[4],gmax[4],gmin[4]; 105 PetscScalar interface_size; 106 PetscReal ratio1=0.,ratio2=0.; 107 Vec counter; 108 109 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use verbose output: %d\n",pcbddc->dbg_flag);CHKERRQ(ierr); 110 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use user-defined CSR: %d\n",!!pcbddc->mat_graph->nvtxs_csr);CHKERRQ(ierr); 111 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use local mat graph: %d\n",pcbddc->use_local_adj);CHKERRQ(ierr); 112 if (pcbddc->mat_graph->twodim) { 113 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Connectivity graph topological dimension: 2\n");CHKERRQ(ierr); 114 } else { 115 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Connectivity graph topological dimension: 3\n");CHKERRQ(ierr); 116 } 117 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use vertices: %d (vertext size %d)\n",pcbddc->use_vertices,pcbddc->vertex_size);CHKERRQ(ierr); 118 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use edges: %d\n",pcbddc->use_edges);CHKERRQ(ierr); 119 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use faces: %d\n",pcbddc->use_faces);CHKERRQ(ierr); 120 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use true near null space: %d\n",pcbddc->use_nnsp_true);CHKERRQ(ierr); 121 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use QR for single constraints on cc: %d\n",pcbddc->use_qr_single);CHKERRQ(ierr); 122 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use change of basis on local edge nodes: %d\n",pcbddc->use_change_of_basis);CHKERRQ(ierr); 123 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use change of basis on local face nodes: %d\n",pcbddc->use_change_on_faces);CHKERRQ(ierr); 124 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Eliminate dirichlet dofs: %d\n",pcbddc->eliminate_dirdofs);CHKERRQ(ierr); 125 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Switch on static condensation ops around the interface preconditioner: %d\n",pcbddc->switch_static);CHKERRQ(ierr); 126 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Multilevel max levels: %d\n",pcbddc->max_levels);CHKERRQ(ierr); 127 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Multilevel coarsening ratio: %d\n",pcbddc->coarsening_ratio);CHKERRQ(ierr); 128 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use estimated eigs for coarse problem: %d\n",pcbddc->use_coarse_estimates);CHKERRQ(ierr); 129 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use deluxe scaling: %d\n",pcbddc->use_deluxe_scaling);CHKERRQ(ierr); 130 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use deluxe zerorows: %d\n",pcbddc->deluxe_zerorows);CHKERRQ(ierr); 131 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Rebuild interface graph for Schur principal minors: %d\n",pcbddc->sub_schurs_rebuild);CHKERRQ(ierr); 132 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Number of dofs' layers for the computation of principal minors: %d\n",pcbddc->sub_schurs_layers);CHKERRQ(ierr); 133 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Use user CSR graph to compute successive layers: %d\n",pcbddc->sub_schurs_use_useradj);CHKERRQ(ierr); 134 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Adaptive constraint selection threshold (active %d, userdefined %d): %g\n",pcbddc->adaptive_threshold,pcbddc->adaptive_selection,pcbddc->adaptive_userdefined);CHKERRQ(ierr); 135 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Min constraints / connected component: %d\n",pcbddc->adaptive_nmin);CHKERRQ(ierr); 136 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Max constraints / connected component: %d\n",pcbddc->adaptive_nmax);CHKERRQ(ierr); 137 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Invert exact Schur complement for adaptive selection: %d\n",pcbddc->sub_schurs_exact_schur);CHKERRQ(ierr); 138 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Symmetric computation of primal basis functions: %d\n",pcbddc->symmetric_primal);CHKERRQ(ierr); 139 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Num. Procs. to map coarse adjacency list: %d\n",pcbddc->coarse_adj_red);CHKERRQ(ierr); 140 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Coarse eqs per proc (significant at the coarsest level): %d\n",pcbddc->coarse_eqs_per_proc);CHKERRQ(ierr); 141 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Detect disconnected: %d\n",pcbddc->detect_disconnected);CHKERRQ(ierr); 142 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Benign subspace trick: %d (algebraic computation of no-net-flux %d)\n",pcbddc->benign_saddle_point,pcbddc->benign_compute_nonetflux);CHKERRQ(ierr); 143 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 144 145 /* compute some number on the domain decomposition */ 146 ierr = VecSet(pcis->vec1_B,1.0);CHKERRQ(ierr); 147 ierr = MatCreateVecs(pc->pmat,&counter,0);CHKERRQ(ierr); 148 ierr = VecSet(counter,0.0);CHKERRQ(ierr); 149 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 150 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,counter,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 151 ierr = VecSum(counter,&interface_size);CHKERRQ(ierr); 152 ierr = VecDestroy(&counter);CHKERRQ(ierr); 153 gsum[0] = 1; 154 gsum[1] = gsum[2] = gsum[3] = 0; 155 loc[0] = !!pcis->n; 156 loc[1] = pcis->n - pcis->n_B; 157 loc[2] = pcis->n_B; 158 loc[3] = pcbddc->local_primal_size; 159 ierr = MPI_Reduce(loc,gsum,4,MPIU_INT64,MPI_SUM,0,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 160 if (!loc[0]) loc[1] = loc[2] = loc[3] = -1; 161 ierr = MPI_Reduce(loc,gmax,4,MPIU_INT64,MPI_MAX,0,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 162 if (!loc[0]) loc[1] = loc[2] = loc[3] = PETSC_MAX_INT; 163 ierr = MPI_Reduce(loc,gmin,4,MPIU_INT64,MPI_MIN,0,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 164 if (pcbddc->coarse_size) { 165 ratio1 = pc->pmat->rmap->N/(1.*pcbddc->coarse_size); 166 ratio2 = PetscRealPart(interface_size)/pcbddc->coarse_size; 167 } 168 ierr = PetscViewerASCIIPrintf(viewer," BDDC: ********************************** STATISTICS AT LEVEL %d **********************************\n",pcbddc->current_level);CHKERRQ(ierr); 169 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Global dofs sizes: all %d interface %d coarse %d\n",pc->pmat->rmap->N,(PetscInt)PetscRealPart(interface_size),pcbddc->coarse_size);CHKERRQ(ierr); 170 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Coarsening ratios: all/coarse %d interface/coarse %d\n",(PetscInt)ratio1,(PetscInt)ratio2);CHKERRQ(ierr); 171 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Active processes : %d\n",(PetscInt)gsum[0]);CHKERRQ(ierr); 172 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Dofs type :\tMIN\tMAX\tMEAN\n",(PetscInt)gmin[1],(PetscInt)gmax[1],(PetscInt)(gsum[1]/gsum[0]));CHKERRQ(ierr); 173 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Interior dofs :\t%d\t%d\t%d\n",(PetscInt)gmin[1],(PetscInt)gmax[1],(PetscInt)(gsum[1]/gsum[0]));CHKERRQ(ierr); 174 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Interface dofs :\t%d\t%d\t%d\n",(PetscInt)gmin[2],(PetscInt)gmax[2],(PetscInt)(gsum[2]/gsum[0]));CHKERRQ(ierr); 175 ierr = PetscViewerASCIIPrintf(viewer," BDDC: Primal dofs :\t%d\t%d\t%d\n",(PetscInt)gmin[3],(PetscInt)gmax[3],(PetscInt)(gsum[3]/gsum[0]));CHKERRQ(ierr); 176 ierr = PetscViewerASCIIPrintf(viewer," BDDC: *******************************************************************************************\n");CHKERRQ(ierr); 177 178 /* the coarse problem can be handled with a different communicator */ 179 if (pcbddc->coarse_ksp) color = 1; 180 else color = 0; 181 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 182 ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)pc),&subcomm);CHKERRQ(ierr); 183 ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr); 184 ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr); 185 ierr = PetscViewerGetSubViewer(viewer,PetscSubcommChild(subcomm),&subviewer);CHKERRQ(ierr); 186 if (color == 1) { 187 ierr = KSPView(pcbddc->coarse_ksp,subviewer);CHKERRQ(ierr); 188 ierr = PetscViewerFlush(subviewer);CHKERRQ(ierr); 189 } 190 ierr = PetscViewerRestoreSubViewer(viewer,PetscSubcommChild(subcomm),&subviewer);CHKERRQ(ierr); 191 ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr); 192 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 193 } 194 PetscFunctionReturn(0); 195 } 196 197 /* -------------------------------------------------------------------------- */ 198 #undef __FUNCT__ 199 #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC" 200 static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change, PetscBool interior) 201 { 202 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 203 PetscErrorCode ierr; 204 205 PetscFunctionBegin; 206 ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr); 207 ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr); 208 pcbddc->user_ChangeOfBasisMatrix = change; 209 pcbddc->change_interior = interior; 210 PetscFunctionReturn(0); 211 } 212 #undef __FUNCT__ 213 #define __FUNCT__ "PCBDDCSetChangeOfBasisMat" 214 /*@ 215 PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs 216 217 Collective on PC 218 219 Input Parameters: 220 + pc - the preconditioning context 221 . change - the change of basis matrix 222 - interior - whether or not the change of basis modifies interior dofs 223 224 Level: intermediate 225 226 Notes: 227 228 .seealso: PCBDDC 229 @*/ 230 PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change, PetscBool interior) 231 { 232 PetscErrorCode ierr; 233 234 PetscFunctionBegin; 235 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 236 PetscValidHeaderSpecific(change,MAT_CLASSID,2); 237 PetscCheckSameComm(pc,1,change,2); 238 if (pc->mat) { 239 PetscInt rows_c,cols_c,rows,cols; 240 ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr); 241 ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr); 242 if (rows_c != rows) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows); 243 if (cols_c != cols) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols); 244 ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr); 245 ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr); 246 if (rows_c != rows) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows); 247 if (cols_c != cols) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols); 248 } 249 ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat,PetscBool),(pc,change,interior));CHKERRQ(ierr); 250 PetscFunctionReturn(0); 251 } 252 /* -------------------------------------------------------------------------- */ 253 #undef __FUNCT__ 254 #define __FUNCT__ "PCBDDCSetPrimalVerticesIS_BDDC" 255 static PetscErrorCode PCBDDCSetPrimalVerticesIS_BDDC(PC pc, IS PrimalVertices) 256 { 257 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 258 PetscBool isequal = PETSC_FALSE; 259 PetscErrorCode ierr; 260 261 PetscFunctionBegin; 262 ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 263 if (pcbddc->user_primal_vertices) { 264 ierr = ISEqual(PrimalVertices,pcbddc->user_primal_vertices,&isequal);CHKERRQ(ierr); 265 } 266 ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 267 ierr = ISDestroy(&pcbddc->user_primal_vertices_local);CHKERRQ(ierr); 268 pcbddc->user_primal_vertices = PrimalVertices; 269 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 270 PetscFunctionReturn(0); 271 } 272 #undef __FUNCT__ 273 #define __FUNCT__ "PCBDDCSetPrimalVerticesIS" 274 /*@ 275 PCBDDCSetPrimalVerticesIS - Set additional user defined primal vertices in PCBDDC 276 277 Collective 278 279 Input Parameters: 280 + pc - the preconditioning context 281 - PrimalVertices - index set of primal vertices in global numbering (can be empty) 282 283 Level: intermediate 284 285 Notes: 286 Any process can list any global node 287 288 .seealso: PCBDDC 289 @*/ 290 PetscErrorCode PCBDDCSetPrimalVerticesIS(PC pc, IS PrimalVertices) 291 { 292 PetscErrorCode ierr; 293 294 PetscFunctionBegin; 295 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 296 PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 297 PetscCheckSameComm(pc,1,PrimalVertices,2); 298 ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 299 PetscFunctionReturn(0); 300 } 301 /* -------------------------------------------------------------------------- */ 302 #undef __FUNCT__ 303 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 304 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 305 { 306 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 307 PetscBool isequal = PETSC_FALSE; 308 PetscErrorCode ierr; 309 310 PetscFunctionBegin; 311 ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 312 if (pcbddc->user_primal_vertices_local) { 313 ierr = ISEqual(PrimalVertices,pcbddc->user_primal_vertices_local,&isequal);CHKERRQ(ierr); 314 } 315 ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 316 ierr = ISDestroy(&pcbddc->user_primal_vertices_local);CHKERRQ(ierr); 317 pcbddc->user_primal_vertices_local = PrimalVertices; 318 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 319 PetscFunctionReturn(0); 320 } 321 #undef __FUNCT__ 322 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 323 /*@ 324 PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 325 326 Collective 327 328 Input Parameters: 329 + pc - the preconditioning context 330 - PrimalVertices - index set of primal vertices in local numbering (can be empty) 331 332 Level: intermediate 333 334 Notes: 335 336 .seealso: PCBDDC 337 @*/ 338 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 339 { 340 PetscErrorCode ierr; 341 342 PetscFunctionBegin; 343 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 344 PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 345 PetscCheckSameComm(pc,1,PrimalVertices,2); 346 ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 347 PetscFunctionReturn(0); 348 } 349 /* -------------------------------------------------------------------------- */ 350 #undef __FUNCT__ 351 #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 352 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 353 { 354 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 355 356 PetscFunctionBegin; 357 pcbddc->coarsening_ratio = k; 358 PetscFunctionReturn(0); 359 } 360 361 #undef __FUNCT__ 362 #define __FUNCT__ "PCBDDCSetCoarseningRatio" 363 /*@ 364 PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 365 366 Logically collective on PC 367 368 Input Parameters: 369 + pc - the preconditioning context 370 - k - coarsening ratio (H/h at the coarser level) 371 372 Options Database Keys: 373 . -pc_bddc_coarsening_ratio 374 375 Level: intermediate 376 377 Notes: 378 Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 379 380 .seealso: PCBDDC, PCBDDCSetLevels() 381 @*/ 382 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 383 { 384 PetscErrorCode ierr; 385 386 PetscFunctionBegin; 387 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 388 PetscValidLogicalCollectiveInt(pc,k,2); 389 ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 390 PetscFunctionReturn(0); 391 } 392 393 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 394 #undef __FUNCT__ 395 #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 396 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 397 { 398 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 399 400 PetscFunctionBegin; 401 pcbddc->use_exact_dirichlet_trick = flg; 402 PetscFunctionReturn(0); 403 } 404 405 #undef __FUNCT__ 406 #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 407 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 408 { 409 PetscErrorCode ierr; 410 411 PetscFunctionBegin; 412 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 413 PetscValidLogicalCollectiveBool(pc,flg,2); 414 ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 415 PetscFunctionReturn(0); 416 } 417 418 #undef __FUNCT__ 419 #define __FUNCT__ "PCBDDCSetLevel_BDDC" 420 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 421 { 422 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 423 424 PetscFunctionBegin; 425 pcbddc->current_level = level; 426 PetscFunctionReturn(0); 427 } 428 429 #undef __FUNCT__ 430 #define __FUNCT__ "PCBDDCSetLevel" 431 PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 432 { 433 PetscErrorCode ierr; 434 435 PetscFunctionBegin; 436 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 437 PetscValidLogicalCollectiveInt(pc,level,2); 438 ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 439 PetscFunctionReturn(0); 440 } 441 442 #undef __FUNCT__ 443 #define __FUNCT__ "PCBDDCSetLevels_BDDC" 444 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 445 { 446 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 447 448 PetscFunctionBegin; 449 pcbddc->max_levels = levels; 450 PetscFunctionReturn(0); 451 } 452 453 #undef __FUNCT__ 454 #define __FUNCT__ "PCBDDCSetLevels" 455 /*@ 456 PCBDDCSetLevels - Sets the maximum number of levels for multilevel 457 458 Logically collective on PC 459 460 Input Parameters: 461 + pc - the preconditioning context 462 - levels - the maximum number of levels (max 9) 463 464 Options Database Keys: 465 . -pc_bddc_levels 466 467 Level: intermediate 468 469 Notes: 470 Default value is 0, i.e. traditional one-level BDDC 471 472 .seealso: PCBDDC, PCBDDCSetCoarseningRatio() 473 @*/ 474 PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 475 { 476 PetscErrorCode ierr; 477 478 PetscFunctionBegin; 479 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 480 PetscValidLogicalCollectiveInt(pc,levels,2); 481 if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n"); 482 ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 483 PetscFunctionReturn(0); 484 } 485 /* -------------------------------------------------------------------------- */ 486 487 #undef __FUNCT__ 488 #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 489 static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 490 { 491 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 492 PetscErrorCode ierr; 493 494 PetscFunctionBegin; 495 ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 496 ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 497 pcbddc->NullSpace = NullSpace; 498 PetscFunctionReturn(0); 499 } 500 501 #undef __FUNCT__ 502 #define __FUNCT__ "PCBDDCSetNullSpace" 503 /*@ 504 PCBDDCSetNullSpace - Set nullspace for BDDC operator 505 506 Logically collective on PC and MatNullSpace 507 508 Input Parameters: 509 + pc - the preconditioning context 510 - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 511 512 Level: intermediate 513 514 Notes: 515 516 .seealso: PCBDDC 517 @*/ 518 PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 519 { 520 PetscErrorCode ierr; 521 522 PetscFunctionBegin; 523 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 524 PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 525 PetscCheckSameComm(pc,1,NullSpace,2); 526 ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 527 PetscFunctionReturn(0); 528 } 529 /* -------------------------------------------------------------------------- */ 530 531 #undef __FUNCT__ 532 #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 533 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 534 { 535 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 536 PetscBool isequal = PETSC_FALSE; 537 PetscErrorCode ierr; 538 539 PetscFunctionBegin; 540 ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 541 if (pcbddc->DirichletBoundaries) { 542 ierr = ISEqual(DirichletBoundaries,pcbddc->DirichletBoundaries,&isequal);CHKERRQ(ierr); 543 } 544 /* last user setting takes precendence -> destroy any other customization */ 545 ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 546 ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 547 pcbddc->DirichletBoundaries = DirichletBoundaries; 548 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 549 PetscFunctionReturn(0); 550 } 551 552 #undef __FUNCT__ 553 #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 554 /*@ 555 PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 556 557 Collective 558 559 Input Parameters: 560 + pc - the preconditioning context 561 - DirichletBoundaries - parallel IS defining the Dirichlet boundaries 562 563 Level: intermediate 564 565 Notes: 566 Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node 567 568 .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns() 569 @*/ 570 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 571 { 572 PetscErrorCode ierr; 573 574 PetscFunctionBegin; 575 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 576 PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 577 PetscCheckSameComm(pc,1,DirichletBoundaries,2); 578 ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 579 PetscFunctionReturn(0); 580 } 581 /* -------------------------------------------------------------------------- */ 582 583 #undef __FUNCT__ 584 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC" 585 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries) 586 { 587 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 588 PetscBool isequal = PETSC_FALSE; 589 PetscErrorCode ierr; 590 591 PetscFunctionBegin; 592 ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 593 if (pcbddc->DirichletBoundariesLocal) { 594 ierr = ISEqual(DirichletBoundaries,pcbddc->DirichletBoundariesLocal,&isequal);CHKERRQ(ierr); 595 } 596 /* last user setting takes precendence -> destroy any other customization */ 597 ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 598 ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 599 pcbddc->DirichletBoundariesLocal = DirichletBoundaries; 600 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 601 PetscFunctionReturn(0); 602 } 603 604 #undef __FUNCT__ 605 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal" 606 /*@ 607 PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering. 608 609 Collective 610 611 Input Parameters: 612 + pc - the preconditioning context 613 - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering) 614 615 Level: intermediate 616 617 Notes: 618 619 .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns() 620 @*/ 621 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries) 622 { 623 PetscErrorCode ierr; 624 625 PetscFunctionBegin; 626 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 627 PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 628 PetscCheckSameComm(pc,1,DirichletBoundaries,2); 629 ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 630 PetscFunctionReturn(0); 631 } 632 /* -------------------------------------------------------------------------- */ 633 634 #undef __FUNCT__ 635 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 636 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 637 { 638 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 639 PetscBool isequal = PETSC_FALSE; 640 PetscErrorCode ierr; 641 642 PetscFunctionBegin; 643 ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 644 if (pcbddc->NeumannBoundaries) { 645 ierr = ISEqual(NeumannBoundaries,pcbddc->NeumannBoundaries,&isequal);CHKERRQ(ierr); 646 } 647 /* last user setting takes precendence -> destroy any other customization */ 648 ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 649 ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 650 pcbddc->NeumannBoundaries = NeumannBoundaries; 651 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 652 PetscFunctionReturn(0); 653 } 654 655 #undef __FUNCT__ 656 #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 657 /*@ 658 PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 659 660 Collective 661 662 Input Parameters: 663 + pc - the preconditioning context 664 - NeumannBoundaries - parallel IS defining the Neumann boundaries 665 666 Level: intermediate 667 668 Notes: 669 Any process can list any global node 670 671 .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal() 672 @*/ 673 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 674 { 675 PetscErrorCode ierr; 676 677 PetscFunctionBegin; 678 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 679 PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 680 PetscCheckSameComm(pc,1,NeumannBoundaries,2); 681 ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 682 PetscFunctionReturn(0); 683 } 684 /* -------------------------------------------------------------------------- */ 685 686 #undef __FUNCT__ 687 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC" 688 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries) 689 { 690 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 691 PetscBool isequal = PETSC_FALSE; 692 PetscErrorCode ierr; 693 694 PetscFunctionBegin; 695 ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 696 if (pcbddc->NeumannBoundariesLocal) { 697 ierr = ISEqual(NeumannBoundaries,pcbddc->NeumannBoundariesLocal,&isequal);CHKERRQ(ierr); 698 } 699 /* last user setting takes precendence -> destroy any other customization */ 700 ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 701 ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 702 pcbddc->NeumannBoundariesLocal = NeumannBoundaries; 703 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 704 PetscFunctionReturn(0); 705 } 706 707 #undef __FUNCT__ 708 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal" 709 /*@ 710 PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering. 711 712 Collective 713 714 Input Parameters: 715 + pc - the preconditioning context 716 - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering) 717 718 Level: intermediate 719 720 Notes: 721 722 .seealso: PCBDDC, PCBDDCSetNeumannBoundaries() 723 @*/ 724 PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries) 725 { 726 PetscErrorCode ierr; 727 728 PetscFunctionBegin; 729 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 730 PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 731 PetscCheckSameComm(pc,1,NeumannBoundaries,2); 732 ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 733 PetscFunctionReturn(0); 734 } 735 /* -------------------------------------------------------------------------- */ 736 737 #undef __FUNCT__ 738 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 739 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 740 { 741 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 742 743 PetscFunctionBegin; 744 *DirichletBoundaries = pcbddc->DirichletBoundaries; 745 PetscFunctionReturn(0); 746 } 747 748 #undef __FUNCT__ 749 #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 750 /*@ 751 PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries 752 753 Collective 754 755 Input Parameters: 756 . pc - the preconditioning context 757 758 Output Parameters: 759 . DirichletBoundaries - index set defining the Dirichlet boundaries 760 761 Level: intermediate 762 763 Notes: 764 The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries 765 766 .seealso: PCBDDC 767 @*/ 768 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 769 { 770 PetscErrorCode ierr; 771 772 PetscFunctionBegin; 773 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 774 ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 775 PetscFunctionReturn(0); 776 } 777 /* -------------------------------------------------------------------------- */ 778 779 #undef __FUNCT__ 780 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC" 781 static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries) 782 { 783 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 784 785 PetscFunctionBegin; 786 *DirichletBoundaries = pcbddc->DirichletBoundariesLocal; 787 PetscFunctionReturn(0); 788 } 789 790 #undef __FUNCT__ 791 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal" 792 /*@ 793 PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering) 794 795 Collective 796 797 Input Parameters: 798 . pc - the preconditioning context 799 800 Output Parameters: 801 . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 802 803 Level: intermediate 804 805 Notes: 806 The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetDirichletBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetDirichletBoundaries). 807 In the latter case, the IS will be available after PCSetUp. 808 809 .seealso: PCBDDC 810 @*/ 811 PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries) 812 { 813 PetscErrorCode ierr; 814 815 PetscFunctionBegin; 816 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 817 ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 818 PetscFunctionReturn(0); 819 } 820 /* -------------------------------------------------------------------------- */ 821 822 #undef __FUNCT__ 823 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 824 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 825 { 826 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 827 828 PetscFunctionBegin; 829 *NeumannBoundaries = pcbddc->NeumannBoundaries; 830 PetscFunctionReturn(0); 831 } 832 833 #undef __FUNCT__ 834 #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 835 /*@ 836 PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries 837 838 Collective 839 840 Input Parameters: 841 . pc - the preconditioning context 842 843 Output Parameters: 844 . NeumannBoundaries - index set defining the Neumann boundaries 845 846 Level: intermediate 847 848 Notes: 849 The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries 850 851 .seealso: PCBDDC 852 @*/ 853 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 854 { 855 PetscErrorCode ierr; 856 857 PetscFunctionBegin; 858 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 859 ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 860 PetscFunctionReturn(0); 861 } 862 /* -------------------------------------------------------------------------- */ 863 864 #undef __FUNCT__ 865 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC" 866 static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries) 867 { 868 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 869 870 PetscFunctionBegin; 871 *NeumannBoundaries = pcbddc->NeumannBoundariesLocal; 872 PetscFunctionReturn(0); 873 } 874 875 #undef __FUNCT__ 876 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal" 877 /*@ 878 PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering) 879 880 Collective 881 882 Input Parameters: 883 . pc - the preconditioning context 884 885 Output Parameters: 886 . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 887 888 Level: intermediate 889 890 Notes: 891 The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetNeumannBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetNeumannBoundaries). 892 In the latter case, the IS will be available after PCSetUp. 893 894 .seealso: PCBDDC 895 @*/ 896 PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries) 897 { 898 PetscErrorCode ierr; 899 900 PetscFunctionBegin; 901 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 902 ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 903 PetscFunctionReturn(0); 904 } 905 /* -------------------------------------------------------------------------- */ 906 907 #undef __FUNCT__ 908 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 909 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 910 { 911 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 912 PCBDDCGraph mat_graph = pcbddc->mat_graph; 913 PetscBool same_data = PETSC_FALSE; 914 PetscErrorCode ierr; 915 916 PetscFunctionBegin; 917 if (!nvtxs) { 918 ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 919 PetscFunctionReturn(0); 920 } 921 if (mat_graph->nvtxs_csr == nvtxs) { 922 if (mat_graph->xadj == xadj && mat_graph->adjncy == adjncy) same_data = PETSC_TRUE; 923 if (!same_data && mat_graph->xadj[nvtxs] == xadj[nvtxs]) { 924 ierr = PetscMemcmp(adjncy,mat_graph->adjncy,nvtxs*sizeof(PetscInt),&same_data);CHKERRQ(ierr); 925 } 926 } 927 if (!same_data) { 928 /* free old CSR */ 929 ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 930 /* get CSR into graph structure */ 931 if (copymode == PETSC_COPY_VALUES) { 932 ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr); 933 ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr); 934 ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 935 ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 936 } else if (copymode == PETSC_OWN_POINTER) { 937 mat_graph->xadj = (PetscInt*)xadj; 938 mat_graph->adjncy = (PetscInt*)adjncy; 939 } 940 mat_graph->nvtxs_csr = nvtxs; 941 pcbddc->recompute_topography = PETSC_TRUE; 942 } 943 PetscFunctionReturn(0); 944 } 945 946 #undef __FUNCT__ 947 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 948 /*@ 949 PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local matrix 950 951 Not collective 952 953 Input Parameters: 954 + pc - the preconditioning context 955 . nvtxs - number of local vertices of the graph (i.e., the size of the local problem) 956 . xadj, adjncy - the CSR graph 957 - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 958 959 Level: intermediate 960 961 Notes: 962 963 .seealso: PCBDDC,PetscCopyMode 964 @*/ 965 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 966 { 967 void (*f)(void) = 0; 968 PetscErrorCode ierr; 969 970 PetscFunctionBegin; 971 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 972 if (nvtxs) { 973 PetscValidIntPointer(xadj,3); 974 PetscValidIntPointer(adjncy,4); 975 } 976 if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d",copymode); 977 ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 978 /* free arrays if PCBDDC is not the PC type */ 979 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 980 if (!f && copymode == PETSC_OWN_POINTER) { 981 ierr = PetscFree(xadj);CHKERRQ(ierr); 982 ierr = PetscFree(adjncy);CHKERRQ(ierr); 983 } 984 PetscFunctionReturn(0); 985 } 986 /* -------------------------------------------------------------------------- */ 987 988 #undef __FUNCT__ 989 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC" 990 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 991 { 992 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 993 PetscInt i; 994 PetscBool isequal = PETSC_FALSE; 995 PetscErrorCode ierr; 996 997 PetscFunctionBegin; 998 if (pcbddc->n_ISForDofsLocal == n_is) { 999 for (i=0;i<n_is;i++) { 1000 PetscBool isequalt; 1001 ierr = ISEqual(ISForDofs[i],pcbddc->ISForDofsLocal[i],&isequalt);CHKERRQ(ierr); 1002 if (!isequalt) break; 1003 } 1004 if (i == n_is) isequal = PETSC_TRUE; 1005 } 1006 for (i=0;i<n_is;i++) { 1007 ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 1008 } 1009 /* Destroy ISes if they were already set */ 1010 for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 1011 ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 1012 } 1013 ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 1014 /* last user setting takes precendence -> destroy any other customization */ 1015 for (i=0;i<pcbddc->n_ISForDofs;i++) { 1016 ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 1017 } 1018 ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 1019 pcbddc->n_ISForDofs = 0; 1020 /* allocate space then set */ 1021 if (n_is) { 1022 ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr); 1023 } 1024 for (i=0;i<n_is;i++) { 1025 pcbddc->ISForDofsLocal[i] = ISForDofs[i]; 1026 } 1027 pcbddc->n_ISForDofsLocal = n_is; 1028 if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 1029 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 1030 PetscFunctionReturn(0); 1031 } 1032 1033 #undef __FUNCT__ 1034 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal" 1035 /*@ 1036 PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix 1037 1038 Collective 1039 1040 Input Parameters: 1041 + pc - the preconditioning context 1042 . n_is - number of index sets defining the fields 1043 - ISForDofs - array of IS describing the fields in local ordering 1044 1045 Level: intermediate 1046 1047 Notes: 1048 n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field. 1049 1050 .seealso: PCBDDC 1051 @*/ 1052 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[]) 1053 { 1054 PetscInt i; 1055 PetscErrorCode ierr; 1056 1057 PetscFunctionBegin; 1058 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1059 PetscValidLogicalCollectiveInt(pc,n_is,2); 1060 for (i=0;i<n_is;i++) { 1061 PetscCheckSameComm(pc,1,ISForDofs[i],3); 1062 PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 1063 } 1064 ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 1065 PetscFunctionReturn(0); 1066 } 1067 /* -------------------------------------------------------------------------- */ 1068 1069 #undef __FUNCT__ 1070 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 1071 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 1072 { 1073 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1074 PetscInt i; 1075 PetscBool isequal = PETSC_FALSE; 1076 PetscErrorCode ierr; 1077 1078 PetscFunctionBegin; 1079 if (pcbddc->n_ISForDofs == n_is) { 1080 for (i=0;i<n_is;i++) { 1081 PetscBool isequalt; 1082 ierr = ISEqual(ISForDofs[i],pcbddc->ISForDofs[i],&isequalt);CHKERRQ(ierr); 1083 if (!isequalt) break; 1084 } 1085 if (i == n_is) isequal = PETSC_TRUE; 1086 } 1087 for (i=0;i<n_is;i++) { 1088 ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 1089 } 1090 /* Destroy ISes if they were already set */ 1091 for (i=0;i<pcbddc->n_ISForDofs;i++) { 1092 ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 1093 } 1094 ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 1095 /* last user setting takes precendence -> destroy any other customization */ 1096 for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 1097 ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 1098 } 1099 ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 1100 pcbddc->n_ISForDofsLocal = 0; 1101 /* allocate space then set */ 1102 if (n_is) { 1103 ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr); 1104 } 1105 for (i=0;i<n_is;i++) { 1106 pcbddc->ISForDofs[i] = ISForDofs[i]; 1107 } 1108 pcbddc->n_ISForDofs = n_is; 1109 if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 1110 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 1111 PetscFunctionReturn(0); 1112 } 1113 1114 #undef __FUNCT__ 1115 #define __FUNCT__ "PCBDDCSetDofsSplitting" 1116 /*@ 1117 PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix 1118 1119 Collective 1120 1121 Input Parameters: 1122 + pc - the preconditioning context 1123 . n_is - number of index sets defining the fields 1124 - ISForDofs - array of IS describing the fields in global ordering 1125 1126 Level: intermediate 1127 1128 Notes: 1129 Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field. 1130 1131 .seealso: PCBDDC 1132 @*/ 1133 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 1134 { 1135 PetscInt i; 1136 PetscErrorCode ierr; 1137 1138 PetscFunctionBegin; 1139 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1140 PetscValidLogicalCollectiveInt(pc,n_is,2); 1141 for (i=0;i<n_is;i++) { 1142 PetscCheckSameComm(pc,1,ISForDofs[i],3); 1143 PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 1144 } 1145 ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 1146 PetscFunctionReturn(0); 1147 } 1148 1149 /* -------------------------------------------------------------------------- */ 1150 #undef __FUNCT__ 1151 #define __FUNCT__ "PCPreSolve_BDDC" 1152 /* -------------------------------------------------------------------------- */ 1153 /* 1154 PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 1155 guess if a transformation of basis approach has been selected. 1156 1157 Input Parameter: 1158 + pc - the preconditioner contex 1159 1160 Application Interface Routine: PCPreSolve() 1161 1162 Notes: 1163 The interface routine PCPreSolve() is not usually called directly by 1164 the user, but instead is called by KSPSolve(). 1165 */ 1166 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1167 { 1168 PetscErrorCode ierr; 1169 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1170 PC_IS *pcis = (PC_IS*)(pc->data); 1171 Vec used_vec; 1172 PetscBool copy_rhs = PETSC_TRUE; 1173 PetscBool iscg = PETSC_FALSE; 1174 1175 PetscFunctionBegin; 1176 /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 1177 if (ksp) { 1178 ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 1179 if (!iscg || pcbddc->switch_static) { 1180 ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 1181 } 1182 } 1183 /* Creates parallel work vectors used in presolve */ 1184 if (!pcbddc->original_rhs) { 1185 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 1186 } 1187 if (!pcbddc->temp_solution) { 1188 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 1189 } 1190 1191 ierr = VecSet(pcbddc->temp_solution,0.);CHKERRQ(ierr); 1192 if (x) { 1193 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 1194 used_vec = x; 1195 } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */ 1196 ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 1197 used_vec = pcbddc->temp_solution; 1198 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 1199 } 1200 1201 /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */ 1202 if (ksp) { 1203 /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */ 1204 ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr); 1205 if (!pcbddc->ksp_guess_nonzero) { 1206 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 1207 } 1208 } 1209 1210 pcbddc->rhs_change = PETSC_FALSE; 1211 /* Take into account zeroed rows -> change rhs and store solution removed */ 1212 if (rhs && pcbddc->eliminate_dirdofs) { 1213 IS dirIS = NULL; 1214 1215 /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */ 1216 ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr); 1217 if (dirIS) { 1218 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 1219 PetscInt dirsize,i,*is_indices; 1220 PetscScalar *array_x; 1221 const PetscScalar *array_diagonal; 1222 1223 ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 1224 ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 1225 ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1226 ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1227 ierr = VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1228 ierr = VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1229 ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr); 1230 ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 1231 ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 1232 ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 1233 for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 1234 ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 1235 ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 1236 ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 1237 ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1238 ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1239 pcbddc->rhs_change = PETSC_TRUE; 1240 ierr = ISDestroy(&dirIS);CHKERRQ(ierr); 1241 } 1242 } 1243 1244 /* remove the computed solution or the initial guess from the rhs */ 1245 if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) { 1246 /* store the original rhs */ 1247 if (copy_rhs) { 1248 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1249 copy_rhs = PETSC_FALSE; 1250 } 1251 pcbddc->rhs_change = PETSC_TRUE; 1252 ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 1253 ierr = MatMultAdd(pc->mat,used_vec,rhs,rhs);CHKERRQ(ierr); 1254 ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 1255 ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 1256 if (ksp) { 1257 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr); 1258 } 1259 } 1260 ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 1261 1262 /* compute initial vector in benign space (TODO: what about FETI-DP?)*/ 1263 if (rhs && pcbddc->benign_saddle_point) { 1264 /* compute u^*_h as in Xuemin Tu's thesis (see Section 4.8.1) */ 1265 if (!pcbddc->benign_vec) { 1266 ierr = VecDuplicate(rhs,&pcbddc->benign_vec);CHKERRQ(ierr); 1267 } 1268 pcbddc->benign_apply_coarse_only = PETSC_TRUE; 1269 ierr = PCApply_BDDC(pc,rhs,pcbddc->benign_vec);CHKERRQ(ierr); 1270 pcbddc->benign_apply_coarse_only = PETSC_FALSE; 1271 /* TODO: what about Stokes? */ 1272 #if 0 1273 ierr = VecNorm(pcbddc->benign_vec,NORM_INFINITY,&norm);CHKERRQ(ierr); 1274 if (norm < PETSC_SMALL) benign_correction_is_zero = PETSC_TRUE; 1275 #endif 1276 } 1277 1278 /* remove non-benign solution from the rhs */ 1279 if (pcbddc->benign_saddle_point) { 1280 /* store the original rhs */ 1281 if (copy_rhs) { 1282 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1283 copy_rhs = PETSC_FALSE; 1284 } 1285 /* this code is disabled, until I test against incompressible Stokes */ 1286 #if 0 1287 if (benign_correction_is_zero) { /* still need to understand why it works great */ 1288 ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 1289 ierr = PCApply_BDDC(pc,rhs,pcbddc->benign_vec);CHKERRQ(ierr); 1290 } 1291 #endif 1292 ierr = VecScale(pcbddc->benign_vec,-1.0);CHKERRQ(ierr); 1293 ierr = MatMultAdd(pc->mat,pcbddc->benign_vec,rhs,rhs);CHKERRQ(ierr); 1294 pcbddc->rhs_change = PETSC_TRUE; 1295 ierr = VecAXPY(pcbddc->temp_solution,-1.0,pcbddc->benign_vec);CHKERRQ(ierr); 1296 } 1297 1298 /* set initial guess if using PCG */ 1299 if (x && pcbddc->use_exact_dirichlet_trick) { 1300 ierr = VecSet(x,0.0);CHKERRQ(ierr); 1301 if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) { 1302 if (pcbddc->benign_saddle_point) { /* we have already saved the changed rhs */ 1303 ierr = VecLockPop(pcis->vec1_global);CHKERRQ(ierr); 1304 } else { 1305 ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,rhs,pcis->vec1_global);CHKERRQ(ierr); 1306 } 1307 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_global,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1308 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec1_global,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1309 } else { 1310 ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1311 ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1312 } 1313 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 1314 if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) { 1315 ierr = VecSet(pcis->vec1_global,0.);CHKERRQ(ierr); 1316 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1317 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1318 ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_global,x);CHKERRQ(ierr); 1319 } else { 1320 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1321 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1322 } 1323 if (ksp) { 1324 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 1325 } 1326 } 1327 PetscFunctionReturn(0); 1328 } 1329 1330 /* -------------------------------------------------------------------------- */ 1331 #undef __FUNCT__ 1332 #define __FUNCT__ "PCPostSolve_BDDC" 1333 /* -------------------------------------------------------------------------- */ 1334 /* 1335 PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 1336 approach has been selected. Also, restores rhs to its original state. 1337 1338 Input Parameter: 1339 + pc - the preconditioner contex 1340 1341 Application Interface Routine: PCPostSolve() 1342 1343 Notes: 1344 The interface routine PCPostSolve() is not usually called directly by 1345 the user, but instead is called by KSPSolve(). 1346 */ 1347 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1348 { 1349 PetscErrorCode ierr; 1350 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1351 1352 PetscFunctionBegin; 1353 /* add solution removed in presolve */ 1354 if (x && pcbddc->rhs_change) { 1355 ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 1356 } 1357 /* restore rhs to its original state */ 1358 if (rhs && pcbddc->rhs_change) { 1359 ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 1360 } 1361 pcbddc->rhs_change = PETSC_FALSE; 1362 /* restore ksp guess state */ 1363 if (ksp) { 1364 ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr); 1365 } 1366 PetscFunctionReturn(0); 1367 } 1368 /* -------------------------------------------------------------------------- */ 1369 #undef __FUNCT__ 1370 #define __FUNCT__ "PCSetUp_BDDC" 1371 /* -------------------------------------------------------------------------- */ 1372 /* 1373 PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 1374 by setting data structures and options. 1375 1376 Input Parameter: 1377 + pc - the preconditioner context 1378 1379 Application Interface Routine: PCSetUp() 1380 1381 Notes: 1382 The interface routine PCSetUp() is not usually called directly by 1383 the user, but instead is called by PCApply() if necessary. 1384 */ 1385 PetscErrorCode PCSetUp_BDDC(PC pc) 1386 { 1387 PetscErrorCode ierr; 1388 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 1389 Mat_IS* matis; 1390 MatNullSpace nearnullspace; 1391 IS zerodiag = NULL; 1392 PetscInt nrows,ncols; 1393 PetscBool computetopography,computesolvers,computesubschurs; 1394 PetscBool computeconstraintsmatrix; 1395 PetscBool new_nearnullspace_provided,ismatis; 1396 1397 PetscFunctionBegin; 1398 ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr); 1399 if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS"); 1400 ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr); 1401 if (nrows != ncols) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix"); 1402 matis = (Mat_IS*)pc->pmat->data; 1403 /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 1404 /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1405 Also, BDDC directly build the Dirichlet problem */ 1406 /* split work */ 1407 if (pc->setupcalled) { 1408 if (pc->flag == SAME_NONZERO_PATTERN) { 1409 computetopography = PETSC_FALSE; 1410 computesolvers = PETSC_TRUE; 1411 } else { /* DIFFERENT_NONZERO_PATTERN */ 1412 computetopography = PETSC_TRUE; 1413 computesolvers = PETSC_TRUE; 1414 } 1415 } else { 1416 computetopography = PETSC_TRUE; 1417 computesolvers = PETSC_TRUE; 1418 } 1419 if (pcbddc->recompute_topography) { 1420 computetopography = PETSC_TRUE; 1421 } 1422 pcbddc->recompute_topography = computetopography; 1423 computeconstraintsmatrix = PETSC_FALSE; 1424 1425 /* check parameters' compatibility */ 1426 if (!pcbddc->use_deluxe_scaling) pcbddc->deluxe_zerorows = PETSC_FALSE; 1427 pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0); 1428 pcbddc->adaptive_userdefined = (PetscBool)(pcbddc->adaptive_selection && pcbddc->adaptive_userdefined); 1429 if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE; 1430 1431 computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling); 1432 if (pcbddc->switch_static) { 1433 PetscBool ismatis; 1434 ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&ismatis);CHKERRQ(ierr); 1435 if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"When the static switch is one, the iteration matrix should be of type MATIS"); 1436 } 1437 1438 /* Get stdout for dbg */ 1439 if (pcbddc->dbg_flag) { 1440 if (!pcbddc->dbg_viewer) { 1441 pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)); 1442 ierr = PetscViewerASCIIPushSynchronized(pcbddc->dbg_viewer);CHKERRQ(ierr); 1443 } 1444 ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1445 } 1446 1447 if (pcbddc->user_ChangeOfBasisMatrix) { 1448 /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */ 1449 pcbddc->use_change_of_basis = PETSC_FALSE; 1450 ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr); 1451 } else { 1452 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 1453 ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 1454 pcbddc->local_mat = matis->A; 1455 } 1456 1457 /* detect local disconnected subdomains if requested and not done before */ 1458 if (pcbddc->detect_disconnected && !pcbddc->n_local_subs) { 1459 ierr = MatDetectDisconnectedComponents(pcbddc->local_mat,PETSC_FALSE,&pcbddc->n_local_subs,&pcbddc->local_subs);CHKERRQ(ierr); 1460 } 1461 1462 /* 1463 Compute change of basis on local pressures (aka zerodiag dofs) 1464 This should come earlier then PCISSetUp for extracting the correct subdomain matrices 1465 */ 1466 ierr = PCBDDCBenignShellMat(pc,PETSC_TRUE);CHKERRQ(ierr); 1467 if (pcbddc->benign_saddle_point) { 1468 PC_IS* pcis = (PC_IS*)pc->data; 1469 1470 if (pcbddc->user_ChangeOfBasisMatrix || pcbddc->use_change_of_basis || !computesubschurs) pcbddc->benign_change_explicit = PETSC_TRUE; 1471 /* detect local saddle point and change the basis in pcbddc->local_mat */ 1472 ierr = PCBDDCBenignDetectSaddlePoint(pc,&zerodiag);CHKERRQ(ierr); 1473 /* pop B0 mat from local mat */ 1474 ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr); 1475 /* give pcis a hint to not reuse submatrices during PCISCreate */ 1476 if (pc->flag == SAME_NONZERO_PATTERN && pcis->reusesubmatrices == PETSC_TRUE) { 1477 if (pcbddc->benign_n && (pcbddc->benign_change_explicit || pcbddc->dbg_flag)) { 1478 pcis->reusesubmatrices = PETSC_FALSE; 1479 } else { 1480 pcis->reusesubmatrices = PETSC_TRUE; 1481 } 1482 } else { 1483 pcis->reusesubmatrices = PETSC_FALSE; 1484 } 1485 } 1486 /* propagate relevant information */ 1487 #if !defined(PETSC_USE_COMPLEX) /* workaround for reals */ 1488 if (matis->A->symmetric_set) { 1489 ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr); 1490 } 1491 #endif 1492 if (matis->A->symmetric_set) { 1493 ierr = MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);CHKERRQ(ierr); 1494 } 1495 if (matis->A->spd_set) { 1496 ierr = MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);CHKERRQ(ierr); 1497 } 1498 1499 /* Set up all the "iterative substructuring" common block without computing solvers */ 1500 { 1501 Mat temp_mat; 1502 1503 temp_mat = matis->A; 1504 matis->A = pcbddc->local_mat; 1505 ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr); 1506 pcbddc->local_mat = matis->A; 1507 matis->A = temp_mat; 1508 } 1509 1510 /* Analyze interface */ 1511 if (computetopography) { 1512 ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 1513 computeconstraintsmatrix = PETSC_TRUE; 1514 } 1515 1516 /* check existence of a divergence free extension, i.e. 1517 b(v_I,p_0) = 0 for all v_I (raise error if not). 1518 Also, check that PCBDDCBenignGetOrSetP0 works */ 1519 #if defined(PETSC_USE_DEBUG) 1520 if (pcbddc->benign_saddle_point) { 1521 ierr = PCBDDCBenignCheck(pc,zerodiag);CHKERRQ(ierr); 1522 } 1523 #endif 1524 ierr = ISDestroy(&zerodiag);CHKERRQ(ierr); 1525 1526 /* Setup local dirichlet solver ksp_D and sub_schurs solvers */ 1527 if (computesolvers) { 1528 PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs; 1529 1530 if (computesubschurs && computetopography) { 1531 ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr); 1532 } 1533 /* SetUp Scaling operator (scaling matrices could be needed in SubSchursSetUp)*/ 1534 if (!pcbddc->use_deluxe_scaling) { 1535 ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 1536 } 1537 if (sub_schurs->schur_explicit) { 1538 if (computesubschurs) { 1539 ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr); 1540 } 1541 ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 1542 } else { 1543 ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 1544 if (computesubschurs) { 1545 ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr); 1546 } 1547 } 1548 if (pcbddc->adaptive_selection) { 1549 ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr); 1550 computeconstraintsmatrix = PETSC_TRUE; 1551 } 1552 } 1553 1554 /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1555 new_nearnullspace_provided = PETSC_FALSE; 1556 ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr); 1557 if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1558 if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1559 new_nearnullspace_provided = PETSC_TRUE; 1560 } else { 1561 /* determine if the two nullspaces are different (should be lightweight) */ 1562 if (nearnullspace != pcbddc->onearnullspace) { 1563 new_nearnullspace_provided = PETSC_TRUE; 1564 } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1565 PetscInt i; 1566 const Vec *nearnullvecs; 1567 PetscObjectState state; 1568 PetscInt nnsp_size; 1569 ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 1570 for (i=0;i<nnsp_size;i++) { 1571 ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr); 1572 if (pcbddc->onearnullvecs_state[i] != state) { 1573 new_nearnullspace_provided = PETSC_TRUE; 1574 break; 1575 } 1576 } 1577 } 1578 } 1579 } else { 1580 if (!nearnullspace) { /* both nearnullspaces are null */ 1581 new_nearnullspace_provided = PETSC_FALSE; 1582 } else { /* nearnullspace attached later */ 1583 new_nearnullspace_provided = PETSC_TRUE; 1584 } 1585 } 1586 1587 /* Setup constraints and related work vectors */ 1588 /* reset primal space flags */ 1589 pcbddc->new_primal_space = PETSC_FALSE; 1590 pcbddc->new_primal_space_local = PETSC_FALSE; 1591 if (computeconstraintsmatrix || new_nearnullspace_provided) { 1592 /* It also sets the primal space flags */ 1593 ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 1594 /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1595 ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr); 1596 } 1597 1598 if (computesolvers || pcbddc->new_primal_space) { 1599 if (pcbddc->use_change_of_basis) { 1600 PC_IS *pcis = (PC_IS*)(pc->data); 1601 1602 ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 1603 if (pcbddc->benign_change) { 1604 ierr = MatDestroy(&pcbddc->benign_B0);CHKERRQ(ierr); 1605 /* pop B0 from pcbddc->local_mat */ 1606 ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr); 1607 } 1608 /* get submatrices */ 1609 ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 1610 ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 1611 ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 1612 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr); 1613 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr); 1614 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr); 1615 /* set flag in pcis to not reuse submatrices during PCISCreate */ 1616 pcis->reusesubmatrices = PETSC_FALSE; 1617 } else if (!pcbddc->user_ChangeOfBasisMatrix && !pcbddc->benign_change) { 1618 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 1619 ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 1620 pcbddc->local_mat = matis->A; 1621 } 1622 /* SetUp coarse and local Neumann solvers */ 1623 ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 1624 /* SetUp Scaling operator */ 1625 if (pcbddc->use_deluxe_scaling) { 1626 ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 1627 } 1628 } 1629 /* mark topography as done */ 1630 pcbddc->recompute_topography = PETSC_FALSE; 1631 1632 /* wrap pcis->A_IB and pcis->A_BI if we did not change explicitly the variables on the pressures */ 1633 ierr = PCBDDCBenignShellMat(pc,PETSC_FALSE);CHKERRQ(ierr); 1634 1635 if (pcbddc->dbg_flag) { 1636 ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1637 } 1638 PetscFunctionReturn(0); 1639 } 1640 1641 /* -------------------------------------------------------------------------- */ 1642 /* 1643 PCApply_BDDC - Applies the BDDC operator to a vector. 1644 1645 Input Parameters: 1646 + pc - the preconditioner context 1647 - r - input vector (global) 1648 1649 Output Parameter: 1650 . z - output vector (global) 1651 1652 Application Interface Routine: PCApply() 1653 */ 1654 #undef __FUNCT__ 1655 #define __FUNCT__ "PCApply_BDDC" 1656 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 1657 { 1658 PC_IS *pcis = (PC_IS*)(pc->data); 1659 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1660 PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B; 1661 PetscErrorCode ierr; 1662 const PetscScalar one = 1.0; 1663 const PetscScalar m_one = -1.0; 1664 const PetscScalar zero = 0.0; 1665 1666 /* This code is similar to that provided in nn.c for PCNN 1667 NN interface preconditioner changed to BDDC 1668 Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */ 1669 1670 PetscFunctionBegin; 1671 if (pcbddc->ChangeOfBasisMatrix) { 1672 Vec swap; 1673 ierr = VecCopy(r,pcbddc->work_change);CHKERRQ(ierr); 1674 swap = pcbddc->work_change; 1675 pcbddc->work_change = r; 1676 r = swap; 1677 ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,r);CHKERRQ(ierr); 1678 /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */ 1679 if (pcbddc->benign_apply_coarse_only && pcbddc->use_exact_dirichlet_trick && pcbddc->change_interior) { 1680 ierr = VecCopy(r,pcis->vec1_global);CHKERRQ(ierr); 1681 ierr = VecLockPush(pcis->vec1_global);CHKERRQ(ierr); 1682 } 1683 } 1684 if (pcbddc->benign_saddle_point) { /* get p0 from r */ 1685 ierr = PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);CHKERRQ(ierr); 1686 } 1687 if (!pcbddc->use_exact_dirichlet_trick) { 1688 ierr = VecCopy(r,z);CHKERRQ(ierr); 1689 /* First Dirichlet solve */ 1690 ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1691 ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1692 /* 1693 Assembling right hand side for BDDC operator 1694 - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE) 1695 - pcis->vec1_B the interface part of the global vector z 1696 */ 1697 if (n_D) { 1698 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 1699 ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 1700 if (pcbddc->switch_static) { 1701 Mat_IS *matis = (Mat_IS*)(pc->mat->data); 1702 1703 ierr = VecSet(pcis->vec1_N,0.);CHKERRQ(ierr); 1704 ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1705 ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1706 if (!pcbddc->switch_static_change) { 1707 ierr = MatMult(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1708 } else { 1709 ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1710 ierr = MatMult(matis->A,pcis->vec2_N,pcis->vec1_N);CHKERRQ(ierr); 1711 ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1712 } 1713 ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1714 ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1715 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1716 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1717 } else { 1718 ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 1719 } 1720 } else { 1721 ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr); 1722 } 1723 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1724 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1725 ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1726 } else { 1727 if (!pcbddc->benign_apply_coarse_only) { 1728 ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1729 } 1730 } 1731 1732 /* Apply interface preconditioner 1733 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1734 ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr); 1735 1736 /* Apply transpose of partition of unity operator */ 1737 ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 1738 1739 /* Second Dirichlet solve and assembling of output */ 1740 ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1741 ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1742 if (n_B) { 1743 if (pcbddc->switch_static) { 1744 Mat_IS *matis = (Mat_IS*)(pc->mat->data); 1745 1746 ierr = VecScatterBegin(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1747 ierr = VecScatterEnd(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1748 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1749 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1750 if (!pcbddc->switch_static_change) { 1751 ierr = MatMult(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1752 } else { 1753 ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1754 ierr = MatMult(matis->A,pcis->vec2_N,pcis->vec1_N);CHKERRQ(ierr); 1755 ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1756 } 1757 ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1758 ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1759 } else { 1760 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 1761 } 1762 } else if (pcbddc->switch_static) { /* n_B is zero */ 1763 Mat_IS *matis = (Mat_IS*)(pc->mat->data); 1764 1765 if (!pcbddc->switch_static_change) { 1766 ierr = MatMult(matis->A,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr); 1767 } else { 1768 ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_D,pcis->vec1_N);CHKERRQ(ierr); 1769 ierr = MatMult(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1770 ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec2_N,pcis->vec3_D);CHKERRQ(ierr); 1771 } 1772 } 1773 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1774 1775 if (!pcbddc->use_exact_dirichlet_trick) { 1776 if (pcbddc->switch_static) { 1777 ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr); 1778 } else { 1779 ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr); 1780 } 1781 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1782 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1783 } else { 1784 if (pcbddc->switch_static) { 1785 ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr); 1786 } else { 1787 ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1788 } 1789 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1790 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1791 } 1792 1793 if (pcbddc->benign_saddle_point) { /* set p0 (computed in PCBDDCApplyInterface) */ 1794 if (pcbddc->benign_apply_coarse_only) { 1795 ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 1796 } 1797 ierr = PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);CHKERRQ(ierr); 1798 } 1799 if (pcbddc->ChangeOfBasisMatrix) { 1800 Vec swap; 1801 1802 swap = r; 1803 r = pcbddc->work_change; 1804 pcbddc->work_change = swap; 1805 ierr = VecCopy(z,pcbddc->work_change);CHKERRQ(ierr); 1806 ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,z);CHKERRQ(ierr); 1807 } 1808 PetscFunctionReturn(0); 1809 } 1810 1811 /* -------------------------------------------------------------------------- */ 1812 /* 1813 PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector. 1814 1815 Input Parameters: 1816 + pc - the preconditioner context 1817 - r - input vector (global) 1818 1819 Output Parameter: 1820 . z - output vector (global) 1821 1822 Application Interface Routine: PCApplyTranspose() 1823 */ 1824 #undef __FUNCT__ 1825 #define __FUNCT__ "PCApplyTranspose_BDDC" 1826 PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z) 1827 { 1828 PC_IS *pcis = (PC_IS*)(pc->data); 1829 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1830 PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B; 1831 PetscErrorCode ierr; 1832 const PetscScalar one = 1.0; 1833 const PetscScalar m_one = -1.0; 1834 const PetscScalar zero = 0.0; 1835 1836 PetscFunctionBegin; 1837 if (pcbddc->ChangeOfBasisMatrix) { 1838 Vec swap; 1839 ierr = VecCopy(r,pcbddc->work_change);CHKERRQ(ierr); 1840 swap = pcbddc->work_change; 1841 pcbddc->work_change = r; 1842 r = swap; 1843 ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,r);CHKERRQ(ierr); 1844 } 1845 if (pcbddc->benign_saddle_point) { /* get p0 from r */ 1846 ierr = PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);CHKERRQ(ierr); 1847 } 1848 if (!pcbddc->use_exact_dirichlet_trick) { 1849 ierr = VecCopy(r,z);CHKERRQ(ierr); 1850 /* First Dirichlet solve */ 1851 ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1852 ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1853 /* 1854 Assembling right hand side for BDDC operator 1855 - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE) 1856 - pcis->vec1_B the interface part of the global vector z 1857 */ 1858 if (n_D) { 1859 ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 1860 ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 1861 if (pcbddc->switch_static) { 1862 Mat_IS *matis = (Mat_IS*)(pc->mat->data); 1863 1864 ierr = VecSet(pcis->vec1_N,0.);CHKERRQ(ierr); 1865 ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1866 ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1867 if (!pcbddc->switch_static_change) { 1868 ierr = MatMultTranspose(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1869 } else { 1870 ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1871 ierr = MatMultTranspose(matis->A,pcis->vec2_N,pcis->vec1_N);CHKERRQ(ierr); 1872 ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1873 } 1874 ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1875 ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1876 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1877 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1878 } else { 1879 ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 1880 } 1881 } else { 1882 ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr); 1883 } 1884 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1885 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1886 ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1887 } else { 1888 ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1889 } 1890 1891 /* Apply interface preconditioner 1892 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1893 ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr); 1894 1895 /* Apply transpose of partition of unity operator */ 1896 ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 1897 1898 /* Second Dirichlet solve and assembling of output */ 1899 ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1900 ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1901 if (n_B) { 1902 if (pcbddc->switch_static) { 1903 Mat_IS *matis = (Mat_IS*)(pc->mat->data); 1904 1905 ierr = VecScatterBegin(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1906 ierr = VecScatterEnd(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1907 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1908 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1909 if (!pcbddc->switch_static_change) { 1910 ierr = MatMultTranspose(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1911 } else { 1912 ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1913 ierr = MatMultTranspose(matis->A,pcis->vec2_N,pcis->vec1_N);CHKERRQ(ierr); 1914 ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1915 } 1916 ierr = VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1917 ierr = VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1918 } else { 1919 ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 1920 } 1921 } else if (pcbddc->switch_static) { /* n_B is zero */ 1922 Mat_IS *matis = (Mat_IS*)(pc->mat->data); 1923 1924 if (!pcbddc->switch_static_change) { 1925 ierr = MatMultTranspose(matis->A,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr); 1926 } else { 1927 ierr = MatMult(pcbddc->switch_static_change,pcis->vec1_D,pcis->vec1_N);CHKERRQ(ierr); 1928 ierr = MatMultTranspose(matis->A,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1929 ierr = MatMultTranspose(pcbddc->switch_static_change,pcis->vec2_N,pcis->vec3_D);CHKERRQ(ierr); 1930 } 1931 } 1932 ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1933 if (!pcbddc->use_exact_dirichlet_trick) { 1934 if (pcbddc->switch_static) { 1935 ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr); 1936 } else { 1937 ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr); 1938 } 1939 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1940 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1941 } else { 1942 if (pcbddc->switch_static) { 1943 ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr); 1944 } else { 1945 ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1946 } 1947 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1948 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1949 } 1950 if (pcbddc->benign_saddle_point) { /* set p0 (computed in PCBDDCApplyInterface) */ 1951 ierr = PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);CHKERRQ(ierr); 1952 } 1953 if (pcbddc->ChangeOfBasisMatrix) { 1954 Vec swap; 1955 1956 swap = r; 1957 r = pcbddc->work_change; 1958 pcbddc->work_change = swap; 1959 ierr = VecCopy(z,pcbddc->work_change);CHKERRQ(ierr); 1960 ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,z);CHKERRQ(ierr); 1961 } 1962 PetscFunctionReturn(0); 1963 } 1964 /* -------------------------------------------------------------------------- */ 1965 1966 #undef __FUNCT__ 1967 #define __FUNCT__ "PCDestroy_BDDC" 1968 PetscErrorCode PCDestroy_BDDC(PC pc) 1969 { 1970 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1971 PetscErrorCode ierr; 1972 1973 PetscFunctionBegin; 1974 /* free BDDC custom data */ 1975 ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 1976 /* destroy objects related to topography */ 1977 ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 1978 /* free allocated graph structure */ 1979 ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 1980 /* free allocated sub schurs structure */ 1981 ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr); 1982 /* destroy objects for scaling operator */ 1983 ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 1984 ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr); 1985 /* free solvers stuff */ 1986 ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 1987 /* free global vectors needed in presolve */ 1988 ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 1989 ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 1990 /* free data created by PCIS */ 1991 ierr = PCISDestroy(pc);CHKERRQ(ierr); 1992 /* remove functions */ 1993 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr); 1994 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 1995 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",NULL);CHKERRQ(ierr); 1996 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 1997 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 1998 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 1999 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 2000 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 2001 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 2002 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 2003 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 2004 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 2005 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 2006 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 2007 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 2008 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 2009 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 2010 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr); 2011 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 2012 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 2013 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 2014 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 2015 /* Free the private data structure */ 2016 ierr = PetscFree(pc->data);CHKERRQ(ierr); 2017 PetscFunctionReturn(0); 2018 } 2019 /* -------------------------------------------------------------------------- */ 2020 2021 #undef __FUNCT__ 2022 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 2023 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 2024 { 2025 FETIDPMat_ctx mat_ctx; 2026 Vec copy_standard_rhs; 2027 PC_IS* pcis; 2028 PC_BDDC* pcbddc; 2029 PetscErrorCode ierr; 2030 2031 PetscFunctionBegin; 2032 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 2033 pcis = (PC_IS*)mat_ctx->pc->data; 2034 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 2035 2036 /* 2037 change of basis for physical rhs if needed 2038 It also changes the rhs in case of dirichlet boundaries 2039 TODO: better management when FETIDP will have its own class 2040 */ 2041 ierr = VecDuplicate(standard_rhs,©_standard_rhs);CHKERRQ(ierr); 2042 ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr); 2043 ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr); 2044 /* store vectors for computation of fetidp final solution */ 2045 ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2046 ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2047 /* scale rhs since it should be unassembled */ 2048 /* TODO use counter scaling? (also below) */ 2049 ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2050 ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2051 /* Apply partition of unity */ 2052 ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 2053 /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 2054 if (!pcbddc->switch_static) { 2055 /* compute partially subassembled Schur complement right-hand side */ 2056 ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 2057 ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 2058 ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 2059 ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr); 2060 ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2061 ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2062 /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 2063 ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2064 ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2065 ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 2066 } 2067 ierr = VecDestroy(©_standard_rhs);CHKERRQ(ierr); 2068 /* BDDC rhs */ 2069 ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 2070 if (pcbddc->switch_static) { 2071 ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 2072 } 2073 /* apply BDDC */ 2074 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 2075 /* Application of B_delta and assembling of rhs for fetidp fluxes */ 2076 ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 2077 ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 2078 ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2079 ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2080 PetscFunctionReturn(0); 2081 } 2082 2083 #undef __FUNCT__ 2084 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 2085 /*@ 2086 PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side 2087 2088 Collective 2089 2090 Input Parameters: 2091 + fetidp_mat - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators 2092 - standard_rhs - the right-hand side of the original linear system 2093 2094 Output Parameters: 2095 . fetidp_flux_rhs - the right-hand side for the FETI-DP linear system 2096 2097 Level: developer 2098 2099 Notes: 2100 2101 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution 2102 @*/ 2103 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 2104 { 2105 FETIDPMat_ctx mat_ctx; 2106 PetscErrorCode ierr; 2107 2108 PetscFunctionBegin; 2109 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 2110 ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 2111 PetscFunctionReturn(0); 2112 } 2113 /* -------------------------------------------------------------------------- */ 2114 2115 #undef __FUNCT__ 2116 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 2117 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 2118 { 2119 FETIDPMat_ctx mat_ctx; 2120 PC_IS* pcis; 2121 PC_BDDC* pcbddc; 2122 PetscErrorCode ierr; 2123 2124 PetscFunctionBegin; 2125 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 2126 pcis = (PC_IS*)mat_ctx->pc->data; 2127 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 2128 2129 /* apply B_delta^T */ 2130 ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2131 ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2132 ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 2133 /* compute rhs for BDDC application */ 2134 ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 2135 if (pcbddc->switch_static) { 2136 ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 2137 } 2138 /* apply BDDC */ 2139 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 2140 /* put values into standard global vector */ 2141 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2142 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2143 if (!pcbddc->switch_static) { 2144 /* compute values into the interior if solved for the partially subassembled Schur complement */ 2145 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 2146 ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 2147 ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 2148 } 2149 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2150 ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2151 /* final change of basis if needed 2152 Is also sums the dirichlet part removed during RHS assembling */ 2153 ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 2154 PetscFunctionReturn(0); 2155 } 2156 2157 #undef __FUNCT__ 2158 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 2159 /*@ 2160 PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system 2161 2162 Collective 2163 2164 Input Parameters: 2165 + fetidp_mat - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators 2166 - fetidp_flux_sol - the solution of the FETI-DP linear system 2167 2168 Output Parameters: 2169 . standard_sol - the solution defined on the physical domain 2170 2171 Level: developer 2172 2173 Notes: 2174 2175 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS 2176 @*/ 2177 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 2178 { 2179 FETIDPMat_ctx mat_ctx; 2180 PetscErrorCode ierr; 2181 2182 PetscFunctionBegin; 2183 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 2184 ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 2185 PetscFunctionReturn(0); 2186 } 2187 /* -------------------------------------------------------------------------- */ 2188 2189 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 2190 extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec); 2191 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 2192 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 2193 extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec); 2194 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 2195 2196 #undef __FUNCT__ 2197 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 2198 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 2199 { 2200 2201 FETIDPMat_ctx fetidpmat_ctx; 2202 Mat newmat; 2203 FETIDPPC_ctx fetidppc_ctx; 2204 PC newpc; 2205 MPI_Comm comm; 2206 PetscErrorCode ierr; 2207 2208 PetscFunctionBegin; 2209 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 2210 /* FETIDP linear matrix */ 2211 ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 2212 ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 2213 ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 2214 ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 2215 ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr); 2216 ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 2217 ierr = MatSetUp(newmat);CHKERRQ(ierr); 2218 /* FETIDP preconditioner */ 2219 ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 2220 ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 2221 ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 2222 ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 2223 ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 2224 ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 2225 ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr); 2226 ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 2227 ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr); 2228 ierr = PCSetUp(newpc);CHKERRQ(ierr); 2229 /* return pointers for objects created */ 2230 *fetidp_mat=newmat; 2231 *fetidp_pc=newpc; 2232 PetscFunctionReturn(0); 2233 } 2234 2235 #undef __FUNCT__ 2236 #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 2237 /*@ 2238 PCBDDCCreateFETIDPOperators - Create FETI-DP operators 2239 2240 Collective 2241 2242 Input Parameters: 2243 . pc - the BDDC preconditioning context (setup should have been called before) 2244 2245 Output Parameters: 2246 + fetidp_mat - shell FETI-DP matrix object 2247 - fetidp_pc - shell Dirichlet preconditioner for FETI-DP matrix 2248 2249 Options Database Keys: 2250 . -fetidp_fullyredundant <false> - use or not a fully redundant set of Lagrange multipliers 2251 2252 Level: developer 2253 2254 Notes: 2255 Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose 2256 2257 .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution 2258 @*/ 2259 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 2260 { 2261 PetscErrorCode ierr; 2262 2263 PetscFunctionBegin; 2264 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2265 if (pc->setupcalled) { 2266 ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 2267 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 2268 PetscFunctionReturn(0); 2269 } 2270 /* -------------------------------------------------------------------------- */ 2271 /*MC 2272 PCBDDC - Balancing Domain Decomposition by Constraints. 2273 2274 An implementation of the BDDC preconditioner based on 2275 2276 .vb 2277 [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 2278 [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 2279 [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 2280 [4] C. Pechstein and C. R. Dohrmann. "Modern domain decomposition methods BDDC, deluxe scaling, and an algebraic approach", Seminar talk, Linz, December 2013, http://people.ricam.oeaw.ac.at/c.pechstein/pechstein-bddc2013.pdf 2281 .ve 2282 2283 The matrix to be preconditioned (Pmat) must be of type MATIS. 2284 2285 Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 2286 2287 It also works with unsymmetric and indefinite problems. 2288 2289 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. 2290 2291 Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace() 2292 2293 Boundary nodes are split in vertices, edges and faces classes 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() 2294 Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesIS() and their local counterparts. 2295 2296 Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD. 2297 2298 Change of basis is performed similarly to [2] when requested. When more than 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. 2299 User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat() 2300 2301 The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object. 2302 2303 Adaptive selection of primal constraints [4] is supported for SPD systems with high-contrast in the coefficients if MUMPS or MKL_PARDISO are present. Future versions of the code will also consider using PASTIX. 2304 2305 An experimental interface to the FETI-DP method is available. FETI-DP operators could be created using PCBDDCCreateFETIDPOperators(). A stand-alone class for the FETI-DP method will be provided in the next releases. 2306 Deluxe scaling is not supported yet for FETI-DP. 2307 2308 Options Database Keys (some of them, run with -h for a complete list): 2309 2310 . -pc_bddc_use_vertices <true> - use or not vertices in primal space 2311 . -pc_bddc_use_edges <true> - use or not edges in primal space 2312 . -pc_bddc_use_faces <false> - use or not faces in primal space 2313 . -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems 2314 . -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only) 2315 . -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested 2316 . -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1]) 2317 . -pc_bddc_levels <0> - maximum number of levels for multilevel 2318 . -pc_bddc_coarsening_ratio <8> - number of subdomains which will be aggregated together at the coarser level (e.g. H/h ratio at the coarser level, significative only in the multilevel case) 2319 . -pc_bddc_redistribute <0> - size of a subset of processors where the coarse problem will be remapped (the value is ignored if not at the coarsest level) 2320 . -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling 2321 . -pc_bddc_schur_layers <-1> - select the economic version of deluxe scaling by specifying the number of layers (-1 corresponds to the original deluxe scaling) 2322 . -pc_bddc_adaptive_threshold <0.0> - when a value greater than one is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS or MKL_PARDISO installed) 2323 - -pc_bddc_check_level <0> - set verbosity level of debugging output 2324 2325 Options for Dirichlet, Neumann or coarse solver can be set with 2326 .vb 2327 -pc_bddc_dirichlet_ 2328 -pc_bddc_neumann_ 2329 -pc_bddc_coarse_ 2330 .ve 2331 e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU. 2332 2333 When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as 2334 .vb 2335 -pc_bddc_dirichlet_lN_ 2336 -pc_bddc_neumann_lN_ 2337 -pc_bddc_coarse_lN_ 2338 .ve 2339 Note that level number ranges from the finest (0) to the coarsest (N). 2340 In order to specify options for the BDDC operators at the coarser levels (and not for the solvers), prepend -pc_bddc_coarse_ or -pc_bddc_coarse_l to the option, e.g. 2341 .vb 2342 -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3 2343 .ve 2344 will use a threshold of 5 for constraints' selection at the first coarse level and will redistribute the coarse problem of the first coarse level on 3 processors 2345 2346 Level: intermediate 2347 2348 Developer notes: 2349 2350 Contributed by Stefano Zampini 2351 2352 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 2353 M*/ 2354 2355 #undef __FUNCT__ 2356 #define __FUNCT__ "PCCreate_BDDC" 2357 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 2358 { 2359 PetscErrorCode ierr; 2360 PC_BDDC *pcbddc; 2361 2362 PetscFunctionBegin; 2363 /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 2364 ierr = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr); 2365 pc->data = (void*)pcbddc; 2366 2367 /* create PCIS data structure */ 2368 ierr = PCISCreate(pc);CHKERRQ(ierr); 2369 2370 /* BDDC customization */ 2371 pcbddc->use_local_adj = PETSC_TRUE; 2372 pcbddc->use_vertices = PETSC_TRUE; 2373 pcbddc->use_edges = PETSC_TRUE; 2374 pcbddc->use_faces = PETSC_FALSE; 2375 pcbddc->use_change_of_basis = PETSC_FALSE; 2376 pcbddc->use_change_on_faces = PETSC_FALSE; 2377 pcbddc->switch_static = PETSC_FALSE; 2378 pcbddc->use_nnsp_true = PETSC_FALSE; 2379 pcbddc->use_qr_single = PETSC_FALSE; 2380 pcbddc->symmetric_primal = PETSC_TRUE; 2381 pcbddc->benign_saddle_point = PETSC_FALSE; 2382 pcbddc->vertex_size = 1; 2383 pcbddc->dbg_flag = 0; 2384 /* private */ 2385 pcbddc->local_primal_size = 0; 2386 pcbddc->local_primal_size_cc = 0; 2387 pcbddc->local_primal_ref_node = 0; 2388 pcbddc->local_primal_ref_mult = 0; 2389 pcbddc->n_vertices = 0; 2390 pcbddc->primal_indices_local_idxs = 0; 2391 pcbddc->recompute_topography = PETSC_FALSE; 2392 pcbddc->coarse_size = -1; 2393 pcbddc->new_primal_space = PETSC_FALSE; 2394 pcbddc->new_primal_space_local = PETSC_FALSE; 2395 pcbddc->global_primal_indices = 0; 2396 pcbddc->onearnullspace = 0; 2397 pcbddc->onearnullvecs_state = 0; 2398 pcbddc->user_primal_vertices = 0; 2399 pcbddc->user_primal_vertices_local = 0; 2400 pcbddc->NullSpace = 0; 2401 pcbddc->temp_solution = 0; 2402 pcbddc->original_rhs = 0; 2403 pcbddc->local_mat = 0; 2404 pcbddc->ChangeOfBasisMatrix = 0; 2405 pcbddc->user_ChangeOfBasisMatrix = 0; 2406 pcbddc->coarse_vec = 0; 2407 pcbddc->coarse_ksp = 0; 2408 pcbddc->coarse_phi_B = 0; 2409 pcbddc->coarse_phi_D = 0; 2410 pcbddc->coarse_psi_B = 0; 2411 pcbddc->coarse_psi_D = 0; 2412 pcbddc->vec1_P = 0; 2413 pcbddc->vec1_R = 0; 2414 pcbddc->vec2_R = 0; 2415 pcbddc->local_auxmat1 = 0; 2416 pcbddc->local_auxmat2 = 0; 2417 pcbddc->R_to_B = 0; 2418 pcbddc->R_to_D = 0; 2419 pcbddc->ksp_D = 0; 2420 pcbddc->ksp_R = 0; 2421 pcbddc->NeumannBoundaries = 0; 2422 pcbddc->NeumannBoundariesLocal = 0; 2423 pcbddc->DirichletBoundaries = 0; 2424 pcbddc->DirichletBoundariesLocal = 0; 2425 pcbddc->user_provided_isfordofs = PETSC_FALSE; 2426 pcbddc->n_ISForDofs = 0; 2427 pcbddc->n_ISForDofsLocal = 0; 2428 pcbddc->ISForDofs = 0; 2429 pcbddc->ISForDofsLocal = 0; 2430 pcbddc->ConstraintMatrix = 0; 2431 pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 2432 pcbddc->coarse_loc_to_glob = 0; 2433 pcbddc->coarsening_ratio = 8; 2434 pcbddc->coarse_adj_red = 0; 2435 pcbddc->current_level = 0; 2436 pcbddc->max_levels = 0; 2437 pcbddc->use_coarse_estimates = PETSC_FALSE; 2438 pcbddc->coarse_eqs_per_proc = 1; 2439 pcbddc->coarse_subassembling = 0; 2440 pcbddc->detect_disconnected = PETSC_FALSE; 2441 pcbddc->n_local_subs = 0; 2442 pcbddc->local_subs = NULL; 2443 2444 /* benign subspace trick */ 2445 pcbddc->benign_change = 0; 2446 pcbddc->benign_vec = 0; 2447 pcbddc->benign_original_mat = 0; 2448 pcbddc->benign_sf = 0; 2449 pcbddc->benign_B0 = 0; 2450 pcbddc->benign_n = 0; 2451 pcbddc->benign_p0 = NULL; 2452 pcbddc->benign_p0_lidx = NULL; 2453 pcbddc->benign_p0_gidx = NULL; 2454 pcbddc->benign_null = PETSC_FALSE; 2455 2456 /* create local graph structure */ 2457 ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 2458 2459 /* scaling */ 2460 pcbddc->work_scaling = 0; 2461 pcbddc->use_deluxe_scaling = PETSC_FALSE; 2462 2463 /* create sub schurs structure */ 2464 ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr); 2465 pcbddc->sub_schurs_rebuild = PETSC_FALSE; 2466 pcbddc->sub_schurs_layers = -1; 2467 pcbddc->sub_schurs_use_useradj = PETSC_FALSE; 2468 2469 pcbddc->computed_rowadj = PETSC_FALSE; 2470 2471 /* adaptivity */ 2472 pcbddc->adaptive_threshold = 0.0; 2473 pcbddc->adaptive_nmax = 0; 2474 pcbddc->adaptive_nmin = 0; 2475 2476 /* function pointers */ 2477 pc->ops->apply = PCApply_BDDC; 2478 pc->ops->applytranspose = PCApplyTranspose_BDDC; 2479 pc->ops->setup = PCSetUp_BDDC; 2480 pc->ops->destroy = PCDestroy_BDDC; 2481 pc->ops->setfromoptions = PCSetFromOptions_BDDC; 2482 pc->ops->view = PCView_BDDC; 2483 pc->ops->applyrichardson = 0; 2484 pc->ops->applysymmetricleft = 0; 2485 pc->ops->applysymmetricright = 0; 2486 pc->ops->presolve = PCPreSolve_BDDC; 2487 pc->ops->postsolve = PCPostSolve_BDDC; 2488 2489 /* composing function */ 2490 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr); 2491 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 2492 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",PCBDDCSetPrimalVerticesIS_BDDC);CHKERRQ(ierr); 2493 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 2494 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 2495 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 2496 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 2497 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 2498 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 2499 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 2500 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 2501 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 2502 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 2503 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 2504 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 2505 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 2506 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 2507 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr); 2508 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 2509 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 2510 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 2511 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 2512 PetscFunctionReturn(0); 2513 } 2514 2515