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 MATIS related operations contained in BDDC code 18 - Provide general case for subassembling 19 20 */ 21 22 #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 23 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 24 #include <petscblaslapack.h> 25 26 static PetscBool PCBDDCPackageInitialized = PETSC_FALSE; 27 28 static PetscBool cited = PETSC_FALSE; 29 static const char citation[] = 30 "@article{ZampiniPCBDDC,\n" 31 "author = {Stefano Zampini},\n" 32 "title = {{PCBDDC}: A Class of Robust Dual-Primal Methods in {PETS}c},\n" 33 "journal = {SIAM Journal on Scientific Computing},\n" 34 "volume = {38},\n" 35 "number = {5},\n" 36 "pages = {S282-S306},\n" 37 "year = {2016},\n" 38 "doi = {10.1137/15M1025785},\n" 39 "URL = {http://dx.doi.org/10.1137/15M1025785},\n" 40 "eprint = {http://dx.doi.org/10.1137/15M1025785}\n" 41 "}\n"; 42 43 PetscLogEvent PC_BDDC_Topology[PETSC_PCBDDC_MAXLEVELS]; 44 PetscLogEvent PC_BDDC_LocalSolvers[PETSC_PCBDDC_MAXLEVELS]; 45 PetscLogEvent PC_BDDC_LocalWork[PETSC_PCBDDC_MAXLEVELS]; 46 PetscLogEvent PC_BDDC_CorrectionSetUp[PETSC_PCBDDC_MAXLEVELS]; 47 PetscLogEvent PC_BDDC_ApproxSetUp[PETSC_PCBDDC_MAXLEVELS]; 48 PetscLogEvent PC_BDDC_ApproxApply[PETSC_PCBDDC_MAXLEVELS]; 49 PetscLogEvent PC_BDDC_CoarseSetUp[PETSC_PCBDDC_MAXLEVELS]; 50 PetscLogEvent PC_BDDC_CoarseSolver[PETSC_PCBDDC_MAXLEVELS]; 51 PetscLogEvent PC_BDDC_AdaptiveSetUp[PETSC_PCBDDC_MAXLEVELS]; 52 PetscLogEvent PC_BDDC_Scaling[PETSC_PCBDDC_MAXLEVELS]; 53 PetscLogEvent PC_BDDC_Schurs[PETSC_PCBDDC_MAXLEVELS]; 54 PetscLogEvent PC_BDDC_Solves[PETSC_PCBDDC_MAXLEVELS][3]; 55 56 const char *const PCBDDCInterfaceExtTypes[] = {"DIRICHLET","LUMP","PCBDDCInterfaceExtType","PC_BDDC_INTERFACE_EXT_",NULL}; 57 58 PetscErrorCode PCApply_BDDC(PC,Vec,Vec); 59 60 PetscErrorCode PCSetFromOptions_BDDC(PetscOptionItems *PetscOptionsObject,PC pc) 61 { 62 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 63 PetscInt nt,i; 64 65 PetscFunctionBegin; 66 PetscCall(PetscOptionsHead(PetscOptionsObject,"BDDC options")); 67 /* Verbose debugging */ 68 PetscCall(PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL)); 69 /* Approximate solvers */ 70 PetscCall(PetscOptionsEnum("-pc_bddc_interface_ext_type","Use DIRICHLET or LUMP to extend interface corrections to interior","PCBDDCSetInterfaceExtType",PCBDDCInterfaceExtTypes,(PetscEnum)pcbddc->interface_extension,(PetscEnum*)&pcbddc->interface_extension,NULL)); 71 if (pcbddc->interface_extension == PC_BDDC_INTERFACE_EXT_DIRICHLET) { 72 PetscCall(PetscOptionsBool("-pc_bddc_dirichlet_approximate","Inform PCBDDC that we are using approximate Dirichlet solvers","none",pcbddc->NullSpace_corr[0],&pcbddc->NullSpace_corr[0],NULL)); 73 PetscCall(PetscOptionsBool("-pc_bddc_dirichlet_approximate_scale","Inform PCBDDC that we need to scale the Dirichlet solve","none",pcbddc->NullSpace_corr[1],&pcbddc->NullSpace_corr[1],NULL)); 74 } else { 75 /* This flag is needed/implied by lumping */ 76 pcbddc->switch_static = PETSC_TRUE; 77 } 78 PetscCall(PetscOptionsBool("-pc_bddc_neumann_approximate","Inform PCBDDC that we are using approximate Neumann solvers","none",pcbddc->NullSpace_corr[2],&pcbddc->NullSpace_corr[2],NULL)); 79 PetscCall(PetscOptionsBool("-pc_bddc_neumann_approximate_scale","Inform PCBDDC that we need to scale the Neumann solve","none",pcbddc->NullSpace_corr[3],&pcbddc->NullSpace_corr[3],NULL)); 80 /* Primal space customization */ 81 PetscCall(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)); 82 PetscCall(PetscOptionsInt("-pc_bddc_graph_maxcount","Maximum number of shared subdomains for a connected component","none",pcbddc->graphmaxcount,&pcbddc->graphmaxcount,NULL)); 83 PetscCall(PetscOptionsBool("-pc_bddc_corner_selection","Activates face-based corner selection","none",pcbddc->corner_selection,&pcbddc->corner_selection,NULL)); 84 PetscCall(PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL)); 85 PetscCall(PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL)); 86 PetscCall(PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL)); 87 PetscCall(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)); 88 PetscCall(PetscOptionsBool("-pc_bddc_use_nnsp","Use near null space attached to the matrix to compute constraints","none",pcbddc->use_nnsp,&pcbddc->use_nnsp,NULL)); 89 PetscCall(PetscOptionsBool("-pc_bddc_use_nnsp_true","Use near null space attached to the matrix to compute constraints as is","none",pcbddc->use_nnsp_true,&pcbddc->use_nnsp_true,NULL)); 90 PetscCall(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)); 91 /* Change of basis */ 92 PetscCall(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)); 93 PetscCall(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)); 94 if (!pcbddc->use_change_of_basis) { 95 pcbddc->use_change_on_faces = PETSC_FALSE; 96 } 97 /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 98 PetscCall(PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL)); 99 PetscCall(PetscOptionsInt("-pc_bddc_coarse_eqs_per_proc","Target 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)); 100 i = pcbddc->coarsening_ratio; 101 PetscCall(PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","PCBDDCSetCoarseningRatio",i,&i,NULL)); 102 PetscCall(PCBDDCSetCoarseningRatio(pc,i)); 103 i = pcbddc->max_levels; 104 PetscCall(PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","PCBDDCSetLevels",i,&i,NULL)); 105 PetscCall(PCBDDCSetLevels(pc,i)); 106 PetscCall(PetscOptionsInt("-pc_bddc_coarse_eqs_limit","Set maximum number of equations on coarsest grid to aim for","none",pcbddc->coarse_eqs_limit,&pcbddc->coarse_eqs_limit,NULL)); 107 PetscCall(PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL)); 108 PetscCall(PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL)); 109 PetscCall(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)); 110 PetscCall(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)); 111 PetscCall(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)); 112 PetscCall(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)); 113 PetscCall(PetscOptionsBool("-pc_bddc_deluxe_zerorows","Zero rows and columns of deluxe operators associated with primal dofs","none",pcbddc->deluxe_zerorows,&pcbddc->deluxe_zerorows,NULL)); 114 PetscCall(PetscOptionsBool("-pc_bddc_deluxe_singlemat","Collapse deluxe operators","none",pcbddc->deluxe_singlemat,&pcbddc->deluxe_singlemat,NULL)); 115 PetscCall(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)); 116 nt = 2; 117 PetscCall(PetscOptionsRealArray("-pc_bddc_adaptive_threshold","Thresholds to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&nt,NULL)); 118 if (nt == 1) pcbddc->adaptive_threshold[1] = pcbddc->adaptive_threshold[0]; 119 PetscCall(PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL)); 120 PetscCall(PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL)); 121 PetscCall(PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL)); 122 PetscCall(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)); 123 PetscCall(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)); 124 PetscCall(PetscOptionsBool("-pc_bddc_benign_change","Compute the pressure change of basis explicitly","none",pcbddc->benign_change_explicit,&pcbddc->benign_change_explicit,NULL)); 125 PetscCall(PetscOptionsBool("-pc_bddc_benign_compute_correction","Compute the benign correction during PreSolve","none",pcbddc->benign_compute_correction,&pcbddc->benign_compute_correction,NULL)); 126 PetscCall(PetscOptionsBool("-pc_bddc_nonetflux","Automatic computation of no-net-flux quadrature weights","none",pcbddc->compute_nonetflux,&pcbddc->compute_nonetflux,NULL)); 127 PetscCall(PetscOptionsBool("-pc_bddc_detect_disconnected","Detects disconnected subdomains","none",pcbddc->detect_disconnected,&pcbddc->detect_disconnected,NULL)); 128 PetscCall(PetscOptionsBool("-pc_bddc_detect_disconnected_filter","Filters out small entries in the local matrix when detecting disconnected subdomains","none",pcbddc->detect_disconnected_filter,&pcbddc->detect_disconnected_filter,NULL)); 129 PetscCall(PetscOptionsBool("-pc_bddc_eliminate_dirichlet","Whether or not we want to eliminate dirichlet dofs during presolve","none",pcbddc->eliminate_dirdofs,&pcbddc->eliminate_dirdofs,NULL)); 130 PetscCall(PetscOptionsTail()); 131 PetscFunctionReturn(0); 132 } 133 134 static PetscErrorCode PCView_BDDC(PC pc,PetscViewer viewer) 135 { 136 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 137 PC_IS *pcis = (PC_IS*)pc->data; 138 PetscBool isascii; 139 PetscSubcomm subcomm; 140 PetscViewer subviewer; 141 142 PetscFunctionBegin; 143 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 144 /* ASCII viewer */ 145 if (isascii) { 146 PetscMPIInt color,rank,size; 147 PetscInt64 loc[7],gsum[6],gmax[6],gmin[6],totbenign; 148 PetscScalar interface_size; 149 PetscReal ratio1=0.,ratio2=0.; 150 Vec counter; 151 152 if (!pc->setupcalled) { 153 PetscCall(PetscViewerASCIIPrintf(viewer," Partial information available: preconditioner has not been setup yet\n")); 154 } 155 PetscCall(PetscViewerASCIIPrintf(viewer," Use verbose output: %D\n",pcbddc->dbg_flag)); 156 PetscCall(PetscViewerASCIIPrintf(viewer," Use user-defined CSR: %d\n",!!pcbddc->mat_graph->nvtxs_csr)); 157 PetscCall(PetscViewerASCIIPrintf(viewer," Use local mat graph: %d\n",pcbddc->use_local_adj && !pcbddc->mat_graph->nvtxs_csr)); 158 if (pcbddc->mat_graph->twodim) { 159 PetscCall(PetscViewerASCIIPrintf(viewer," Connectivity graph topological dimension: 2\n")); 160 } else { 161 PetscCall(PetscViewerASCIIPrintf(viewer," Connectivity graph topological dimension: 3\n")); 162 } 163 if (pcbddc->graphmaxcount != PETSC_MAX_INT) { 164 PetscCall(PetscViewerASCIIPrintf(viewer," Graph max count: %D\n",pcbddc->graphmaxcount)); 165 } 166 PetscCall(PetscViewerASCIIPrintf(viewer," Corner selection: %d (selected %d)\n",pcbddc->corner_selection,pcbddc->corner_selected)); 167 PetscCall(PetscViewerASCIIPrintf(viewer," Use vertices: %d (vertex size %D)\n",pcbddc->use_vertices,pcbddc->vertex_size)); 168 PetscCall(PetscViewerASCIIPrintf(viewer," Use edges: %d\n",pcbddc->use_edges)); 169 PetscCall(PetscViewerASCIIPrintf(viewer," Use faces: %d\n",pcbddc->use_faces)); 170 PetscCall(PetscViewerASCIIPrintf(viewer," Use true near null space: %d\n",pcbddc->use_nnsp_true)); 171 PetscCall(PetscViewerASCIIPrintf(viewer," Use QR for single constraints on cc: %d\n",pcbddc->use_qr_single)); 172 PetscCall(PetscViewerASCIIPrintf(viewer," Use change of basis on local edge nodes: %d\n",pcbddc->use_change_of_basis)); 173 PetscCall(PetscViewerASCIIPrintf(viewer," Use change of basis on local face nodes: %d\n",pcbddc->use_change_on_faces)); 174 PetscCall(PetscViewerASCIIPrintf(viewer," User defined change of basis matrix: %d\n",!!pcbddc->user_ChangeOfBasisMatrix)); 175 PetscCall(PetscViewerASCIIPrintf(viewer," Has change of basis matrix: %d\n",!!pcbddc->ChangeOfBasisMatrix)); 176 PetscCall(PetscViewerASCIIPrintf(viewer," Eliminate dirichlet boundary dofs: %d\n",pcbddc->eliminate_dirdofs)); 177 PetscCall(PetscViewerASCIIPrintf(viewer," Switch on static condensation ops around the interface preconditioner: %d\n",pcbddc->switch_static)); 178 PetscCall(PetscViewerASCIIPrintf(viewer," Use exact dirichlet trick: %d\n",pcbddc->use_exact_dirichlet_trick)); 179 PetscCall(PetscViewerASCIIPrintf(viewer," Interface extension: %s\n",PCBDDCInterfaceExtTypes[pcbddc->interface_extension])); 180 PetscCall(PetscViewerASCIIPrintf(viewer," Multilevel max levels: %D\n",pcbddc->max_levels)); 181 PetscCall(PetscViewerASCIIPrintf(viewer," Multilevel coarsening ratio: %D\n",pcbddc->coarsening_ratio)); 182 PetscCall(PetscViewerASCIIPrintf(viewer," Use estimated eigs for coarse problem: %d\n",pcbddc->use_coarse_estimates)); 183 PetscCall(PetscViewerASCIIPrintf(viewer," Use deluxe scaling: %d\n",pcbddc->use_deluxe_scaling)); 184 PetscCall(PetscViewerASCIIPrintf(viewer," Use deluxe zerorows: %d\n",pcbddc->deluxe_zerorows)); 185 PetscCall(PetscViewerASCIIPrintf(viewer," Use deluxe singlemat: %d\n",pcbddc->deluxe_singlemat)); 186 PetscCall(PetscViewerASCIIPrintf(viewer," Rebuild interface graph for Schur principal minors: %d\n",pcbddc->sub_schurs_rebuild)); 187 PetscCall(PetscViewerASCIIPrintf(viewer," Number of dofs' layers for the computation of principal minors: %D\n",pcbddc->sub_schurs_layers)); 188 PetscCall(PetscViewerASCIIPrintf(viewer," Use user CSR graph to compute successive layers: %d\n",pcbddc->sub_schurs_use_useradj)); 189 if (pcbddc->adaptive_threshold[1] != pcbddc->adaptive_threshold[0]) { 190 PetscCall(PetscViewerASCIIPrintf(viewer," Adaptive constraint selection thresholds (active %d, userdefined %d): %g,%g\n",pcbddc->adaptive_selection,pcbddc->adaptive_userdefined,pcbddc->adaptive_threshold[0],pcbddc->adaptive_threshold[1])); 191 } else { 192 PetscCall(PetscViewerASCIIPrintf(viewer," Adaptive constraint selection threshold (active %d, userdefined %d): %g\n",pcbddc->adaptive_selection,pcbddc->adaptive_userdefined,pcbddc->adaptive_threshold[0])); 193 } 194 PetscCall(PetscViewerASCIIPrintf(viewer," Min constraints / connected component: %D\n",pcbddc->adaptive_nmin)); 195 PetscCall(PetscViewerASCIIPrintf(viewer," Max constraints / connected component: %D\n",pcbddc->adaptive_nmax)); 196 PetscCall(PetscViewerASCIIPrintf(viewer," Invert exact Schur complement for adaptive selection: %d\n",pcbddc->sub_schurs_exact_schur)); 197 PetscCall(PetscViewerASCIIPrintf(viewer," Symmetric computation of primal basis functions: %d\n",pcbddc->symmetric_primal)); 198 PetscCall(PetscViewerASCIIPrintf(viewer," Num. Procs. to map coarse adjacency list: %D\n",pcbddc->coarse_adj_red)); 199 PetscCall(PetscViewerASCIIPrintf(viewer," Coarse eqs per proc (significant at the coarsest level): %D\n",pcbddc->coarse_eqs_per_proc)); 200 PetscCall(PetscViewerASCIIPrintf(viewer," Detect disconnected: %d (filter %d)\n",pcbddc->detect_disconnected,pcbddc->detect_disconnected_filter)); 201 PetscCall(PetscViewerASCIIPrintf(viewer," Benign subspace trick: %d (change explicit %d)\n",pcbddc->benign_saddle_point,pcbddc->benign_change_explicit)); 202 PetscCall(PetscViewerASCIIPrintf(viewer," Benign subspace trick is active: %d\n",pcbddc->benign_have_null)); 203 PetscCall(PetscViewerASCIIPrintf(viewer," Algebraic computation of no-net-flux: %d\n",pcbddc->compute_nonetflux)); 204 if (!pc->setupcalled) PetscFunctionReturn(0); 205 206 /* compute interface size */ 207 PetscCall(VecSet(pcis->vec1_B,1.0)); 208 PetscCall(MatCreateVecs(pc->pmat,&counter,NULL)); 209 PetscCall(VecSet(counter,0.0)); 210 PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,INSERT_VALUES,SCATTER_REVERSE)); 211 PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_B,counter,INSERT_VALUES,SCATTER_REVERSE)); 212 PetscCall(VecSum(counter,&interface_size)); 213 PetscCall(VecDestroy(&counter)); 214 215 /* compute some statistics on the domain decomposition */ 216 gsum[0] = 1; 217 gsum[1] = gsum[2] = gsum[3] = gsum[4] = gsum[5] = 0; 218 loc[0] = !!pcis->n; 219 loc[1] = pcis->n - pcis->n_B; 220 loc[2] = pcis->n_B; 221 loc[3] = pcbddc->local_primal_size; 222 loc[4] = pcis->n; 223 loc[5] = pcbddc->n_local_subs > 0 ? pcbddc->n_local_subs : (pcis->n ? 1 : 0); 224 loc[6] = pcbddc->benign_n; 225 PetscCallMPI(MPI_Reduce(loc,gsum,6,MPIU_INT64,MPI_SUM,0,PetscObjectComm((PetscObject)pc))); 226 if (!loc[0]) loc[1] = loc[2] = loc[3] = loc[4] = loc[5] = -1; 227 PetscCallMPI(MPI_Reduce(loc,gmax,6,MPIU_INT64,MPI_MAX,0,PetscObjectComm((PetscObject)pc))); 228 if (!loc[0]) loc[1] = loc[2] = loc[3] = loc[4] = loc[5] = PETSC_MAX_INT; 229 PetscCallMPI(MPI_Reduce(loc,gmin,6,MPIU_INT64,MPI_MIN,0,PetscObjectComm((PetscObject)pc))); 230 PetscCallMPI(MPI_Reduce(&loc[6],&totbenign,1,MPIU_INT64,MPI_SUM,0,PetscObjectComm((PetscObject)pc))); 231 if (pcbddc->coarse_size) { 232 ratio1 = pc->pmat->rmap->N/(1.*pcbddc->coarse_size); 233 ratio2 = PetscRealPart(interface_size)/pcbddc->coarse_size; 234 } 235 PetscCall(PetscViewerASCIIPrintf(viewer,"********************************** STATISTICS AT LEVEL %d **********************************\n",pcbddc->current_level)); 236 PetscCall(PetscViewerASCIIPrintf(viewer," Global dofs sizes: all %D interface %D coarse %D\n",pc->pmat->rmap->N,(PetscInt)PetscRealPart(interface_size),pcbddc->coarse_size)); 237 PetscCall(PetscViewerASCIIPrintf(viewer," Coarsening ratios: all/coarse %D interface/coarse %D\n",(PetscInt)ratio1,(PetscInt)ratio2)); 238 PetscCall(PetscViewerASCIIPrintf(viewer," Active processes : %D\n",(PetscInt)gsum[0])); 239 PetscCall(PetscViewerASCIIPrintf(viewer," Total subdomains : %D\n",(PetscInt)gsum[5])); 240 if (pcbddc->benign_have_null) { 241 PetscCall(PetscViewerASCIIPrintf(viewer," Benign subs : %D\n",(PetscInt)totbenign)); 242 } 243 PetscCall(PetscViewerASCIIPrintf(viewer," Dofs type :\tMIN\tMAX\tMEAN\n")); 244 PetscCall(PetscViewerASCIIPrintf(viewer," Interior dofs :\t%D\t%D\t%D\n",(PetscInt)gmin[1],(PetscInt)gmax[1],(PetscInt)(gsum[1]/gsum[0]))); 245 PetscCall(PetscViewerASCIIPrintf(viewer," Interface dofs :\t%D\t%D\t%D\n",(PetscInt)gmin[2],(PetscInt)gmax[2],(PetscInt)(gsum[2]/gsum[0]))); 246 PetscCall(PetscViewerASCIIPrintf(viewer," Primal dofs :\t%D\t%D\t%D\n",(PetscInt)gmin[3],(PetscInt)gmax[3],(PetscInt)(gsum[3]/gsum[0]))); 247 PetscCall(PetscViewerASCIIPrintf(viewer," Local dofs :\t%D\t%D\t%D\n",(PetscInt)gmin[4],(PetscInt)gmax[4],(PetscInt)(gsum[4]/gsum[0]))); 248 PetscCall(PetscViewerASCIIPrintf(viewer," Local subs :\t%D\t%D\n" ,(PetscInt)gmin[5],(PetscInt)gmax[5])); 249 PetscCall(PetscViewerFlush(viewer)); 250 251 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 252 253 /* local solvers */ 254 PetscCall(PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pcbddc->ksp_D),&subviewer)); 255 if (rank == 0) { 256 PetscCall(PetscViewerASCIIPrintf(subviewer,"--- Interior solver (rank 0)\n")); 257 PetscCall(PetscViewerASCIIPushTab(subviewer)); 258 PetscCall(KSPView(pcbddc->ksp_D,subviewer)); 259 PetscCall(PetscViewerASCIIPopTab(subviewer)); 260 PetscCall(PetscViewerASCIIPrintf(subviewer,"--- Correction solver (rank 0)\n")); 261 PetscCall(PetscViewerASCIIPushTab(subviewer)); 262 PetscCall(KSPView(pcbddc->ksp_R,subviewer)); 263 PetscCall(PetscViewerASCIIPopTab(subviewer)); 264 PetscCall(PetscViewerFlush(subviewer)); 265 } 266 PetscCall(PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pcbddc->ksp_D),&subviewer)); 267 PetscCall(PetscViewerFlush(viewer)); 268 269 /* the coarse problem can be handled by a different communicator */ 270 if (pcbddc->coarse_ksp) color = 1; 271 else color = 0; 272 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 273 PetscCall(PetscSubcommCreate(PetscObjectComm((PetscObject)pc),&subcomm)); 274 PetscCall(PetscSubcommSetNumber(subcomm,PetscMin(size,2))); 275 PetscCall(PetscSubcommSetTypeGeneral(subcomm,color,rank)); 276 PetscCall(PetscViewerGetSubViewer(viewer,PetscSubcommChild(subcomm),&subviewer)); 277 if (color == 1) { 278 PetscCall(PetscViewerASCIIPrintf(subviewer,"--- Coarse solver\n")); 279 PetscCall(PetscViewerASCIIPushTab(subviewer)); 280 PetscCall(KSPView(pcbddc->coarse_ksp,subviewer)); 281 PetscCall(PetscViewerASCIIPopTab(subviewer)); 282 PetscCall(PetscViewerFlush(subviewer)); 283 } 284 PetscCall(PetscViewerRestoreSubViewer(viewer,PetscSubcommChild(subcomm),&subviewer)); 285 PetscCall(PetscSubcommDestroy(&subcomm)); 286 PetscCall(PetscViewerFlush(viewer)); 287 } 288 PetscFunctionReturn(0); 289 } 290 291 static PetscErrorCode PCBDDCSetDiscreteGradient_BDDC(PC pc, Mat G, PetscInt order, PetscInt field, PetscBool global, PetscBool conforming) 292 { 293 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 294 295 PetscFunctionBegin; 296 PetscCall(PetscObjectReference((PetscObject)G)); 297 PetscCall(MatDestroy(&pcbddc->discretegradient)); 298 pcbddc->discretegradient = G; 299 pcbddc->nedorder = order > 0 ? order : -order; 300 pcbddc->nedfield = field; 301 pcbddc->nedglobal = global; 302 pcbddc->conforming = conforming; 303 PetscFunctionReturn(0); 304 } 305 306 /*@ 307 PCBDDCSetDiscreteGradient - Sets the discrete gradient 308 309 Collective on PC 310 311 Input Parameters: 312 + pc - the preconditioning context 313 . G - the discrete gradient matrix (should be in AIJ format) 314 . order - the order of the Nedelec space (1 for the lowest order) 315 . field - the field id of the Nedelec dofs (not used if the fields have not been specified) 316 . global - the type of global ordering for the rows of G 317 - conforming - whether the mesh is conforming or not 318 319 Level: advanced 320 321 Notes: 322 The discrete gradient matrix G is used to analyze the subdomain edges, and it should not contain any zero entry. 323 For variable order spaces, the order should be set to zero. 324 If global is true, the rows of G should be given in global ordering for the whole dofs; 325 if false, the ordering should be global for the Nedelec field. 326 In the latter case, it should hold gid[i] < gid[j] iff geid[i] < geid[j], with gid the global orderding for all the dofs 327 and geid the one for the Nedelec field. 328 329 .seealso: PCBDDC,PCBDDCSetDofsSplitting(),PCBDDCSetDofsSplittingLocal() 330 @*/ 331 PetscErrorCode PCBDDCSetDiscreteGradient(PC pc, Mat G, PetscInt order, PetscInt field, PetscBool global, PetscBool conforming) 332 { 333 PetscFunctionBegin; 334 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 335 PetscValidHeaderSpecific(G,MAT_CLASSID,2); 336 PetscValidLogicalCollectiveInt(pc,order,3); 337 PetscValidLogicalCollectiveInt(pc,field,4); 338 PetscValidLogicalCollectiveBool(pc,global,5); 339 PetscValidLogicalCollectiveBool(pc,conforming,6); 340 PetscCheckSameComm(pc,1,G,2); 341 PetscTryMethod(pc,"PCBDDCSetDiscreteGradient_C",(PC,Mat,PetscInt,PetscInt,PetscBool,PetscBool),(pc,G,order,field,global,conforming)); 342 PetscFunctionReturn(0); 343 } 344 345 static PetscErrorCode PCBDDCSetDivergenceMat_BDDC(PC pc, Mat divudotp, PetscBool trans, IS vl2l) 346 { 347 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 348 349 PetscFunctionBegin; 350 PetscCall(PetscObjectReference((PetscObject)divudotp)); 351 PetscCall(MatDestroy(&pcbddc->divudotp)); 352 pcbddc->divudotp = divudotp; 353 pcbddc->divudotp_trans = trans; 354 pcbddc->compute_nonetflux = PETSC_TRUE; 355 if (vl2l) { 356 PetscCall(PetscObjectReference((PetscObject)vl2l)); 357 PetscCall(ISDestroy(&pcbddc->divudotp_vl2l)); 358 pcbddc->divudotp_vl2l = vl2l; 359 } 360 PetscFunctionReturn(0); 361 } 362 363 /*@ 364 PCBDDCSetDivergenceMat - Sets the linear operator representing \int_\Omega \div {\bf u} \cdot p dx 365 366 Collective on PC 367 368 Input Parameters: 369 + pc - the preconditioning context 370 . divudotp - the matrix (must be of type MATIS) 371 . trans - if trans if false (resp. true), then pressures are in the test (trial) space and velocities are in the trial (test) space. 372 - vl2l - optional index set describing the local (wrt the local matrix in divudotp) to local (wrt the local matrix in the preconditioning matrix) map for the velocities 373 374 Level: advanced 375 376 Notes: 377 This auxiliary matrix is used to compute quadrature weights representing the net-flux across subdomain boundaries 378 If vl2l is NULL, the local ordering for velocities in divudotp should match that of the preconditioning matrix 379 380 .seealso: PCBDDC 381 @*/ 382 PetscErrorCode PCBDDCSetDivergenceMat(PC pc, Mat divudotp, PetscBool trans, IS vl2l) 383 { 384 PetscBool ismatis; 385 386 PetscFunctionBegin; 387 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 388 PetscValidHeaderSpecific(divudotp,MAT_CLASSID,2); 389 PetscCheckSameComm(pc,1,divudotp,2); 390 PetscValidLogicalCollectiveBool(pc,trans,3); 391 if (vl2l) PetscValidHeaderSpecific(vl2l,IS_CLASSID,4); 392 PetscCall(PetscObjectTypeCompare((PetscObject)divudotp,MATIS,&ismatis)); 393 PetscCheck(ismatis,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Divergence matrix needs to be of type MATIS"); 394 PetscTryMethod(pc,"PCBDDCSetDivergenceMat_C",(PC,Mat,PetscBool,IS),(pc,divudotp,trans,vl2l)); 395 PetscFunctionReturn(0); 396 } 397 398 static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change, PetscBool interior) 399 { 400 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 401 402 PetscFunctionBegin; 403 PetscCall(PetscObjectReference((PetscObject)change)); 404 PetscCall(MatDestroy(&pcbddc->user_ChangeOfBasisMatrix)); 405 pcbddc->user_ChangeOfBasisMatrix = change; 406 pcbddc->change_interior = interior; 407 PetscFunctionReturn(0); 408 } 409 /*@ 410 PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs 411 412 Collective on PC 413 414 Input Parameters: 415 + pc - the preconditioning context 416 . change - the change of basis matrix 417 - interior - whether or not the change of basis modifies interior dofs 418 419 Level: intermediate 420 421 Notes: 422 423 .seealso: PCBDDC 424 @*/ 425 PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change, PetscBool interior) 426 { 427 PetscFunctionBegin; 428 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 429 PetscValidHeaderSpecific(change,MAT_CLASSID,2); 430 PetscCheckSameComm(pc,1,change,2); 431 if (pc->mat) { 432 PetscInt rows_c,cols_c,rows,cols; 433 PetscCall(MatGetSize(pc->mat,&rows,&cols)); 434 PetscCall(MatGetSize(change,&rows_c,&cols_c)); 435 PetscCheck(rows_c == rows,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %D != %D",rows_c,rows); 436 PetscCheck(cols_c == cols,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %D != %D",cols_c,cols); 437 PetscCall(MatGetLocalSize(pc->mat,&rows,&cols)); 438 PetscCall(MatGetLocalSize(change,&rows_c,&cols_c)); 439 PetscCheck(rows_c == rows,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %D != %D",rows_c,rows); 440 PetscCheck(cols_c == cols,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %D != %D",cols_c,cols); 441 } 442 PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat,PetscBool),(pc,change,interior)); 443 PetscFunctionReturn(0); 444 } 445 446 static PetscErrorCode PCBDDCSetPrimalVerticesIS_BDDC(PC pc, IS PrimalVertices) 447 { 448 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 449 PetscBool isequal = PETSC_FALSE; 450 451 PetscFunctionBegin; 452 PetscCall(PetscObjectReference((PetscObject)PrimalVertices)); 453 if (pcbddc->user_primal_vertices) { 454 PetscCall(ISEqual(PrimalVertices,pcbddc->user_primal_vertices,&isequal)); 455 } 456 PetscCall(ISDestroy(&pcbddc->user_primal_vertices)); 457 PetscCall(ISDestroy(&pcbddc->user_primal_vertices_local)); 458 pcbddc->user_primal_vertices = PrimalVertices; 459 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 460 PetscFunctionReturn(0); 461 } 462 463 /*@ 464 PCBDDCSetPrimalVerticesIS - Set additional user defined primal vertices in PCBDDC 465 466 Collective 467 468 Input Parameters: 469 + pc - the preconditioning context 470 - PrimalVertices - index set of primal vertices in global numbering (can be empty) 471 472 Level: intermediate 473 474 Notes: 475 Any process can list any global node 476 477 .seealso: PCBDDC, PCBDDCGetPrimalVerticesIS(), PCBDDCSetPrimalVerticesLocalIS(), PCBDDCGetPrimalVerticesLocalIS() 478 @*/ 479 PetscErrorCode PCBDDCSetPrimalVerticesIS(PC pc, IS PrimalVertices) 480 { 481 PetscFunctionBegin; 482 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 483 PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 484 PetscCheckSameComm(pc,1,PrimalVertices,2); 485 PetscTryMethod(pc,"PCBDDCSetPrimalVerticesIS_C",(PC,IS),(pc,PrimalVertices)); 486 PetscFunctionReturn(0); 487 } 488 489 static PetscErrorCode PCBDDCGetPrimalVerticesIS_BDDC(PC pc, IS *is) 490 { 491 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 492 493 PetscFunctionBegin; 494 *is = pcbddc->user_primal_vertices; 495 PetscFunctionReturn(0); 496 } 497 498 /*@ 499 PCBDDCGetPrimalVerticesIS - Get user defined primal vertices set with PCBDDCSetPrimalVerticesIS() 500 501 Collective 502 503 Input Parameters: 504 . pc - the preconditioning context 505 506 Output Parameters: 507 . is - index set of primal vertices in global numbering (NULL if not set) 508 509 Level: intermediate 510 511 Notes: 512 513 .seealso: PCBDDC, PCBDDCSetPrimalVerticesIS(), PCBDDCSetPrimalVerticesLocalIS(), PCBDDCGetPrimalVerticesLocalIS() 514 @*/ 515 PetscErrorCode PCBDDCGetPrimalVerticesIS(PC pc, IS *is) 516 { 517 PetscFunctionBegin; 518 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 519 PetscValidPointer(is,2); 520 PetscUseMethod(pc,"PCBDDCGetPrimalVerticesIS_C",(PC,IS*),(pc,is)); 521 PetscFunctionReturn(0); 522 } 523 524 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 525 { 526 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 527 PetscBool isequal = PETSC_FALSE; 528 529 PetscFunctionBegin; 530 PetscCall(PetscObjectReference((PetscObject)PrimalVertices)); 531 if (pcbddc->user_primal_vertices_local) { 532 PetscCall(ISEqual(PrimalVertices,pcbddc->user_primal_vertices_local,&isequal)); 533 } 534 PetscCall(ISDestroy(&pcbddc->user_primal_vertices)); 535 PetscCall(ISDestroy(&pcbddc->user_primal_vertices_local)); 536 pcbddc->user_primal_vertices_local = PrimalVertices; 537 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 538 PetscFunctionReturn(0); 539 } 540 541 /*@ 542 PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 543 544 Collective 545 546 Input Parameters: 547 + pc - the preconditioning context 548 - PrimalVertices - index set of primal vertices in local numbering (can be empty) 549 550 Level: intermediate 551 552 Notes: 553 554 .seealso: PCBDDC, PCBDDCSetPrimalVerticesIS(), PCBDDCGetPrimalVerticesIS(), PCBDDCGetPrimalVerticesLocalIS() 555 @*/ 556 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 557 { 558 PetscFunctionBegin; 559 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 560 PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 561 PetscCheckSameComm(pc,1,PrimalVertices,2); 562 PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices)); 563 PetscFunctionReturn(0); 564 } 565 566 static PetscErrorCode PCBDDCGetPrimalVerticesLocalIS_BDDC(PC pc, IS *is) 567 { 568 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 569 570 PetscFunctionBegin; 571 *is = pcbddc->user_primal_vertices_local; 572 PetscFunctionReturn(0); 573 } 574 575 /*@ 576 PCBDDCGetPrimalVerticesLocalIS - Get user defined primal vertices set with PCBDDCSetPrimalVerticesLocalIS() 577 578 Collective 579 580 Input Parameters: 581 . pc - the preconditioning context 582 583 Output Parameters: 584 . is - index set of primal vertices in local numbering (NULL if not set) 585 586 Level: intermediate 587 588 Notes: 589 590 .seealso: PCBDDC, PCBDDCSetPrimalVerticesIS(), PCBDDCGetPrimalVerticesIS(), PCBDDCSetPrimalVerticesLocalIS() 591 @*/ 592 PetscErrorCode PCBDDCGetPrimalVerticesLocalIS(PC pc, IS *is) 593 { 594 PetscFunctionBegin; 595 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 596 PetscValidPointer(is,2); 597 PetscUseMethod(pc,"PCBDDCGetPrimalVerticesLocalIS_C",(PC,IS*),(pc,is)); 598 PetscFunctionReturn(0); 599 } 600 601 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 602 { 603 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 604 605 PetscFunctionBegin; 606 pcbddc->coarsening_ratio = k; 607 PetscFunctionReturn(0); 608 } 609 610 /*@ 611 PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 612 613 Logically collective on PC 614 615 Input Parameters: 616 + pc - the preconditioning context 617 - k - coarsening ratio (H/h at the coarser level) 618 619 Options Database Keys: 620 . -pc_bddc_coarsening_ratio <int> - Set coarsening ratio used in multilevel coarsening 621 622 Level: intermediate 623 624 Notes: 625 Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 626 627 .seealso: PCBDDC, PCBDDCSetLevels() 628 @*/ 629 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 630 { 631 PetscFunctionBegin; 632 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 633 PetscValidLogicalCollectiveInt(pc,k,2); 634 PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k)); 635 PetscFunctionReturn(0); 636 } 637 638 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 639 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 640 { 641 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 642 643 PetscFunctionBegin; 644 pcbddc->use_exact_dirichlet_trick = flg; 645 PetscFunctionReturn(0); 646 } 647 648 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 649 { 650 PetscFunctionBegin; 651 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 652 PetscValidLogicalCollectiveBool(pc,flg,2); 653 PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg)); 654 PetscFunctionReturn(0); 655 } 656 657 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 658 { 659 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 660 661 PetscFunctionBegin; 662 pcbddc->current_level = level; 663 PetscFunctionReturn(0); 664 } 665 666 PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 667 { 668 PetscFunctionBegin; 669 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 670 PetscValidLogicalCollectiveInt(pc,level,2); 671 PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level)); 672 PetscFunctionReturn(0); 673 } 674 675 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 676 { 677 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 678 679 PetscFunctionBegin; 680 PetscCheck(levels < PETSC_PCBDDC_MAXLEVELS,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of additional levels for BDDC is %d",PETSC_PCBDDC_MAXLEVELS-1); 681 pcbddc->max_levels = levels; 682 PetscFunctionReturn(0); 683 } 684 685 /*@ 686 PCBDDCSetLevels - Sets the maximum number of additional levels allowed for multilevel BDDC 687 688 Logically collective on PC 689 690 Input Parameters: 691 + pc - the preconditioning context 692 - levels - the maximum number of levels 693 694 Options Database Keys: 695 . -pc_bddc_levels <int> - Set maximum number of levels for multilevel 696 697 Level: intermediate 698 699 Notes: 700 The default value is 0, that gives the classical two-levels BDDC 701 702 .seealso: PCBDDC, PCBDDCSetCoarseningRatio() 703 @*/ 704 PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 705 { 706 PetscFunctionBegin; 707 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 708 PetscValidLogicalCollectiveInt(pc,levels,2); 709 PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels)); 710 PetscFunctionReturn(0); 711 } 712 713 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 714 { 715 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 716 PetscBool isequal = PETSC_FALSE; 717 718 PetscFunctionBegin; 719 PetscCall(PetscObjectReference((PetscObject)DirichletBoundaries)); 720 if (pcbddc->DirichletBoundaries) { 721 PetscCall(ISEqual(DirichletBoundaries,pcbddc->DirichletBoundaries,&isequal)); 722 } 723 /* last user setting takes precedence -> destroy any other customization */ 724 PetscCall(ISDestroy(&pcbddc->DirichletBoundariesLocal)); 725 PetscCall(ISDestroy(&pcbddc->DirichletBoundaries)); 726 pcbddc->DirichletBoundaries = DirichletBoundaries; 727 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 728 PetscFunctionReturn(0); 729 } 730 731 /*@ 732 PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 733 734 Collective 735 736 Input Parameters: 737 + pc - the preconditioning context 738 - DirichletBoundaries - parallel IS defining the Dirichlet boundaries 739 740 Level: intermediate 741 742 Notes: 743 Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node 744 745 .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns() 746 @*/ 747 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 748 { 749 PetscFunctionBegin; 750 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 751 PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 752 PetscCheckSameComm(pc,1,DirichletBoundaries,2); 753 PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries)); 754 PetscFunctionReturn(0); 755 } 756 757 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries) 758 { 759 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 760 PetscBool isequal = PETSC_FALSE; 761 762 PetscFunctionBegin; 763 PetscCall(PetscObjectReference((PetscObject)DirichletBoundaries)); 764 if (pcbddc->DirichletBoundariesLocal) { 765 PetscCall(ISEqual(DirichletBoundaries,pcbddc->DirichletBoundariesLocal,&isequal)); 766 } 767 /* last user setting takes precedence -> destroy any other customization */ 768 PetscCall(ISDestroy(&pcbddc->DirichletBoundariesLocal)); 769 PetscCall(ISDestroy(&pcbddc->DirichletBoundaries)); 770 pcbddc->DirichletBoundariesLocal = DirichletBoundaries; 771 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 772 PetscFunctionReturn(0); 773 } 774 775 /*@ 776 PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering. 777 778 Collective 779 780 Input Parameters: 781 + pc - the preconditioning context 782 - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering) 783 784 Level: intermediate 785 786 Notes: 787 788 .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns() 789 @*/ 790 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries) 791 { 792 PetscFunctionBegin; 793 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 794 PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 795 PetscCheckSameComm(pc,1,DirichletBoundaries,2); 796 PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries)); 797 PetscFunctionReturn(0); 798 } 799 800 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 801 { 802 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 803 PetscBool isequal = PETSC_FALSE; 804 805 PetscFunctionBegin; 806 PetscCall(PetscObjectReference((PetscObject)NeumannBoundaries)); 807 if (pcbddc->NeumannBoundaries) { 808 PetscCall(ISEqual(NeumannBoundaries,pcbddc->NeumannBoundaries,&isequal)); 809 } 810 /* last user setting takes precedence -> destroy any other customization */ 811 PetscCall(ISDestroy(&pcbddc->NeumannBoundariesLocal)); 812 PetscCall(ISDestroy(&pcbddc->NeumannBoundaries)); 813 pcbddc->NeumannBoundaries = NeumannBoundaries; 814 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 815 PetscFunctionReturn(0); 816 } 817 818 /*@ 819 PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 820 821 Collective 822 823 Input Parameters: 824 + pc - the preconditioning context 825 - NeumannBoundaries - parallel IS defining the Neumann boundaries 826 827 Level: intermediate 828 829 Notes: 830 Any process can list any global node 831 832 .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal() 833 @*/ 834 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 835 { 836 PetscFunctionBegin; 837 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 838 PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 839 PetscCheckSameComm(pc,1,NeumannBoundaries,2); 840 PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries)); 841 PetscFunctionReturn(0); 842 } 843 844 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries) 845 { 846 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 847 PetscBool isequal = PETSC_FALSE; 848 849 PetscFunctionBegin; 850 PetscCall(PetscObjectReference((PetscObject)NeumannBoundaries)); 851 if (pcbddc->NeumannBoundariesLocal) { 852 PetscCall(ISEqual(NeumannBoundaries,pcbddc->NeumannBoundariesLocal,&isequal)); 853 } 854 /* last user setting takes precedence -> destroy any other customization */ 855 PetscCall(ISDestroy(&pcbddc->NeumannBoundariesLocal)); 856 PetscCall(ISDestroy(&pcbddc->NeumannBoundaries)); 857 pcbddc->NeumannBoundariesLocal = NeumannBoundaries; 858 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 859 PetscFunctionReturn(0); 860 } 861 862 /*@ 863 PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering. 864 865 Collective 866 867 Input Parameters: 868 + pc - the preconditioning context 869 - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering) 870 871 Level: intermediate 872 873 Notes: 874 875 .seealso: PCBDDC, PCBDDCSetNeumannBoundaries() 876 @*/ 877 PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries) 878 { 879 PetscFunctionBegin; 880 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 881 PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 882 PetscCheckSameComm(pc,1,NeumannBoundaries,2); 883 PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries)); 884 PetscFunctionReturn(0); 885 } 886 887 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 888 { 889 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 890 891 PetscFunctionBegin; 892 *DirichletBoundaries = pcbddc->DirichletBoundaries; 893 PetscFunctionReturn(0); 894 } 895 896 /*@ 897 PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries 898 899 Collective 900 901 Input Parameters: 902 . pc - the preconditioning context 903 904 Output Parameters: 905 . DirichletBoundaries - index set defining the Dirichlet boundaries 906 907 Level: intermediate 908 909 Notes: 910 The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries 911 912 .seealso: PCBDDC 913 @*/ 914 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 915 { 916 PetscFunctionBegin; 917 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 918 PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries)); 919 PetscFunctionReturn(0); 920 } 921 922 static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries) 923 { 924 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 925 926 PetscFunctionBegin; 927 *DirichletBoundaries = pcbddc->DirichletBoundariesLocal; 928 PetscFunctionReturn(0); 929 } 930 931 /*@ 932 PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering) 933 934 Collective 935 936 Input Parameters: 937 . pc - the preconditioning context 938 939 Output Parameters: 940 . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 941 942 Level: intermediate 943 944 Notes: 945 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). 946 In the latter case, the IS will be available after PCSetUp. 947 948 .seealso: PCBDDC 949 @*/ 950 PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries) 951 { 952 PetscFunctionBegin; 953 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 954 PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries)); 955 PetscFunctionReturn(0); 956 } 957 958 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 959 { 960 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 961 962 PetscFunctionBegin; 963 *NeumannBoundaries = pcbddc->NeumannBoundaries; 964 PetscFunctionReturn(0); 965 } 966 967 /*@ 968 PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries 969 970 Collective 971 972 Input Parameters: 973 . pc - the preconditioning context 974 975 Output Parameters: 976 . NeumannBoundaries - index set defining the Neumann boundaries 977 978 Level: intermediate 979 980 Notes: 981 The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries 982 983 .seealso: PCBDDC 984 @*/ 985 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 986 { 987 PetscFunctionBegin; 988 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 989 PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries)); 990 PetscFunctionReturn(0); 991 } 992 993 static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries) 994 { 995 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 996 997 PetscFunctionBegin; 998 *NeumannBoundaries = pcbddc->NeumannBoundariesLocal; 999 PetscFunctionReturn(0); 1000 } 1001 1002 /*@ 1003 PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering) 1004 1005 Collective 1006 1007 Input Parameters: 1008 . pc - the preconditioning context 1009 1010 Output Parameters: 1011 . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 1012 1013 Level: intermediate 1014 1015 Notes: 1016 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). 1017 In the latter case, the IS will be available after PCSetUp. 1018 1019 .seealso: PCBDDC 1020 @*/ 1021 PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries) 1022 { 1023 PetscFunctionBegin; 1024 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1025 PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries)); 1026 PetscFunctionReturn(0); 1027 } 1028 1029 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 1030 { 1031 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1032 PCBDDCGraph mat_graph = pcbddc->mat_graph; 1033 PetscBool same_data = PETSC_FALSE; 1034 1035 PetscFunctionBegin; 1036 if (!nvtxs) { 1037 if (copymode == PETSC_OWN_POINTER) { 1038 PetscCall(PetscFree(xadj)); 1039 PetscCall(PetscFree(adjncy)); 1040 } 1041 PetscCall(PCBDDCGraphResetCSR(mat_graph)); 1042 PetscFunctionReturn(0); 1043 } 1044 if (mat_graph->nvtxs == nvtxs && mat_graph->freecsr) { /* we own the data */ 1045 if (mat_graph->xadj == xadj && mat_graph->adjncy == adjncy) same_data = PETSC_TRUE; 1046 if (!same_data && mat_graph->xadj[nvtxs] == xadj[nvtxs]) { 1047 PetscCall(PetscArraycmp(xadj,mat_graph->xadj,nvtxs+1,&same_data)); 1048 if (same_data) { 1049 PetscCall(PetscArraycmp(adjncy,mat_graph->adjncy,xadj[nvtxs],&same_data)); 1050 } 1051 } 1052 } 1053 if (!same_data) { 1054 /* free old CSR */ 1055 PetscCall(PCBDDCGraphResetCSR(mat_graph)); 1056 /* get CSR into graph structure */ 1057 if (copymode == PETSC_COPY_VALUES) { 1058 PetscCall(PetscMalloc1(nvtxs+1,&mat_graph->xadj)); 1059 PetscCall(PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy)); 1060 PetscCall(PetscArraycpy(mat_graph->xadj,xadj,nvtxs+1)); 1061 PetscCall(PetscArraycpy(mat_graph->adjncy,adjncy,xadj[nvtxs])); 1062 mat_graph->freecsr = PETSC_TRUE; 1063 } else if (copymode == PETSC_OWN_POINTER) { 1064 mat_graph->xadj = (PetscInt*)xadj; 1065 mat_graph->adjncy = (PetscInt*)adjncy; 1066 mat_graph->freecsr = PETSC_TRUE; 1067 } else if (copymode == PETSC_USE_POINTER) { 1068 mat_graph->xadj = (PetscInt*)xadj; 1069 mat_graph->adjncy = (PetscInt*)adjncy; 1070 mat_graph->freecsr = PETSC_FALSE; 1071 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %D",copymode); 1072 mat_graph->nvtxs_csr = nvtxs; 1073 pcbddc->recompute_topography = PETSC_TRUE; 1074 } 1075 PetscFunctionReturn(0); 1076 } 1077 1078 /*@ 1079 PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local degrees of freedom. 1080 1081 Not collective 1082 1083 Input Parameters: 1084 + pc - the preconditioning context. 1085 . nvtxs - number of local vertices of the graph (i.e., the number of local dofs). 1086 . xadj, adjncy - the connectivity of the dofs in CSR format. 1087 - copymode - supported modes are PETSC_COPY_VALUES, PETSC_USE_POINTER or PETSC_OWN_POINTER. 1088 1089 Level: intermediate 1090 1091 Notes: 1092 A dof is considered connected with all local dofs if xadj[dof+1]-xadj[dof] == 1 and adjncy[xadj[dof]] is negative. 1093 1094 .seealso: PCBDDC,PetscCopyMode 1095 @*/ 1096 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 1097 { 1098 void (*f)(void) = NULL; 1099 1100 PetscFunctionBegin; 1101 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1102 if (nvtxs) { 1103 PetscValidIntPointer(xadj,3); 1104 if (xadj[nvtxs]) PetscValidIntPointer(adjncy,4); 1105 } 1106 PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode)); 1107 /* free arrays if PCBDDC is not the PC type */ 1108 PetscCall(PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f)); 1109 if (!f && copymode == PETSC_OWN_POINTER) { 1110 PetscCall(PetscFree(xadj)); 1111 PetscCall(PetscFree(adjncy)); 1112 } 1113 PetscFunctionReturn(0); 1114 } 1115 1116 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 1117 { 1118 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1119 PetscInt i; 1120 PetscBool isequal = PETSC_FALSE; 1121 1122 PetscFunctionBegin; 1123 if (pcbddc->n_ISForDofsLocal == n_is) { 1124 for (i=0;i<n_is;i++) { 1125 PetscBool isequalt; 1126 PetscCall(ISEqual(ISForDofs[i],pcbddc->ISForDofsLocal[i],&isequalt)); 1127 if (!isequalt) break; 1128 } 1129 if (i == n_is) isequal = PETSC_TRUE; 1130 } 1131 for (i=0;i<n_is;i++) { 1132 PetscCall(PetscObjectReference((PetscObject)ISForDofs[i])); 1133 } 1134 /* Destroy ISes if they were already set */ 1135 for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 1136 PetscCall(ISDestroy(&pcbddc->ISForDofsLocal[i])); 1137 } 1138 PetscCall(PetscFree(pcbddc->ISForDofsLocal)); 1139 /* last user setting takes precedence -> destroy any other customization */ 1140 for (i=0;i<pcbddc->n_ISForDofs;i++) { 1141 PetscCall(ISDestroy(&pcbddc->ISForDofs[i])); 1142 } 1143 PetscCall(PetscFree(pcbddc->ISForDofs)); 1144 pcbddc->n_ISForDofs = 0; 1145 /* allocate space then set */ 1146 if (n_is) { 1147 PetscCall(PetscMalloc1(n_is,&pcbddc->ISForDofsLocal)); 1148 } 1149 for (i=0;i<n_is;i++) { 1150 pcbddc->ISForDofsLocal[i] = ISForDofs[i]; 1151 } 1152 pcbddc->n_ISForDofsLocal = n_is; 1153 if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 1154 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 1155 PetscFunctionReturn(0); 1156 } 1157 1158 /*@ 1159 PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix 1160 1161 Collective 1162 1163 Input Parameters: 1164 + pc - the preconditioning context 1165 . n_is - number of index sets defining the fields 1166 - ISForDofs - array of IS describing the fields in local ordering 1167 1168 Level: intermediate 1169 1170 Notes: 1171 n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field. 1172 1173 .seealso: PCBDDC 1174 @*/ 1175 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[]) 1176 { 1177 PetscInt i; 1178 1179 PetscFunctionBegin; 1180 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1181 PetscValidLogicalCollectiveInt(pc,n_is,2); 1182 for (i=0;i<n_is;i++) { 1183 PetscCheckSameComm(pc,1,ISForDofs[i],3); 1184 PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 1185 } 1186 PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs)); 1187 PetscFunctionReturn(0); 1188 } 1189 1190 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 1191 { 1192 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1193 PetscInt i; 1194 PetscBool isequal = PETSC_FALSE; 1195 1196 PetscFunctionBegin; 1197 if (pcbddc->n_ISForDofs == n_is) { 1198 for (i=0;i<n_is;i++) { 1199 PetscBool isequalt; 1200 PetscCall(ISEqual(ISForDofs[i],pcbddc->ISForDofs[i],&isequalt)); 1201 if (!isequalt) break; 1202 } 1203 if (i == n_is) isequal = PETSC_TRUE; 1204 } 1205 for (i=0;i<n_is;i++) { 1206 PetscCall(PetscObjectReference((PetscObject)ISForDofs[i])); 1207 } 1208 /* Destroy ISes if they were already set */ 1209 for (i=0;i<pcbddc->n_ISForDofs;i++) { 1210 PetscCall(ISDestroy(&pcbddc->ISForDofs[i])); 1211 } 1212 PetscCall(PetscFree(pcbddc->ISForDofs)); 1213 /* last user setting takes precedence -> destroy any other customization */ 1214 for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 1215 PetscCall(ISDestroy(&pcbddc->ISForDofsLocal[i])); 1216 } 1217 PetscCall(PetscFree(pcbddc->ISForDofsLocal)); 1218 pcbddc->n_ISForDofsLocal = 0; 1219 /* allocate space then set */ 1220 if (n_is) { 1221 PetscCall(PetscMalloc1(n_is,&pcbddc->ISForDofs)); 1222 } 1223 for (i=0;i<n_is;i++) { 1224 pcbddc->ISForDofs[i] = ISForDofs[i]; 1225 } 1226 pcbddc->n_ISForDofs = n_is; 1227 if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 1228 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 1229 PetscFunctionReturn(0); 1230 } 1231 1232 /*@ 1233 PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix 1234 1235 Collective 1236 1237 Input Parameters: 1238 + pc - the preconditioning context 1239 . n_is - number of index sets defining the fields 1240 - ISForDofs - array of IS describing the fields in global ordering 1241 1242 Level: intermediate 1243 1244 Notes: 1245 Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field. 1246 1247 .seealso: PCBDDC 1248 @*/ 1249 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 1250 { 1251 PetscInt i; 1252 1253 PetscFunctionBegin; 1254 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1255 PetscValidLogicalCollectiveInt(pc,n_is,2); 1256 for (i=0;i<n_is;i++) { 1257 PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 1258 PetscCheckSameComm(pc,1,ISForDofs[i],3); 1259 } 1260 PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs)); 1261 PetscFunctionReturn(0); 1262 } 1263 1264 /* 1265 PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 1266 guess if a transformation of basis approach has been selected. 1267 1268 Input Parameter: 1269 + pc - the preconditioner context 1270 1271 Application Interface Routine: PCPreSolve() 1272 1273 Notes: 1274 The interface routine PCPreSolve() is not usually called directly by 1275 the user, but instead is called by KSPSolve(). 1276 */ 1277 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1278 { 1279 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1280 PC_IS *pcis = (PC_IS*)(pc->data); 1281 Vec used_vec; 1282 PetscBool iscg, save_rhs = PETSC_TRUE, benign_correction_computed; 1283 1284 PetscFunctionBegin; 1285 /* if we are working with CG, one dirichlet solve can be avoided during Krylov iterations */ 1286 if (ksp) { 1287 PetscCall(PetscObjectTypeCompareAny((PetscObject)ksp,&iscg,KSPCG,KSPGROPPCG,KSPPIPECG,KSPPIPELCG,KSPPIPECGRR,"")); 1288 if (pcbddc->benign_apply_coarse_only || pcbddc->switch_static || !iscg || pc->mat != pc->pmat) { 1289 PetscCall(PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE)); 1290 } 1291 } 1292 if (pcbddc->benign_apply_coarse_only || pcbddc->switch_static || pc->mat != pc->pmat) { 1293 PetscCall(PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE)); 1294 } 1295 1296 /* Creates parallel work vectors used in presolve */ 1297 if (!pcbddc->original_rhs) { 1298 PetscCall(VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs)); 1299 } 1300 if (!pcbddc->temp_solution) { 1301 PetscCall(VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution)); 1302 } 1303 1304 pcbddc->temp_solution_used = PETSC_FALSE; 1305 if (x) { 1306 PetscCall(PetscObjectReference((PetscObject)x)); 1307 used_vec = x; 1308 } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */ 1309 PetscCall(PetscObjectReference((PetscObject)pcbddc->temp_solution)); 1310 used_vec = pcbddc->temp_solution; 1311 PetscCall(VecSet(used_vec,0.0)); 1312 pcbddc->temp_solution_used = PETSC_TRUE; 1313 PetscCall(VecCopy(rhs,pcbddc->original_rhs)); 1314 save_rhs = PETSC_FALSE; 1315 pcbddc->eliminate_dirdofs = PETSC_TRUE; 1316 } 1317 1318 /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */ 1319 if (ksp) { 1320 /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */ 1321 PetscCall(KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero)); 1322 if (!pcbddc->ksp_guess_nonzero) { 1323 PetscCall(VecSet(used_vec,0.0)); 1324 } 1325 } 1326 1327 pcbddc->rhs_change = PETSC_FALSE; 1328 /* Take into account zeroed rows -> change rhs and store solution removed */ 1329 if (rhs && pcbddc->eliminate_dirdofs) { 1330 IS dirIS = NULL; 1331 1332 /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */ 1333 PetscCall(PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS)); 1334 if (dirIS) { 1335 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 1336 PetscInt dirsize,i,*is_indices; 1337 PetscScalar *array_x; 1338 const PetscScalar *array_diagonal; 1339 1340 PetscCall(MatGetDiagonal(pc->pmat,pcis->vec1_global)); 1341 PetscCall(VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global)); 1342 PetscCall(VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD)); 1343 PetscCall(VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD)); 1344 PetscCall(VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD)); 1345 PetscCall(VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD)); 1346 PetscCall(ISGetLocalSize(dirIS,&dirsize)); 1347 PetscCall(VecGetArray(pcis->vec1_N,&array_x)); 1348 PetscCall(VecGetArrayRead(pcis->vec2_N,&array_diagonal)); 1349 PetscCall(ISGetIndices(dirIS,(const PetscInt**)&is_indices)); 1350 for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 1351 PetscCall(ISRestoreIndices(dirIS,(const PetscInt**)&is_indices)); 1352 PetscCall(VecRestoreArrayRead(pcis->vec2_N,&array_diagonal)); 1353 PetscCall(VecRestoreArray(pcis->vec1_N,&array_x)); 1354 PetscCall(VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE)); 1355 PetscCall(VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE)); 1356 pcbddc->rhs_change = PETSC_TRUE; 1357 PetscCall(ISDestroy(&dirIS)); 1358 } 1359 } 1360 1361 /* remove the computed solution or the initial guess from the rhs */ 1362 if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero)) { 1363 /* save the original rhs */ 1364 if (save_rhs) { 1365 PetscCall(VecSwap(rhs,pcbddc->original_rhs)); 1366 save_rhs = PETSC_FALSE; 1367 } 1368 pcbddc->rhs_change = PETSC_TRUE; 1369 PetscCall(VecScale(used_vec,-1.0)); 1370 PetscCall(MatMultAdd(pc->mat,used_vec,pcbddc->original_rhs,rhs)); 1371 PetscCall(VecScale(used_vec,-1.0)); 1372 PetscCall(VecCopy(used_vec,pcbddc->temp_solution)); 1373 pcbddc->temp_solution_used = PETSC_TRUE; 1374 if (ksp) { 1375 PetscCall(KSPSetInitialGuessNonzero(ksp,PETSC_FALSE)); 1376 } 1377 } 1378 PetscCall(VecDestroy(&used_vec)); 1379 1380 /* compute initial vector in benign space if needed 1381 and remove non-benign solution from the rhs */ 1382 benign_correction_computed = PETSC_FALSE; 1383 if (rhs && pcbddc->benign_compute_correction && (pcbddc->benign_have_null || pcbddc->benign_apply_coarse_only)) { 1384 /* compute u^*_h using ideas similar to those in Xuemin Tu's PhD thesis (see Section 4.8.1) 1385 Recursively apply BDDC in the multilevel case */ 1386 if (!pcbddc->benign_vec) { 1387 PetscCall(VecDuplicate(rhs,&pcbddc->benign_vec)); 1388 } 1389 /* keep applying coarse solver unless we no longer have benign subdomains */ 1390 pcbddc->benign_apply_coarse_only = pcbddc->benign_have_null ? PETSC_TRUE : PETSC_FALSE; 1391 if (!pcbddc->benign_skip_correction) { 1392 PetscCall(PCApply_BDDC(pc,rhs,pcbddc->benign_vec)); 1393 benign_correction_computed = PETSC_TRUE; 1394 if (pcbddc->temp_solution_used) { 1395 PetscCall(VecAXPY(pcbddc->temp_solution,1.0,pcbddc->benign_vec)); 1396 } 1397 PetscCall(VecScale(pcbddc->benign_vec,-1.0)); 1398 /* store the original rhs if not done earlier */ 1399 if (save_rhs) { 1400 PetscCall(VecSwap(rhs,pcbddc->original_rhs)); 1401 } 1402 if (pcbddc->rhs_change) { 1403 PetscCall(MatMultAdd(pc->mat,pcbddc->benign_vec,rhs,rhs)); 1404 } else { 1405 PetscCall(MatMultAdd(pc->mat,pcbddc->benign_vec,pcbddc->original_rhs,rhs)); 1406 } 1407 pcbddc->rhs_change = PETSC_TRUE; 1408 } 1409 pcbddc->benign_apply_coarse_only = PETSC_FALSE; 1410 } else { 1411 PetscCall(VecDestroy(&pcbddc->benign_vec)); 1412 } 1413 1414 /* dbg output */ 1415 if (pcbddc->dbg_flag && benign_correction_computed) { 1416 Vec v; 1417 1418 PetscCall(VecDuplicate(pcis->vec1_global,&v)); 1419 if (pcbddc->ChangeOfBasisMatrix) { 1420 PetscCall(MatMultTranspose(pcbddc->ChangeOfBasisMatrix,rhs,v)); 1421 } else { 1422 PetscCall(VecCopy(rhs,v)); 1423 } 1424 PetscCall(PCBDDCBenignGetOrSetP0(pc,v,PETSC_TRUE)); 1425 PetscCall(PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"LEVEL %D: is the correction benign?\n",pcbddc->current_level)); 1426 PetscCall(PetscScalarView(pcbddc->benign_n,pcbddc->benign_p0,pcbddc->dbg_viewer)); 1427 PetscCall(PetscViewerFlush(pcbddc->dbg_viewer)); 1428 PetscCall(VecDestroy(&v)); 1429 } 1430 1431 /* set initial guess if using PCG */ 1432 pcbddc->exact_dirichlet_trick_app = PETSC_FALSE; 1433 if (x && pcbddc->use_exact_dirichlet_trick) { 1434 PetscCall(VecSet(x,0.0)); 1435 if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) { 1436 if (benign_correction_computed) { /* we have already saved the changed rhs */ 1437 PetscCall(VecLockReadPop(pcis->vec1_global)); 1438 } else { 1439 PetscCall(MatMultTranspose(pcbddc->ChangeOfBasisMatrix,rhs,pcis->vec1_global)); 1440 } 1441 PetscCall(VecScatterBegin(pcis->global_to_D,pcis->vec1_global,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD)); 1442 PetscCall(VecScatterEnd(pcis->global_to_D,pcis->vec1_global,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD)); 1443 } else { 1444 PetscCall(VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD)); 1445 PetscCall(VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD)); 1446 } 1447 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0],pc,0,0,0)); 1448 PetscCall(KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D)); 1449 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0],pc,0,0,0)); 1450 PetscCall(KSPCheckSolve(pcbddc->ksp_D,pc,pcis->vec2_D)); 1451 if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) { 1452 PetscCall(VecSet(pcis->vec1_global,0.)); 1453 PetscCall(VecScatterBegin(pcis->global_to_D,pcis->vec2_D,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE)); 1454 PetscCall(VecScatterEnd(pcis->global_to_D,pcis->vec2_D,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE)); 1455 PetscCall(MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_global,x)); 1456 } else { 1457 PetscCall(VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE)); 1458 PetscCall(VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE)); 1459 } 1460 if (ksp) { 1461 PetscCall(KSPSetInitialGuessNonzero(ksp,PETSC_TRUE)); 1462 } 1463 pcbddc->exact_dirichlet_trick_app = PETSC_TRUE; 1464 } else if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior && benign_correction_computed && pcbddc->use_exact_dirichlet_trick) { 1465 PetscCall(VecLockReadPop(pcis->vec1_global)); 1466 } 1467 PetscFunctionReturn(0); 1468 } 1469 1470 /* 1471 PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 1472 approach has been selected. Also, restores rhs to its original state. 1473 1474 Input Parameter: 1475 + pc - the preconditioner context 1476 1477 Application Interface Routine: PCPostSolve() 1478 1479 Notes: 1480 The interface routine PCPostSolve() is not usually called directly by 1481 the user, but instead is called by KSPSolve(). 1482 */ 1483 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1484 { 1485 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1486 1487 PetscFunctionBegin; 1488 /* add solution removed in presolve */ 1489 if (x && pcbddc->rhs_change) { 1490 if (pcbddc->temp_solution_used) { 1491 PetscCall(VecAXPY(x,1.0,pcbddc->temp_solution)); 1492 } else if (pcbddc->benign_compute_correction && pcbddc->benign_vec) { 1493 PetscCall(VecAXPY(x,-1.0,pcbddc->benign_vec)); 1494 } 1495 /* restore to original state (not for FETI-DP) */ 1496 if (ksp) pcbddc->temp_solution_used = PETSC_FALSE; 1497 } 1498 1499 /* restore rhs to its original state (not needed for FETI-DP) */ 1500 if (rhs && pcbddc->rhs_change) { 1501 PetscCall(VecSwap(rhs,pcbddc->original_rhs)); 1502 pcbddc->rhs_change = PETSC_FALSE; 1503 } 1504 /* restore ksp guess state */ 1505 if (ksp) { 1506 PetscCall(KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero)); 1507 /* reset flag for exact dirichlet trick */ 1508 pcbddc->exact_dirichlet_trick_app = PETSC_FALSE; 1509 } 1510 PetscFunctionReturn(0); 1511 } 1512 1513 /* 1514 PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 1515 by setting data structures and options. 1516 1517 Input Parameter: 1518 + pc - the preconditioner context 1519 1520 Application Interface Routine: PCSetUp() 1521 1522 Notes: 1523 The interface routine PCSetUp() is not usually called directly by 1524 the user, but instead is called by PCApply() if necessary. 1525 */ 1526 PetscErrorCode PCSetUp_BDDC(PC pc) 1527 { 1528 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 1529 PCBDDCSubSchurs sub_schurs; 1530 Mat_IS* matis; 1531 MatNullSpace nearnullspace; 1532 Mat lA; 1533 IS lP,zerodiag = NULL; 1534 PetscInt nrows,ncols; 1535 PetscMPIInt size; 1536 PetscBool computesubschurs; 1537 PetscBool computeconstraintsmatrix; 1538 PetscBool new_nearnullspace_provided,ismatis,rl; 1539 1540 PetscFunctionBegin; 1541 PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis)); 1542 PetscCheck(ismatis,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS"); 1543 PetscCall(MatGetSize(pc->pmat,&nrows,&ncols)); 1544 PetscCheck(nrows == ncols,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix"); 1545 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 1546 1547 matis = (Mat_IS*)pc->pmat->data; 1548 /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 1549 /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1550 Also, BDDC builds its own KSP for the Dirichlet problem */ 1551 rl = pcbddc->recompute_topography; 1552 if (!pc->setupcalled || pc->flag == DIFFERENT_NONZERO_PATTERN) rl = PETSC_TRUE; 1553 PetscCall(MPIU_Allreduce(&rl,&pcbddc->recompute_topography,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc))); 1554 if (pcbddc->recompute_topography) { 1555 pcbddc->graphanalyzed = PETSC_FALSE; 1556 computeconstraintsmatrix = PETSC_TRUE; 1557 } else { 1558 computeconstraintsmatrix = PETSC_FALSE; 1559 } 1560 1561 /* check parameters' compatibility */ 1562 if (!pcbddc->use_deluxe_scaling) pcbddc->deluxe_zerorows = PETSC_FALSE; 1563 pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold[0] != 0.0 || pcbddc->adaptive_threshold[1] != 0.0); 1564 pcbddc->use_deluxe_scaling = (PetscBool)(pcbddc->use_deluxe_scaling && size > 1); 1565 pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_selection && size > 1); 1566 pcbddc->adaptive_userdefined = (PetscBool)(pcbddc->adaptive_selection && pcbddc->adaptive_userdefined); 1567 if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE; 1568 1569 computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling); 1570 1571 /* activate all connected components if the netflux has been requested */ 1572 if (pcbddc->compute_nonetflux) { 1573 pcbddc->use_vertices = PETSC_TRUE; 1574 pcbddc->use_edges = PETSC_TRUE; 1575 pcbddc->use_faces = PETSC_TRUE; 1576 } 1577 1578 /* Get stdout for dbg */ 1579 if (pcbddc->dbg_flag) { 1580 if (!pcbddc->dbg_viewer) { 1581 pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)); 1582 } 1583 PetscCall(PetscViewerASCIIPushSynchronized(pcbddc->dbg_viewer)); 1584 PetscCall(PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level)); 1585 } 1586 1587 /* process topology information */ 1588 PetscCall(PetscLogEventBegin(PC_BDDC_Topology[pcbddc->current_level],pc,0,0,0)); 1589 if (pcbddc->recompute_topography) { 1590 PetscCall(PCBDDCComputeLocalTopologyInfo(pc)); 1591 if (pcbddc->discretegradient) { 1592 PetscCall(PCBDDCNedelecSupport(pc)); 1593 } 1594 } 1595 if (pcbddc->corner_selected) pcbddc->use_vertices = PETSC_TRUE; 1596 1597 /* change basis if requested by the user */ 1598 if (pcbddc->user_ChangeOfBasisMatrix) { 1599 /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */ 1600 pcbddc->use_change_of_basis = PETSC_FALSE; 1601 PetscCall(PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix)); 1602 } else { 1603 PetscCall(MatDestroy(&pcbddc->local_mat)); 1604 PetscCall(PetscObjectReference((PetscObject)matis->A)); 1605 pcbddc->local_mat = matis->A; 1606 } 1607 1608 /* 1609 Compute change of basis on local pressures (aka zerodiag dofs) with the benign trick 1610 This should come earlier then PCISSetUp for extracting the correct subdomain matrices 1611 */ 1612 PetscCall(PCBDDCBenignShellMat(pc,PETSC_TRUE)); 1613 if (pcbddc->benign_saddle_point) { 1614 PC_IS* pcis = (PC_IS*)pc->data; 1615 1616 if (pcbddc->user_ChangeOfBasisMatrix || pcbddc->use_change_of_basis || !computesubschurs) pcbddc->benign_change_explicit = PETSC_TRUE; 1617 /* detect local saddle point and change the basis in pcbddc->local_mat */ 1618 PetscCall(PCBDDCBenignDetectSaddlePoint(pc,(PetscBool)(!pcbddc->recompute_topography),&zerodiag)); 1619 /* pop B0 mat from local mat */ 1620 PetscCall(PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE)); 1621 /* give pcis a hint to not reuse submatrices during PCISCreate */ 1622 if (pc->flag == SAME_NONZERO_PATTERN && pcis->reusesubmatrices == PETSC_TRUE) { 1623 if (pcbddc->benign_n && (pcbddc->benign_change_explicit || pcbddc->dbg_flag)) { 1624 pcis->reusesubmatrices = PETSC_FALSE; 1625 } else { 1626 pcis->reusesubmatrices = PETSC_TRUE; 1627 } 1628 } else { 1629 pcis->reusesubmatrices = PETSC_FALSE; 1630 } 1631 } 1632 1633 /* propagate relevant information */ 1634 if (matis->A->symmetric_set) { 1635 PetscCall(MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric)); 1636 } 1637 if (matis->A->spd_set) { 1638 PetscCall(MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd)); 1639 } 1640 1641 /* Set up all the "iterative substructuring" common block without computing solvers */ 1642 { 1643 Mat temp_mat; 1644 1645 temp_mat = matis->A; 1646 matis->A = pcbddc->local_mat; 1647 PetscCall(PCISSetUp(pc,PETSC_TRUE,PETSC_FALSE)); 1648 pcbddc->local_mat = matis->A; 1649 matis->A = temp_mat; 1650 } 1651 1652 /* Analyze interface */ 1653 if (!pcbddc->graphanalyzed) { 1654 PetscCall(PCBDDCAnalyzeInterface(pc)); 1655 computeconstraintsmatrix = PETSC_TRUE; 1656 if (pcbddc->adaptive_selection && !pcbddc->use_deluxe_scaling && !pcbddc->mat_graph->twodim) { 1657 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute the adaptive primal space for a problem with 3D edges without deluxe scaling"); 1658 } 1659 if (pcbddc->compute_nonetflux) { 1660 MatNullSpace nnfnnsp; 1661 1662 PetscCheck(pcbddc->divudotp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Missing divudotp operator"); 1663 PetscCall(PCBDDCComputeNoNetFlux(pc->pmat,pcbddc->divudotp,pcbddc->divudotp_trans,pcbddc->divudotp_vl2l,pcbddc->mat_graph,&nnfnnsp)); 1664 /* TODO what if a nearnullspace is already attached? */ 1665 if (nnfnnsp) { 1666 PetscCall(MatSetNearNullSpace(pc->pmat,nnfnnsp)); 1667 PetscCall(MatNullSpaceDestroy(&nnfnnsp)); 1668 } 1669 } 1670 } 1671 PetscCall(PetscLogEventEnd(PC_BDDC_Topology[pcbddc->current_level],pc,0,0,0)); 1672 1673 /* check existence of a divergence free extension, i.e. 1674 b(v_I,p_0) = 0 for all v_I (raise error if not). 1675 Also, check that PCBDDCBenignGetOrSetP0 works */ 1676 if (pcbddc->benign_saddle_point && pcbddc->dbg_flag > 1) { 1677 PetscCall(PCBDDCBenignCheck(pc,zerodiag)); 1678 } 1679 PetscCall(ISDestroy(&zerodiag)); 1680 1681 /* Setup local dirichlet solver ksp_D and sub_schurs solvers */ 1682 if (computesubschurs && pcbddc->recompute_topography) { 1683 PetscCall(PCBDDCInitSubSchurs(pc)); 1684 } 1685 /* SetUp Scaling operator (scaling matrices could be needed in SubSchursSetUp)*/ 1686 if (!pcbddc->use_deluxe_scaling) { 1687 PetscCall(PCBDDCScalingSetUp(pc)); 1688 } 1689 1690 /* finish setup solvers and do adaptive selection of constraints */ 1691 sub_schurs = pcbddc->sub_schurs; 1692 if (sub_schurs && sub_schurs->schur_explicit) { 1693 if (computesubschurs) { 1694 PetscCall(PCBDDCSetUpSubSchurs(pc)); 1695 } 1696 PetscCall(PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE)); 1697 } else { 1698 PetscCall(PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE)); 1699 if (computesubschurs) { 1700 PetscCall(PCBDDCSetUpSubSchurs(pc)); 1701 } 1702 } 1703 if (pcbddc->adaptive_selection) { 1704 PetscCall(PCBDDCAdaptiveSelection(pc)); 1705 computeconstraintsmatrix = PETSC_TRUE; 1706 } 1707 1708 /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1709 new_nearnullspace_provided = PETSC_FALSE; 1710 PetscCall(MatGetNearNullSpace(pc->pmat,&nearnullspace)); 1711 if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1712 if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1713 new_nearnullspace_provided = PETSC_TRUE; 1714 } else { 1715 /* determine if the two nullspaces are different (should be lightweight) */ 1716 if (nearnullspace != pcbddc->onearnullspace) { 1717 new_nearnullspace_provided = PETSC_TRUE; 1718 } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1719 PetscInt i; 1720 const Vec *nearnullvecs; 1721 PetscObjectState state; 1722 PetscInt nnsp_size; 1723 PetscCall(MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs)); 1724 for (i=0;i<nnsp_size;i++) { 1725 PetscCall(PetscObjectStateGet((PetscObject)nearnullvecs[i],&state)); 1726 if (pcbddc->onearnullvecs_state[i] != state) { 1727 new_nearnullspace_provided = PETSC_TRUE; 1728 break; 1729 } 1730 } 1731 } 1732 } 1733 } else { 1734 if (!nearnullspace) { /* both nearnullspaces are null */ 1735 new_nearnullspace_provided = PETSC_FALSE; 1736 } else { /* nearnullspace attached later */ 1737 new_nearnullspace_provided = PETSC_TRUE; 1738 } 1739 } 1740 1741 /* Setup constraints and related work vectors */ 1742 /* reset primal space flags */ 1743 PetscCall(PetscLogEventBegin(PC_BDDC_LocalWork[pcbddc->current_level],pc,0,0,0)); 1744 pcbddc->new_primal_space = PETSC_FALSE; 1745 pcbddc->new_primal_space_local = PETSC_FALSE; 1746 if (computeconstraintsmatrix || new_nearnullspace_provided) { 1747 /* It also sets the primal space flags */ 1748 PetscCall(PCBDDCConstraintsSetUp(pc)); 1749 } 1750 /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1751 PetscCall(PCBDDCSetUpLocalWorkVectors(pc)); 1752 1753 if (pcbddc->use_change_of_basis) { 1754 PC_IS *pcis = (PC_IS*)(pc->data); 1755 1756 PetscCall(PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix)); 1757 if (pcbddc->benign_change) { 1758 PetscCall(MatDestroy(&pcbddc->benign_B0)); 1759 /* pop B0 from pcbddc->local_mat */ 1760 PetscCall(PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE)); 1761 } 1762 /* get submatrices */ 1763 PetscCall(MatDestroy(&pcis->A_IB)); 1764 PetscCall(MatDestroy(&pcis->A_BI)); 1765 PetscCall(MatDestroy(&pcis->A_BB)); 1766 PetscCall(MatCreateSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB)); 1767 PetscCall(MatCreateSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB)); 1768 PetscCall(MatCreateSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI)); 1769 /* set flag in pcis to not reuse submatrices during PCISCreate */ 1770 pcis->reusesubmatrices = PETSC_FALSE; 1771 } else if (!pcbddc->user_ChangeOfBasisMatrix && !pcbddc->benign_change) { 1772 PetscCall(MatDestroy(&pcbddc->local_mat)); 1773 PetscCall(PetscObjectReference((PetscObject)matis->A)); 1774 pcbddc->local_mat = matis->A; 1775 } 1776 1777 /* interface pressure block row for B_C */ 1778 PetscCall(PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_lP" ,(PetscObject*)&lP)); 1779 PetscCall(PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_lA" ,(PetscObject*)&lA)); 1780 if (lA && lP) { 1781 PC_IS* pcis = (PC_IS*)pc->data; 1782 Mat B_BI,B_BB,Bt_BI,Bt_BB; 1783 PetscBool issym; 1784 PetscCall(MatIsSymmetric(lA,PETSC_SMALL,&issym)); 1785 if (issym) { 1786 PetscCall(MatCreateSubMatrix(lA,lP,pcis->is_I_local,MAT_INITIAL_MATRIX,&B_BI)); 1787 PetscCall(MatCreateSubMatrix(lA,lP,pcis->is_B_local,MAT_INITIAL_MATRIX,&B_BB)); 1788 PetscCall(MatCreateTranspose(B_BI,&Bt_BI)); 1789 PetscCall(MatCreateTranspose(B_BB,&Bt_BB)); 1790 } else { 1791 PetscCall(MatCreateSubMatrix(lA,lP,pcis->is_I_local,MAT_INITIAL_MATRIX,&B_BI)); 1792 PetscCall(MatCreateSubMatrix(lA,lP,pcis->is_B_local,MAT_INITIAL_MATRIX,&B_BB)); 1793 PetscCall(MatCreateSubMatrix(lA,pcis->is_I_local,lP,MAT_INITIAL_MATRIX,&Bt_BI)); 1794 PetscCall(MatCreateSubMatrix(lA,pcis->is_B_local,lP,MAT_INITIAL_MATRIX,&Bt_BB)); 1795 } 1796 PetscCall(PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_B_BI",(PetscObject)B_BI)); 1797 PetscCall(PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_B_BB",(PetscObject)B_BB)); 1798 PetscCall(PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_Bt_BI",(PetscObject)Bt_BI)); 1799 PetscCall(PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_Bt_BB",(PetscObject)Bt_BB)); 1800 PetscCall(MatDestroy(&B_BI)); 1801 PetscCall(MatDestroy(&B_BB)); 1802 PetscCall(MatDestroy(&Bt_BI)); 1803 PetscCall(MatDestroy(&Bt_BB)); 1804 } 1805 PetscCall(PetscLogEventEnd(PC_BDDC_LocalWork[pcbddc->current_level],pc,0,0,0)); 1806 1807 /* SetUp coarse and local Neumann solvers */ 1808 PetscCall(PCBDDCSetUpSolvers(pc)); 1809 /* SetUp Scaling operator */ 1810 if (pcbddc->use_deluxe_scaling) { 1811 PetscCall(PCBDDCScalingSetUp(pc)); 1812 } 1813 1814 /* mark topography as done */ 1815 pcbddc->recompute_topography = PETSC_FALSE; 1816 1817 /* wrap pcis->A_IB and pcis->A_BI if we did not change explicitly the variables on the pressures */ 1818 PetscCall(PCBDDCBenignShellMat(pc,PETSC_FALSE)); 1819 1820 if (pcbddc->dbg_flag) { 1821 PetscCall(PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level)); 1822 PetscCall(PetscViewerASCIIPopSynchronized(pcbddc->dbg_viewer)); 1823 } 1824 PetscFunctionReturn(0); 1825 } 1826 1827 /* 1828 PCApply_BDDC - Applies the BDDC operator to a vector. 1829 1830 Input Parameters: 1831 + pc - the preconditioner context 1832 - r - input vector (global) 1833 1834 Output Parameter: 1835 . z - output vector (global) 1836 1837 Application Interface Routine: PCApply() 1838 */ 1839 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 1840 { 1841 PC_IS *pcis = (PC_IS*)(pc->data); 1842 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1843 Mat lA = NULL; 1844 PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B; 1845 const PetscScalar one = 1.0; 1846 const PetscScalar m_one = -1.0; 1847 const PetscScalar zero = 0.0; 1848 /* This code is similar to that provided in nn.c for PCNN 1849 NN interface preconditioner changed to BDDC 1850 Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */ 1851 1852 PetscFunctionBegin; 1853 PetscCall(PetscCitationsRegister(citation,&cited)); 1854 if (pcbddc->switch_static) { 1855 PetscCall(MatISGetLocalMat(pc->useAmat ? pc->mat : pc->pmat,&lA)); 1856 } 1857 1858 if (pcbddc->ChangeOfBasisMatrix) { 1859 Vec swap; 1860 1861 PetscCall(MatMultTranspose(pcbddc->ChangeOfBasisMatrix,r,pcbddc->work_change)); 1862 swap = pcbddc->work_change; 1863 pcbddc->work_change = r; 1864 r = swap; 1865 /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */ 1866 if (pcbddc->benign_apply_coarse_only && pcbddc->use_exact_dirichlet_trick && pcbddc->change_interior) { 1867 PetscCall(VecCopy(r,pcis->vec1_global)); 1868 PetscCall(VecLockReadPush(pcis->vec1_global)); 1869 } 1870 } 1871 if (pcbddc->benign_have_null) { /* get p0 from r */ 1872 PetscCall(PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE)); 1873 } 1874 if (pcbddc->interface_extension == PC_BDDC_INTERFACE_EXT_DIRICHLET && !pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) { 1875 PetscCall(VecCopy(r,z)); 1876 /* First Dirichlet solve */ 1877 PetscCall(VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD)); 1878 PetscCall(VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD)); 1879 /* 1880 Assembling right hand side for BDDC operator 1881 - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE) 1882 - pcis->vec1_B the interface part of the global vector z 1883 */ 1884 if (n_D) { 1885 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0],pc,0,0,0)); 1886 PetscCall(KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D)); 1887 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0],pc,0,0,0)); 1888 PetscCall(KSPCheckSolve(pcbddc->ksp_D,pc,pcis->vec2_D)); 1889 PetscCall(VecScale(pcis->vec2_D,m_one)); 1890 if (pcbddc->switch_static) { 1891 PetscCall(VecSet(pcis->vec1_N,0.)); 1892 PetscCall(VecScatterBegin(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 1893 PetscCall(VecScatterEnd(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 1894 if (!pcbddc->switch_static_change) { 1895 PetscCall(MatMult(lA,pcis->vec1_N,pcis->vec2_N)); 1896 } else { 1897 PetscCall(MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N)); 1898 PetscCall(MatMult(lA,pcis->vec2_N,pcis->vec1_N)); 1899 PetscCall(MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N)); 1900 } 1901 PetscCall(VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD)); 1902 PetscCall(VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD)); 1903 PetscCall(VecScatterBegin(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD)); 1904 PetscCall(VecScatterEnd(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD)); 1905 } else { 1906 PetscCall(MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B)); 1907 } 1908 } else { 1909 PetscCall(VecSet(pcis->vec1_B,zero)); 1910 } 1911 PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE)); 1912 PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE)); 1913 PetscCall(PCBDDCScalingRestriction(pc,z,pcis->vec1_B)); 1914 } else { 1915 if (!pcbddc->benign_apply_coarse_only) { 1916 PetscCall(PCBDDCScalingRestriction(pc,r,pcis->vec1_B)); 1917 } 1918 } 1919 if (pcbddc->interface_extension == PC_BDDC_INTERFACE_EXT_LUMP) { 1920 PetscCheck(pcbddc->switch_static,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"You forgot to pass -pc_bddc_switch_static"); 1921 PetscCall(VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD)); 1922 PetscCall(VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD)); 1923 } 1924 1925 /* Apply interface preconditioner 1926 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1927 PetscCall(PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE)); 1928 1929 /* Apply transpose of partition of unity operator */ 1930 PetscCall(PCBDDCScalingExtension(pc,pcis->vec1_B,z)); 1931 if (pcbddc->interface_extension == PC_BDDC_INTERFACE_EXT_LUMP) { 1932 PetscCall(VecScatterBegin(pcis->global_to_D,pcis->vec1_D,z,INSERT_VALUES,SCATTER_REVERSE)); 1933 PetscCall(VecScatterEnd(pcis->global_to_D,pcis->vec1_D,z,INSERT_VALUES,SCATTER_REVERSE)); 1934 PetscFunctionReturn(0); 1935 } 1936 /* Second Dirichlet solve and assembling of output */ 1937 PetscCall(VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD)); 1938 PetscCall(VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD)); 1939 if (n_B) { 1940 if (pcbddc->switch_static) { 1941 PetscCall(VecScatterBegin(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 1942 PetscCall(VecScatterEnd(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 1943 PetscCall(VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 1944 PetscCall(VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 1945 if (!pcbddc->switch_static_change) { 1946 PetscCall(MatMult(lA,pcis->vec1_N,pcis->vec2_N)); 1947 } else { 1948 PetscCall(MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N)); 1949 PetscCall(MatMult(lA,pcis->vec2_N,pcis->vec1_N)); 1950 PetscCall(MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N)); 1951 } 1952 PetscCall(VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD)); 1953 PetscCall(VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD)); 1954 } else { 1955 PetscCall(MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D)); 1956 } 1957 } else if (pcbddc->switch_static) { /* n_B is zero */ 1958 if (!pcbddc->switch_static_change) { 1959 PetscCall(MatMult(lA,pcis->vec1_D,pcis->vec3_D)); 1960 } else { 1961 PetscCall(MatMult(pcbddc->switch_static_change,pcis->vec1_D,pcis->vec1_N)); 1962 PetscCall(MatMult(lA,pcis->vec1_N,pcis->vec2_N)); 1963 PetscCall(MatMultTranspose(pcbddc->switch_static_change,pcis->vec2_N,pcis->vec3_D)); 1964 } 1965 } 1966 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0],pc,0,0,0)); 1967 PetscCall(KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D)); 1968 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0],pc,0,0,0)); 1969 PetscCall(KSPCheckSolve(pcbddc->ksp_D,pc,pcis->vec4_D)); 1970 1971 if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) { 1972 if (pcbddc->switch_static) { 1973 PetscCall(VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D)); 1974 } else { 1975 PetscCall(VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D)); 1976 } 1977 PetscCall(VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE)); 1978 PetscCall(VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE)); 1979 } else { 1980 if (pcbddc->switch_static) { 1981 PetscCall(VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D)); 1982 } else { 1983 PetscCall(VecScale(pcis->vec4_D,m_one)); 1984 } 1985 PetscCall(VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE)); 1986 PetscCall(VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE)); 1987 } 1988 if (pcbddc->benign_have_null) { /* set p0 (computed in PCBDDCApplyInterface) */ 1989 if (pcbddc->benign_apply_coarse_only) { 1990 PetscCall(PetscArrayzero(pcbddc->benign_p0,pcbddc->benign_n)); 1991 } 1992 PetscCall(PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE)); 1993 } 1994 1995 if (pcbddc->ChangeOfBasisMatrix) { 1996 pcbddc->work_change = r; 1997 PetscCall(VecCopy(z,pcbddc->work_change)); 1998 PetscCall(MatMult(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,z)); 1999 } 2000 PetscFunctionReturn(0); 2001 } 2002 2003 /* 2004 PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector. 2005 2006 Input Parameters: 2007 + pc - the preconditioner context 2008 - r - input vector (global) 2009 2010 Output Parameter: 2011 . z - output vector (global) 2012 2013 Application Interface Routine: PCApplyTranspose() 2014 */ 2015 PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z) 2016 { 2017 PC_IS *pcis = (PC_IS*)(pc->data); 2018 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 2019 Mat lA = NULL; 2020 PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B; 2021 const PetscScalar one = 1.0; 2022 const PetscScalar m_one = -1.0; 2023 const PetscScalar zero = 0.0; 2024 2025 PetscFunctionBegin; 2026 PetscCall(PetscCitationsRegister(citation,&cited)); 2027 if (pcbddc->switch_static) { 2028 PetscCall(MatISGetLocalMat(pc->useAmat ? pc->mat : pc->pmat,&lA)); 2029 } 2030 if (pcbddc->ChangeOfBasisMatrix) { 2031 Vec swap; 2032 2033 PetscCall(MatMultTranspose(pcbddc->ChangeOfBasisMatrix,r,pcbddc->work_change)); 2034 swap = pcbddc->work_change; 2035 pcbddc->work_change = r; 2036 r = swap; 2037 /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */ 2038 if (pcbddc->benign_apply_coarse_only && pcbddc->exact_dirichlet_trick_app && pcbddc->change_interior) { 2039 PetscCall(VecCopy(r,pcis->vec1_global)); 2040 PetscCall(VecLockReadPush(pcis->vec1_global)); 2041 } 2042 } 2043 if (pcbddc->benign_have_null) { /* get p0 from r */ 2044 PetscCall(PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE)); 2045 } 2046 if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) { 2047 PetscCall(VecCopy(r,z)); 2048 /* First Dirichlet solve */ 2049 PetscCall(VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD)); 2050 PetscCall(VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD)); 2051 /* 2052 Assembling right hand side for BDDC operator 2053 - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE) 2054 - pcis->vec1_B the interface part of the global vector z 2055 */ 2056 if (n_D) { 2057 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0],pc,0,0,0)); 2058 PetscCall(KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D)); 2059 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0],pc,0,0,0)); 2060 PetscCall(KSPCheckSolve(pcbddc->ksp_D,pc,pcis->vec2_D)); 2061 PetscCall(VecScale(pcis->vec2_D,m_one)); 2062 if (pcbddc->switch_static) { 2063 PetscCall(VecSet(pcis->vec1_N,0.)); 2064 PetscCall(VecScatterBegin(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 2065 PetscCall(VecScatterEnd(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 2066 if (!pcbddc->switch_static_change) { 2067 PetscCall(MatMultTranspose(lA,pcis->vec1_N,pcis->vec2_N)); 2068 } else { 2069 PetscCall(MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N)); 2070 PetscCall(MatMultTranspose(lA,pcis->vec2_N,pcis->vec1_N)); 2071 PetscCall(MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N)); 2072 } 2073 PetscCall(VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD)); 2074 PetscCall(VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD)); 2075 PetscCall(VecScatterBegin(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD)); 2076 PetscCall(VecScatterEnd(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD)); 2077 } else { 2078 PetscCall(MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B)); 2079 } 2080 } else { 2081 PetscCall(VecSet(pcis->vec1_B,zero)); 2082 } 2083 PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE)); 2084 PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE)); 2085 PetscCall(PCBDDCScalingRestriction(pc,z,pcis->vec1_B)); 2086 } else { 2087 PetscCall(PCBDDCScalingRestriction(pc,r,pcis->vec1_B)); 2088 } 2089 2090 /* Apply interface preconditioner 2091 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 2092 PetscCall(PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE)); 2093 2094 /* Apply transpose of partition of unity operator */ 2095 PetscCall(PCBDDCScalingExtension(pc,pcis->vec1_B,z)); 2096 2097 /* Second Dirichlet solve and assembling of output */ 2098 PetscCall(VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD)); 2099 PetscCall(VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD)); 2100 if (n_B) { 2101 if (pcbddc->switch_static) { 2102 PetscCall(VecScatterBegin(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 2103 PetscCall(VecScatterEnd(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 2104 PetscCall(VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 2105 PetscCall(VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 2106 if (!pcbddc->switch_static_change) { 2107 PetscCall(MatMultTranspose(lA,pcis->vec1_N,pcis->vec2_N)); 2108 } else { 2109 PetscCall(MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N)); 2110 PetscCall(MatMultTranspose(lA,pcis->vec2_N,pcis->vec1_N)); 2111 PetscCall(MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N)); 2112 } 2113 PetscCall(VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD)); 2114 PetscCall(VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD)); 2115 } else { 2116 PetscCall(MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D)); 2117 } 2118 } else if (pcbddc->switch_static) { /* n_B is zero */ 2119 if (!pcbddc->switch_static_change) { 2120 PetscCall(MatMultTranspose(lA,pcis->vec1_D,pcis->vec3_D)); 2121 } else { 2122 PetscCall(MatMult(pcbddc->switch_static_change,pcis->vec1_D,pcis->vec1_N)); 2123 PetscCall(MatMultTranspose(lA,pcis->vec1_N,pcis->vec2_N)); 2124 PetscCall(MatMultTranspose(pcbddc->switch_static_change,pcis->vec2_N,pcis->vec3_D)); 2125 } 2126 } 2127 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0],pc,0,0,0)); 2128 PetscCall(KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D)); 2129 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0],pc,0,0,0)); 2130 PetscCall(KSPCheckSolve(pcbddc->ksp_D,pc,pcis->vec4_D)); 2131 if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) { 2132 if (pcbddc->switch_static) { 2133 PetscCall(VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D)); 2134 } else { 2135 PetscCall(VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D)); 2136 } 2137 PetscCall(VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE)); 2138 PetscCall(VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE)); 2139 } else { 2140 if (pcbddc->switch_static) { 2141 PetscCall(VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D)); 2142 } else { 2143 PetscCall(VecScale(pcis->vec4_D,m_one)); 2144 } 2145 PetscCall(VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE)); 2146 PetscCall(VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE)); 2147 } 2148 if (pcbddc->benign_have_null) { /* set p0 (computed in PCBDDCApplyInterface) */ 2149 PetscCall(PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE)); 2150 } 2151 if (pcbddc->ChangeOfBasisMatrix) { 2152 pcbddc->work_change = r; 2153 PetscCall(VecCopy(z,pcbddc->work_change)); 2154 PetscCall(MatMult(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,z)); 2155 } 2156 PetscFunctionReturn(0); 2157 } 2158 2159 PetscErrorCode PCReset_BDDC(PC pc) 2160 { 2161 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2162 PC_IS *pcis = (PC_IS*)pc->data; 2163 KSP kspD,kspR,kspC; 2164 2165 PetscFunctionBegin; 2166 /* free BDDC custom data */ 2167 PetscCall(PCBDDCResetCustomization(pc)); 2168 /* destroy objects related to topography */ 2169 PetscCall(PCBDDCResetTopography(pc)); 2170 /* destroy objects for scaling operator */ 2171 PetscCall(PCBDDCScalingDestroy(pc)); 2172 /* free solvers stuff */ 2173 PetscCall(PCBDDCResetSolvers(pc)); 2174 /* free global vectors needed in presolve */ 2175 PetscCall(VecDestroy(&pcbddc->temp_solution)); 2176 PetscCall(VecDestroy(&pcbddc->original_rhs)); 2177 /* free data created by PCIS */ 2178 PetscCall(PCISDestroy(pc)); 2179 2180 /* restore defaults */ 2181 kspD = pcbddc->ksp_D; 2182 kspR = pcbddc->ksp_R; 2183 kspC = pcbddc->coarse_ksp; 2184 PetscCall(PetscMemzero(pc->data,sizeof(*pcbddc))); 2185 pcis->n_neigh = -1; 2186 pcis->scaling_factor = 1.0; 2187 pcis->reusesubmatrices = PETSC_TRUE; 2188 pcbddc->use_local_adj = PETSC_TRUE; 2189 pcbddc->use_vertices = PETSC_TRUE; 2190 pcbddc->use_edges = PETSC_TRUE; 2191 pcbddc->symmetric_primal = PETSC_TRUE; 2192 pcbddc->vertex_size = 1; 2193 pcbddc->recompute_topography = PETSC_TRUE; 2194 pcbddc->coarse_size = -1; 2195 pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 2196 pcbddc->coarsening_ratio = 8; 2197 pcbddc->coarse_eqs_per_proc = 1; 2198 pcbddc->benign_compute_correction = PETSC_TRUE; 2199 pcbddc->nedfield = -1; 2200 pcbddc->nedglobal = PETSC_TRUE; 2201 pcbddc->graphmaxcount = PETSC_MAX_INT; 2202 pcbddc->sub_schurs_layers = -1; 2203 pcbddc->ksp_D = kspD; 2204 pcbddc->ksp_R = kspR; 2205 pcbddc->coarse_ksp = kspC; 2206 PetscFunctionReturn(0); 2207 } 2208 2209 PetscErrorCode PCDestroy_BDDC(PC pc) 2210 { 2211 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2212 2213 PetscFunctionBegin; 2214 PetscCall(PCReset_BDDC(pc)); 2215 PetscCall(KSPDestroy(&pcbddc->ksp_D)); 2216 PetscCall(KSPDestroy(&pcbddc->ksp_R)); 2217 PetscCall(KSPDestroy(&pcbddc->coarse_ksp)); 2218 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDiscreteGradient_C",NULL)); 2219 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDivergenceMat_C",NULL)); 2220 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL)); 2221 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL)); 2222 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",NULL)); 2223 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL)); 2224 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL)); 2225 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL)); 2226 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL)); 2227 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL)); 2228 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL)); 2229 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL)); 2230 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL)); 2231 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL)); 2232 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL)); 2233 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL)); 2234 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL)); 2235 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL)); 2236 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL)); 2237 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL)); 2238 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL)); 2239 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL)); 2240 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL)); 2241 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL)); 2242 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL)); 2243 PetscCall(PetscFree(pc->data)); 2244 PetscFunctionReturn(0); 2245 } 2246 2247 static PetscErrorCode PCSetCoordinates_BDDC(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords) 2248 { 2249 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 2250 PCBDDCGraph mat_graph = pcbddc->mat_graph; 2251 2252 PetscFunctionBegin; 2253 PetscCall(PetscFree(mat_graph->coords)); 2254 PetscCall(PetscMalloc1(nloc*dim,&mat_graph->coords)); 2255 PetscCall(PetscArraycpy(mat_graph->coords,coords,nloc*dim)); 2256 mat_graph->cnloc = nloc; 2257 mat_graph->cdim = dim; 2258 mat_graph->cloc = PETSC_FALSE; 2259 /* flg setup */ 2260 pcbddc->recompute_topography = PETSC_TRUE; 2261 pcbddc->corner_selected = PETSC_FALSE; 2262 PetscFunctionReturn(0); 2263 } 2264 2265 static PetscErrorCode PCPreSolveChangeRHS_BDDC(PC pc, PetscBool* change) 2266 { 2267 PetscFunctionBegin; 2268 *change = PETSC_TRUE; 2269 PetscFunctionReturn(0); 2270 } 2271 2272 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 2273 { 2274 FETIDPMat_ctx mat_ctx; 2275 Vec work; 2276 PC_IS* pcis; 2277 PC_BDDC* pcbddc; 2278 2279 PetscFunctionBegin; 2280 PetscCall(MatShellGetContext(fetidp_mat,&mat_ctx)); 2281 pcis = (PC_IS*)mat_ctx->pc->data; 2282 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 2283 2284 PetscCall(VecSet(fetidp_flux_rhs,0.0)); 2285 /* copy rhs since we may change it during PCPreSolve_BDDC */ 2286 if (!pcbddc->original_rhs) { 2287 PetscCall(VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs)); 2288 } 2289 if (mat_ctx->rhs_flip) { 2290 PetscCall(VecPointwiseMult(pcbddc->original_rhs,standard_rhs,mat_ctx->rhs_flip)); 2291 } else { 2292 PetscCall(VecCopy(standard_rhs,pcbddc->original_rhs)); 2293 } 2294 if (mat_ctx->g2g_p) { 2295 /* interface pressure rhs */ 2296 PetscCall(VecScatterBegin(mat_ctx->g2g_p,fetidp_flux_rhs,pcbddc->original_rhs,INSERT_VALUES,SCATTER_REVERSE)); 2297 PetscCall(VecScatterEnd(mat_ctx->g2g_p,fetidp_flux_rhs,pcbddc->original_rhs,INSERT_VALUES,SCATTER_REVERSE)); 2298 PetscCall(VecScatterBegin(mat_ctx->g2g_p,standard_rhs,fetidp_flux_rhs,INSERT_VALUES,SCATTER_FORWARD)); 2299 PetscCall(VecScatterEnd(mat_ctx->g2g_p,standard_rhs,fetidp_flux_rhs,INSERT_VALUES,SCATTER_FORWARD)); 2300 if (!mat_ctx->rhs_flip) { 2301 PetscCall(VecScale(fetidp_flux_rhs,-1.)); 2302 } 2303 } 2304 /* 2305 change of basis for physical rhs if needed 2306 It also changes the rhs in case of dirichlet boundaries 2307 */ 2308 PetscCall(PCPreSolve_BDDC(mat_ctx->pc,NULL,pcbddc->original_rhs,NULL)); 2309 if (pcbddc->ChangeOfBasisMatrix) { 2310 PetscCall(MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcbddc->original_rhs,pcbddc->work_change)); 2311 work = pcbddc->work_change; 2312 } else { 2313 work = pcbddc->original_rhs; 2314 } 2315 /* store vectors for computation of fetidp final solution */ 2316 PetscCall(VecScatterBegin(pcis->global_to_D,work,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD)); 2317 PetscCall(VecScatterEnd(pcis->global_to_D,work,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD)); 2318 /* scale rhs since it should be unassembled */ 2319 /* TODO use counter scaling? (also below) */ 2320 PetscCall(VecScatterBegin(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD)); 2321 PetscCall(VecScatterEnd(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD)); 2322 /* Apply partition of unity */ 2323 PetscCall(VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B)); 2324 /* PetscCall(PCBDDCScalingRestriction(mat_ctx->pc,work,mat_ctx->temp_solution_B)); */ 2325 if (!pcbddc->switch_static) { 2326 /* compute partially subassembled Schur complement right-hand side */ 2327 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0],mat_ctx->pc,0,0,0)); 2328 PetscCall(KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D)); 2329 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0],mat_ctx->pc,0,0,0)); 2330 /* Cannot propagate up error in KSPSolve() because there is no access to the PC */ 2331 PetscCall(MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B)); 2332 PetscCall(VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B)); 2333 PetscCall(VecSet(work,0.0)); 2334 PetscCall(VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,work,ADD_VALUES,SCATTER_REVERSE)); 2335 PetscCall(VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,work,ADD_VALUES,SCATTER_REVERSE)); 2336 /* PetscCall(PCBDDCScalingRestriction(mat_ctx->pc,work,mat_ctx->temp_solution_B)); */ 2337 PetscCall(VecScatterBegin(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD)); 2338 PetscCall(VecScatterEnd(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD)); 2339 PetscCall(VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B)); 2340 } 2341 /* BDDC rhs */ 2342 PetscCall(VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B)); 2343 if (pcbddc->switch_static) { 2344 PetscCall(VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D)); 2345 } 2346 /* apply BDDC */ 2347 PetscCall(PetscArrayzero(pcbddc->benign_p0,pcbddc->benign_n)); 2348 PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE)); 2349 PetscCall(PetscArrayzero(pcbddc->benign_p0,pcbddc->benign_n)); 2350 2351 /* Application of B_delta and assembling of rhs for fetidp fluxes */ 2352 PetscCall(MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local)); 2353 PetscCall(VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD)); 2354 PetscCall(VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD)); 2355 /* Add contribution to interface pressures */ 2356 if (mat_ctx->l2g_p) { 2357 PetscCall(MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP)); 2358 if (pcbddc->switch_static) { 2359 PetscCall(MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP)); 2360 } 2361 PetscCall(VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD)); 2362 PetscCall(VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD)); 2363 } 2364 PetscFunctionReturn(0); 2365 } 2366 2367 /*@ 2368 PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side 2369 2370 Collective 2371 2372 Input Parameters: 2373 + fetidp_mat - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators 2374 - standard_rhs - the right-hand side of the original linear system 2375 2376 Output Parameters: 2377 . fetidp_flux_rhs - the right-hand side for the FETI-DP linear system 2378 2379 Level: developer 2380 2381 Notes: 2382 2383 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution 2384 @*/ 2385 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 2386 { 2387 FETIDPMat_ctx mat_ctx; 2388 2389 PetscFunctionBegin; 2390 PetscValidHeaderSpecific(fetidp_mat,MAT_CLASSID,1); 2391 PetscValidHeaderSpecific(standard_rhs,VEC_CLASSID,2); 2392 PetscValidHeaderSpecific(fetidp_flux_rhs,VEC_CLASSID,3); 2393 PetscCall(MatShellGetContext(fetidp_mat,&mat_ctx)); 2394 PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs)); 2395 PetscFunctionReturn(0); 2396 } 2397 2398 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 2399 { 2400 FETIDPMat_ctx mat_ctx; 2401 PC_IS* pcis; 2402 PC_BDDC* pcbddc; 2403 Vec work; 2404 2405 PetscFunctionBegin; 2406 PetscCall(MatShellGetContext(fetidp_mat,&mat_ctx)); 2407 pcis = (PC_IS*)mat_ctx->pc->data; 2408 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 2409 2410 /* apply B_delta^T */ 2411 PetscCall(VecSet(pcis->vec1_B,0.)); 2412 PetscCall(VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE)); 2413 PetscCall(VecScatterEnd(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE)); 2414 PetscCall(MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B)); 2415 if (mat_ctx->l2g_p) { 2416 PetscCall(VecScatterBegin(mat_ctx->l2g_p,fetidp_flux_sol,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE)); 2417 PetscCall(VecScatterEnd(mat_ctx->l2g_p,fetidp_flux_sol,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE)); 2418 PetscCall(MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B)); 2419 } 2420 2421 /* compute rhs for BDDC application */ 2422 PetscCall(VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B)); 2423 if (pcbddc->switch_static) { 2424 PetscCall(VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D)); 2425 if (mat_ctx->l2g_p) { 2426 PetscCall(VecScale(mat_ctx->vP,-1.)); 2427 PetscCall(MatMultAdd(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D,pcis->vec1_D)); 2428 } 2429 } 2430 2431 /* apply BDDC */ 2432 PetscCall(PetscArrayzero(pcbddc->benign_p0,pcbddc->benign_n)); 2433 PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE)); 2434 2435 /* put values into global vector */ 2436 if (pcbddc->ChangeOfBasisMatrix) work = pcbddc->work_change; 2437 else work = standard_sol; 2438 PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_B,work,INSERT_VALUES,SCATTER_REVERSE)); 2439 PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_B,work,INSERT_VALUES,SCATTER_REVERSE)); 2440 if (!pcbddc->switch_static) { 2441 /* compute values into the interior if solved for the partially subassembled Schur complement */ 2442 PetscCall(MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D)); 2443 PetscCall(VecAYPX(pcis->vec1_D,-1.0,mat_ctx->temp_solution_D)); 2444 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0],mat_ctx->pc,0,0,0)); 2445 PetscCall(KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec1_D)); 2446 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0],mat_ctx->pc,0,0,0)); 2447 /* Cannot propagate up error in KSPSolve() because there is no access to the PC */ 2448 } 2449 2450 PetscCall(VecScatterBegin(pcis->global_to_D,pcis->vec1_D,work,INSERT_VALUES,SCATTER_REVERSE)); 2451 PetscCall(VecScatterEnd(pcis->global_to_D,pcis->vec1_D,work,INSERT_VALUES,SCATTER_REVERSE)); 2452 /* add p0 solution to final solution */ 2453 PetscCall(PCBDDCBenignGetOrSetP0(mat_ctx->pc,work,PETSC_FALSE)); 2454 if (pcbddc->ChangeOfBasisMatrix) { 2455 PetscCall(MatMult(pcbddc->ChangeOfBasisMatrix,work,standard_sol)); 2456 } 2457 PetscCall(PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol)); 2458 if (mat_ctx->g2g_p) { 2459 PetscCall(VecScatterBegin(mat_ctx->g2g_p,fetidp_flux_sol,standard_sol,INSERT_VALUES,SCATTER_REVERSE)); 2460 PetscCall(VecScatterEnd(mat_ctx->g2g_p,fetidp_flux_sol,standard_sol,INSERT_VALUES,SCATTER_REVERSE)); 2461 } 2462 PetscFunctionReturn(0); 2463 } 2464 2465 static PetscErrorCode PCView_BDDCIPC(PC pc, PetscViewer viewer) 2466 { 2467 BDDCIPC_ctx bddcipc_ctx; 2468 PetscBool isascii; 2469 2470 PetscFunctionBegin; 2471 PetscCall(PCShellGetContext(pc,&bddcipc_ctx)); 2472 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 2473 if (isascii) { 2474 PetscCall(PetscViewerASCIIPrintf(viewer,"BDDC interface preconditioner\n")); 2475 } 2476 PetscCall(PetscViewerASCIIPushTab(viewer)); 2477 PetscCall(PCView(bddcipc_ctx->bddc,viewer)); 2478 PetscCall(PetscViewerASCIIPopTab(viewer)); 2479 PetscFunctionReturn(0); 2480 } 2481 2482 static PetscErrorCode PCSetUp_BDDCIPC(PC pc) 2483 { 2484 BDDCIPC_ctx bddcipc_ctx; 2485 PetscBool isbddc; 2486 Vec vv; 2487 IS is; 2488 PC_IS *pcis; 2489 2490 PetscFunctionBegin; 2491 PetscCall(PCShellGetContext(pc,&bddcipc_ctx)); 2492 PetscCall(PetscObjectTypeCompare((PetscObject)bddcipc_ctx->bddc,PCBDDC,&isbddc)); 2493 PetscCheck(isbddc,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid type %s. Must be of type bddc",((PetscObject)bddcipc_ctx->bddc)->type_name); 2494 PetscCall(PCSetUp(bddcipc_ctx->bddc)); 2495 2496 /* create interface scatter */ 2497 pcis = (PC_IS*)(bddcipc_ctx->bddc->data); 2498 PetscCall(VecScatterDestroy(&bddcipc_ctx->g2l)); 2499 PetscCall(MatCreateVecs(pc->pmat,&vv,NULL)); 2500 PetscCall(ISRenumber(pcis->is_B_global,NULL,NULL,&is)); 2501 PetscCall(VecScatterCreate(vv,is,pcis->vec1_B,NULL,&bddcipc_ctx->g2l)); 2502 PetscCall(ISDestroy(&is)); 2503 PetscCall(VecDestroy(&vv)); 2504 PetscFunctionReturn(0); 2505 } 2506 2507 static PetscErrorCode PCApply_BDDCIPC(PC pc, Vec r, Vec x) 2508 { 2509 BDDCIPC_ctx bddcipc_ctx; 2510 PC_IS *pcis; 2511 VecScatter tmps; 2512 2513 PetscFunctionBegin; 2514 PetscCall(PCShellGetContext(pc,&bddcipc_ctx)); 2515 pcis = (PC_IS*)(bddcipc_ctx->bddc->data); 2516 tmps = pcis->global_to_B; 2517 pcis->global_to_B = bddcipc_ctx->g2l; 2518 PetscCall(PCBDDCScalingRestriction(bddcipc_ctx->bddc,r,pcis->vec1_B)); 2519 PetscCall(PCBDDCApplyInterfacePreconditioner(bddcipc_ctx->bddc,PETSC_FALSE)); 2520 PetscCall(PCBDDCScalingExtension(bddcipc_ctx->bddc,pcis->vec1_B,x)); 2521 pcis->global_to_B = tmps; 2522 PetscFunctionReturn(0); 2523 } 2524 2525 static PetscErrorCode PCApplyTranspose_BDDCIPC(PC pc, Vec r, Vec x) 2526 { 2527 BDDCIPC_ctx bddcipc_ctx; 2528 PC_IS *pcis; 2529 VecScatter tmps; 2530 2531 PetscFunctionBegin; 2532 PetscCall(PCShellGetContext(pc,&bddcipc_ctx)); 2533 pcis = (PC_IS*)(bddcipc_ctx->bddc->data); 2534 tmps = pcis->global_to_B; 2535 pcis->global_to_B = bddcipc_ctx->g2l; 2536 PetscCall(PCBDDCScalingRestriction(bddcipc_ctx->bddc,r,pcis->vec1_B)); 2537 PetscCall(PCBDDCApplyInterfacePreconditioner(bddcipc_ctx->bddc,PETSC_TRUE)); 2538 PetscCall(PCBDDCScalingExtension(bddcipc_ctx->bddc,pcis->vec1_B,x)); 2539 pcis->global_to_B = tmps; 2540 PetscFunctionReturn(0); 2541 } 2542 2543 static PetscErrorCode PCDestroy_BDDCIPC(PC pc) 2544 { 2545 BDDCIPC_ctx bddcipc_ctx; 2546 2547 PetscFunctionBegin; 2548 PetscCall(PCShellGetContext(pc,&bddcipc_ctx)); 2549 PetscCall(PCDestroy(&bddcipc_ctx->bddc)); 2550 PetscCall(VecScatterDestroy(&bddcipc_ctx->g2l)); 2551 PetscCall(PetscFree(bddcipc_ctx)); 2552 PetscFunctionReturn(0); 2553 } 2554 2555 /*@ 2556 PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system 2557 2558 Collective 2559 2560 Input Parameters: 2561 + fetidp_mat - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators 2562 - fetidp_flux_sol - the solution of the FETI-DP linear system 2563 2564 Output Parameters: 2565 . standard_sol - the solution defined on the physical domain 2566 2567 Level: developer 2568 2569 Notes: 2570 2571 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS 2572 @*/ 2573 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 2574 { 2575 FETIDPMat_ctx mat_ctx; 2576 2577 PetscFunctionBegin; 2578 PetscValidHeaderSpecific(fetidp_mat,MAT_CLASSID,1); 2579 PetscValidHeaderSpecific(fetidp_flux_sol,VEC_CLASSID,2); 2580 PetscValidHeaderSpecific(standard_sol,VEC_CLASSID,3); 2581 PetscCall(MatShellGetContext(fetidp_mat,&mat_ctx)); 2582 PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol)); 2583 PetscFunctionReturn(0); 2584 } 2585 2586 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, PetscBool fully_redundant, const char* prefix, Mat *fetidp_mat, PC *fetidp_pc) 2587 { 2588 2589 FETIDPMat_ctx fetidpmat_ctx; 2590 Mat newmat; 2591 FETIDPPC_ctx fetidppc_ctx; 2592 PC newpc; 2593 MPI_Comm comm; 2594 2595 PetscFunctionBegin; 2596 PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 2597 /* FETI-DP matrix */ 2598 PetscCall(PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx)); 2599 fetidpmat_ctx->fully_redundant = fully_redundant; 2600 PetscCall(PCBDDCSetupFETIDPMatContext(fetidpmat_ctx)); 2601 PetscCall(MatCreateShell(comm,fetidpmat_ctx->n,fetidpmat_ctx->n,fetidpmat_ctx->N,fetidpmat_ctx->N,fetidpmat_ctx,&newmat)); 2602 PetscCall(PetscObjectSetName((PetscObject)newmat,!fetidpmat_ctx->l2g_lambda_only ? "F" : "G")); 2603 PetscCall(MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult)); 2604 PetscCall(MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose)); 2605 PetscCall(MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat)); 2606 /* propagate MatOptions */ 2607 { 2608 PC_BDDC *pcbddc = (PC_BDDC*)fetidpmat_ctx->pc->data; 2609 PetscBool issym; 2610 2611 PetscCall(MatGetOption(pc->mat,MAT_SYMMETRIC,&issym)); 2612 if (issym || pcbddc->symmetric_primal) { 2613 PetscCall(MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE)); 2614 } 2615 } 2616 PetscCall(MatSetOptionsPrefix(newmat,prefix)); 2617 PetscCall(MatAppendOptionsPrefix(newmat,"fetidp_")); 2618 PetscCall(MatSetUp(newmat)); 2619 /* FETI-DP preconditioner */ 2620 PetscCall(PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx)); 2621 PetscCall(PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx)); 2622 PetscCall(PCCreate(comm,&newpc)); 2623 PetscCall(PCSetOperators(newpc,newmat,newmat)); 2624 PetscCall(PCSetOptionsPrefix(newpc,prefix)); 2625 PetscCall(PCAppendOptionsPrefix(newpc,"fetidp_")); 2626 PetscCall(PCSetErrorIfFailure(newpc,pc->erroriffailure)); 2627 if (!fetidpmat_ctx->l2g_lambda_only) { /* standard FETI-DP */ 2628 PetscCall(PCSetType(newpc,PCSHELL)); 2629 PetscCall(PCShellSetName(newpc,"FETI-DP multipliers")); 2630 PetscCall(PCShellSetContext(newpc,fetidppc_ctx)); 2631 PetscCall(PCShellSetApply(newpc,FETIDPPCApply)); 2632 PetscCall(PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose)); 2633 PetscCall(PCShellSetView(newpc,FETIDPPCView)); 2634 PetscCall(PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC)); 2635 } else { /* saddle-point FETI-DP */ 2636 Mat M; 2637 PetscInt psize; 2638 PetscBool fake = PETSC_FALSE, isfieldsplit; 2639 2640 PetscCall(ISViewFromOptions(fetidpmat_ctx->lagrange,NULL,"-lag_view")); 2641 PetscCall(ISViewFromOptions(fetidpmat_ctx->pressure,NULL,"-press_view")); 2642 PetscCall(PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_PPmat",(PetscObject*)&M)); 2643 PetscCall(PCSetType(newpc,PCFIELDSPLIT)); 2644 PetscCall(PCFieldSplitSetIS(newpc,"lag",fetidpmat_ctx->lagrange)); 2645 PetscCall(PCFieldSplitSetIS(newpc,"p",fetidpmat_ctx->pressure)); 2646 PetscCall(PCFieldSplitSetType(newpc,PC_COMPOSITE_SCHUR)); 2647 PetscCall(PCFieldSplitSetSchurFactType(newpc,PC_FIELDSPLIT_SCHUR_FACT_DIAG)); 2648 PetscCall(ISGetSize(fetidpmat_ctx->pressure,&psize)); 2649 if (psize != M->rmap->N) { 2650 Mat M2; 2651 PetscInt lpsize; 2652 2653 fake = PETSC_TRUE; 2654 PetscCall(ISGetLocalSize(fetidpmat_ctx->pressure,&lpsize)); 2655 PetscCall(MatCreate(comm,&M2)); 2656 PetscCall(MatSetType(M2,MATAIJ)); 2657 PetscCall(MatSetSizes(M2,lpsize,lpsize,psize,psize)); 2658 PetscCall(MatSetUp(M2)); 2659 PetscCall(MatAssemblyBegin(M2,MAT_FINAL_ASSEMBLY)); 2660 PetscCall(MatAssemblyEnd(M2,MAT_FINAL_ASSEMBLY)); 2661 PetscCall(PCFieldSplitSetSchurPre(newpc,PC_FIELDSPLIT_SCHUR_PRE_USER,M2)); 2662 PetscCall(MatDestroy(&M2)); 2663 } else { 2664 PetscCall(PCFieldSplitSetSchurPre(newpc,PC_FIELDSPLIT_SCHUR_PRE_USER,M)); 2665 } 2666 PetscCall(PCFieldSplitSetSchurScale(newpc,1.0)); 2667 2668 /* we need to setfromoptions and setup here to access the blocks */ 2669 PetscCall(PCSetFromOptions(newpc)); 2670 PetscCall(PCSetUp(newpc)); 2671 2672 /* user may have changed the type (e.g. -fetidp_pc_type none) */ 2673 PetscCall(PetscObjectTypeCompare((PetscObject)newpc,PCFIELDSPLIT,&isfieldsplit)); 2674 if (isfieldsplit) { 2675 KSP *ksps; 2676 PC ppc,lagpc; 2677 PetscInt nn; 2678 PetscBool ismatis,matisok = PETSC_FALSE,check = PETSC_FALSE; 2679 2680 /* set the solver for the (0,0) block */ 2681 PetscCall(PCFieldSplitSchurGetSubKSP(newpc,&nn,&ksps)); 2682 if (!nn) { /* not of type PC_COMPOSITE_SCHUR */ 2683 PetscCall(PCFieldSplitGetSubKSP(newpc,&nn,&ksps)); 2684 if (!fake) { /* pass pmat to the pressure solver */ 2685 Mat F; 2686 2687 PetscCall(KSPGetOperators(ksps[1],&F,NULL)); 2688 PetscCall(KSPSetOperators(ksps[1],F,M)); 2689 } 2690 } else { 2691 PetscBool issym; 2692 Mat S; 2693 2694 PetscCall(PCFieldSplitSchurGetS(newpc,&S)); 2695 2696 PetscCall(MatGetOption(newmat,MAT_SYMMETRIC,&issym)); 2697 if (issym) { 2698 PetscCall(MatSetOption(S,MAT_SYMMETRIC,PETSC_TRUE)); 2699 } 2700 } 2701 PetscCall(KSPGetPC(ksps[0],&lagpc)); 2702 PetscCall(PCSetType(lagpc,PCSHELL)); 2703 PetscCall(PCShellSetName(lagpc,"FETI-DP multipliers")); 2704 PetscCall(PCShellSetContext(lagpc,fetidppc_ctx)); 2705 PetscCall(PCShellSetApply(lagpc,FETIDPPCApply)); 2706 PetscCall(PCShellSetApplyTranspose(lagpc,FETIDPPCApplyTranspose)); 2707 PetscCall(PCShellSetView(lagpc,FETIDPPCView)); 2708 PetscCall(PCShellSetDestroy(lagpc,PCBDDCDestroyFETIDPPC)); 2709 2710 /* Olof's idea: interface Schur complement preconditioner for the mass matrix */ 2711 PetscCall(KSPGetPC(ksps[1],&ppc)); 2712 if (fake) { 2713 BDDCIPC_ctx bddcipc_ctx; 2714 PetscContainer c; 2715 2716 matisok = PETSC_TRUE; 2717 2718 /* create inner BDDC solver */ 2719 PetscCall(PetscNew(&bddcipc_ctx)); 2720 PetscCall(PCCreate(comm,&bddcipc_ctx->bddc)); 2721 PetscCall(PCSetType(bddcipc_ctx->bddc,PCBDDC)); 2722 PetscCall(PCSetOperators(bddcipc_ctx->bddc,M,M)); 2723 PetscCall(PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_pCSR",(PetscObject*)&c)); 2724 PetscCall(PetscObjectTypeCompare((PetscObject)M,MATIS,&ismatis)); 2725 if (c && ismatis) { 2726 Mat lM; 2727 PetscInt *csr,n; 2728 2729 PetscCall(MatISGetLocalMat(M,&lM)); 2730 PetscCall(MatGetSize(lM,&n,NULL)); 2731 PetscCall(PetscContainerGetPointer(c,(void**)&csr)); 2732 PetscCall(PCBDDCSetLocalAdjacencyGraph(bddcipc_ctx->bddc,n,csr,csr + (n + 1),PETSC_COPY_VALUES)); 2733 PetscCall(MatISRestoreLocalMat(M,&lM)); 2734 } 2735 PetscCall(PCSetOptionsPrefix(bddcipc_ctx->bddc,((PetscObject)ksps[1])->prefix)); 2736 PetscCall(PCSetErrorIfFailure(bddcipc_ctx->bddc,pc->erroriffailure)); 2737 PetscCall(PCSetFromOptions(bddcipc_ctx->bddc)); 2738 2739 /* wrap the interface application */ 2740 PetscCall(PCSetType(ppc,PCSHELL)); 2741 PetscCall(PCShellSetName(ppc,"FETI-DP pressure")); 2742 PetscCall(PCShellSetContext(ppc,bddcipc_ctx)); 2743 PetscCall(PCShellSetSetUp(ppc,PCSetUp_BDDCIPC)); 2744 PetscCall(PCShellSetApply(ppc,PCApply_BDDCIPC)); 2745 PetscCall(PCShellSetApplyTranspose(ppc,PCApplyTranspose_BDDCIPC)); 2746 PetscCall(PCShellSetView(ppc,PCView_BDDCIPC)); 2747 PetscCall(PCShellSetDestroy(ppc,PCDestroy_BDDCIPC)); 2748 } 2749 2750 /* determine if we need to assemble M to construct a preconditioner */ 2751 if (!matisok) { 2752 PetscCall(PetscObjectTypeCompare((PetscObject)M,MATIS,&ismatis)); 2753 PetscCall(PetscObjectTypeCompareAny((PetscObject)ppc,&matisok,PCBDDC,PCJACOBI,PCNONE,PCMG,"")); 2754 if (ismatis && !matisok) { 2755 PetscCall(MatConvert(M,MATAIJ,MAT_INPLACE_MATRIX,&M)); 2756 } 2757 } 2758 2759 /* run the subproblems to check convergence */ 2760 PetscCall(PetscOptionsGetBool(NULL,((PetscObject)newmat)->prefix,"-check_saddlepoint",&check,NULL)); 2761 if (check) { 2762 PetscInt i; 2763 2764 for (i=0;i<nn;i++) { 2765 KSP kspC; 2766 PC pc; 2767 Mat F,pF; 2768 Vec x,y; 2769 PetscBool isschur,prec = PETSC_TRUE; 2770 2771 PetscCall(KSPCreate(PetscObjectComm((PetscObject)ksps[i]),&kspC)); 2772 PetscCall(KSPSetOptionsPrefix(kspC,((PetscObject)ksps[i])->prefix)); 2773 PetscCall(KSPAppendOptionsPrefix(kspC,"check_")); 2774 PetscCall(KSPGetOperators(ksps[i],&F,&pF)); 2775 PetscCall(PetscObjectTypeCompare((PetscObject)F,MATSCHURCOMPLEMENT,&isschur)); 2776 if (isschur) { 2777 KSP kspS,kspS2; 2778 Mat A00,pA00,A10,A01,A11; 2779 char prefix[256]; 2780 2781 PetscCall(MatSchurComplementGetKSP(F,&kspS)); 2782 PetscCall(MatSchurComplementGetSubMatrices(F,&A00,&pA00,&A01,&A10,&A11)); 2783 PetscCall(MatCreateSchurComplement(A00,pA00,A01,A10,A11,&F)); 2784 PetscCall(MatSchurComplementGetKSP(F,&kspS2)); 2785 PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%sschur_",((PetscObject)kspC)->prefix)); 2786 PetscCall(KSPSetOptionsPrefix(kspS2,prefix)); 2787 PetscCall(KSPGetPC(kspS2,&pc)); 2788 PetscCall(PCSetType(pc,PCKSP)); 2789 PetscCall(PCKSPSetKSP(pc,kspS)); 2790 PetscCall(KSPSetFromOptions(kspS2)); 2791 PetscCall(KSPGetPC(kspS2,&pc)); 2792 PetscCall(PCSetUseAmat(pc,PETSC_TRUE)); 2793 } else { 2794 PetscCall(PetscObjectReference((PetscObject)F)); 2795 } 2796 PetscCall(KSPSetFromOptions(kspC)); 2797 PetscCall(PetscOptionsGetBool(NULL,((PetscObject)kspC)->prefix,"-preconditioned",&prec,NULL)); 2798 if (prec) { 2799 PetscCall(KSPGetPC(ksps[i],&pc)); 2800 PetscCall(KSPSetPC(kspC,pc)); 2801 } 2802 PetscCall(KSPSetOperators(kspC,F,pF)); 2803 PetscCall(MatCreateVecs(F,&x,&y)); 2804 PetscCall(VecSetRandom(x,NULL)); 2805 PetscCall(MatMult(F,x,y)); 2806 PetscCall(KSPSolve(kspC,y,x)); 2807 PetscCall(KSPCheckSolve(kspC,pc,x)); 2808 PetscCall(KSPDestroy(&kspC)); 2809 PetscCall(MatDestroy(&F)); 2810 PetscCall(VecDestroy(&x)); 2811 PetscCall(VecDestroy(&y)); 2812 } 2813 } 2814 PetscCall(PetscFree(ksps)); 2815 } 2816 } 2817 /* return pointers for objects created */ 2818 *fetidp_mat = newmat; 2819 *fetidp_pc = newpc; 2820 PetscFunctionReturn(0); 2821 } 2822 2823 /*@C 2824 PCBDDCCreateFETIDPOperators - Create FETI-DP operators 2825 2826 Collective 2827 2828 Input Parameters: 2829 + pc - the BDDC preconditioning context (setup should have been called before) 2830 . fully_redundant - true for a fully redundant set of Lagrange multipliers 2831 - prefix - optional options database prefix for the objects to be created (can be NULL) 2832 2833 Output Parameters: 2834 + fetidp_mat - shell FETI-DP matrix object 2835 - fetidp_pc - shell Dirichlet preconditioner for FETI-DP matrix 2836 2837 Level: developer 2838 2839 Notes: 2840 Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose 2841 2842 .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution 2843 @*/ 2844 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, PetscBool fully_redundant, const char *prefix, Mat *fetidp_mat, PC *fetidp_pc) 2845 { 2846 PetscFunctionBegin; 2847 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2848 if (pc->setupcalled) { 2849 PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,PetscBool,const char*,Mat*,PC*),(pc,fully_redundant,prefix,fetidp_mat,fetidp_pc)); 2850 } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"You must call PCSetup_BDDC() first"); 2851 PetscFunctionReturn(0); 2852 } 2853 /* -------------------------------------------------------------------------- */ 2854 /*MC 2855 PCBDDC - Balancing Domain Decomposition by Constraints. 2856 2857 An implementation of the BDDC preconditioner based on the bibliography found below. 2858 2859 The matrix to be preconditioned (Pmat) must be of type MATIS. 2860 2861 Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 2862 2863 It also works with unsymmetric and indefinite problems. 2864 2865 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. 2866 2867 Approximate local solvers are automatically adapted (see [1]) if the user has attached a nullspace object to the subdomain matrices, and informed BDDC of using approximate solvers (via the command line). 2868 2869 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() 2870 Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesIS() and their local counterparts. 2871 2872 Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD. 2873 2874 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. 2875 User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat() 2876 2877 The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object. 2878 2879 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. 2880 2881 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. 2882 2883 Options Database Keys (some of them, run with -help for a complete list): 2884 2885 + -pc_bddc_use_vertices <true> - use or not vertices in primal space 2886 . -pc_bddc_use_edges <true> - use or not edges in primal space 2887 . -pc_bddc_use_faces <false> - use or not faces in primal space 2888 . -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems 2889 . -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only) 2890 . -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested 2891 . -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1]) 2892 . -pc_bddc_levels <0> - maximum number of levels for multilevel 2893 . -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) 2894 . -pc_bddc_coarse_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) 2895 . -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling 2896 . -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) 2897 . -pc_bddc_adaptive_threshold <0.0> - when a value different than zero is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS or MKL_PARDISO installed) 2898 - -pc_bddc_check_level <0> - set verbosity level of debugging output 2899 2900 Options for Dirichlet, Neumann or coarse solver can be set with 2901 .vb 2902 -pc_bddc_dirichlet_ 2903 -pc_bddc_neumann_ 2904 -pc_bddc_coarse_ 2905 .ve 2906 e.g. -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KSPPREONLY and PCLU. 2907 2908 When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as 2909 .vb 2910 -pc_bddc_dirichlet_lN_ 2911 -pc_bddc_neumann_lN_ 2912 -pc_bddc_coarse_lN_ 2913 .ve 2914 Note that level number ranges from the finest (0) to the coarsest (N). 2915 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. 2916 .vb 2917 -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3 2918 .ve 2919 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 2920 2921 References: 2922 + * - C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149--168, March 2007 2923 . * - A. Klawonn and O. B. Widlund. "Dual-Primal FETI Methods for Linear Elasticity", Communications on Pure and Applied Mathematics Volume 59, Issue 11, pages 1523--1572, November 2006 2924 . * - J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", Computing Volume 83, Issue 2--3, pages 55--85, November 2008 2925 - * - 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 2926 2927 Level: intermediate 2928 2929 Developer Notes: 2930 2931 Contributed by Stefano Zampini 2932 2933 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 2934 M*/ 2935 2936 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 2937 { 2938 PC_BDDC *pcbddc; 2939 2940 PetscFunctionBegin; 2941 PetscCall(PetscNewLog(pc,&pcbddc)); 2942 pc->data = pcbddc; 2943 2944 /* create PCIS data structure */ 2945 PetscCall(PCISCreate(pc)); 2946 2947 /* create local graph structure */ 2948 PetscCall(PCBDDCGraphCreate(&pcbddc->mat_graph)); 2949 2950 /* BDDC nonzero defaults */ 2951 pcbddc->use_nnsp = PETSC_TRUE; 2952 pcbddc->use_local_adj = PETSC_TRUE; 2953 pcbddc->use_vertices = PETSC_TRUE; 2954 pcbddc->use_edges = PETSC_TRUE; 2955 pcbddc->symmetric_primal = PETSC_TRUE; 2956 pcbddc->vertex_size = 1; 2957 pcbddc->recompute_topography = PETSC_TRUE; 2958 pcbddc->coarse_size = -1; 2959 pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 2960 pcbddc->coarsening_ratio = 8; 2961 pcbddc->coarse_eqs_per_proc = 1; 2962 pcbddc->benign_compute_correction = PETSC_TRUE; 2963 pcbddc->nedfield = -1; 2964 pcbddc->nedglobal = PETSC_TRUE; 2965 pcbddc->graphmaxcount = PETSC_MAX_INT; 2966 pcbddc->sub_schurs_layers = -1; 2967 pcbddc->adaptive_threshold[0] = 0.0; 2968 pcbddc->adaptive_threshold[1] = 0.0; 2969 2970 /* function pointers */ 2971 pc->ops->apply = PCApply_BDDC; 2972 pc->ops->applytranspose = PCApplyTranspose_BDDC; 2973 pc->ops->setup = PCSetUp_BDDC; 2974 pc->ops->destroy = PCDestroy_BDDC; 2975 pc->ops->setfromoptions = PCSetFromOptions_BDDC; 2976 pc->ops->view = PCView_BDDC; 2977 pc->ops->applyrichardson = NULL; 2978 pc->ops->applysymmetricleft = NULL; 2979 pc->ops->applysymmetricright = NULL; 2980 pc->ops->presolve = PCPreSolve_BDDC; 2981 pc->ops->postsolve = PCPostSolve_BDDC; 2982 pc->ops->reset = PCReset_BDDC; 2983 2984 /* composing function */ 2985 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDiscreteGradient_C",PCBDDCSetDiscreteGradient_BDDC)); 2986 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDivergenceMat_C",PCBDDCSetDivergenceMat_BDDC)); 2987 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC)); 2988 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC)); 2989 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",PCBDDCSetPrimalVerticesIS_BDDC)); 2990 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetPrimalVerticesLocalIS_C",PCBDDCGetPrimalVerticesLocalIS_BDDC)); 2991 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetPrimalVerticesIS_C",PCBDDCGetPrimalVerticesIS_BDDC)); 2992 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC)); 2993 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC)); 2994 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC)); 2995 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC)); 2996 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC)); 2997 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC)); 2998 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC)); 2999 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC)); 3000 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC)); 3001 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC)); 3002 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC)); 3003 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC)); 3004 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC)); 3005 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC)); 3006 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC)); 3007 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC)); 3008 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC)); 3009 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC)); 3010 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_BDDC)); 3011 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_BDDC)); 3012 PetscFunctionReturn(0); 3013 } 3014 3015 /*@C 3016 PCBDDCInitializePackage - This function initializes everything in the PCBDDC package. It is called 3017 from PCInitializePackage(). 3018 3019 Level: developer 3020 3021 .seealso: PetscInitialize() 3022 @*/ 3023 PetscErrorCode PCBDDCInitializePackage(void) 3024 { 3025 int i; 3026 3027 PetscFunctionBegin; 3028 if (PCBDDCPackageInitialized) PetscFunctionReturn(0); 3029 PCBDDCPackageInitialized = PETSC_TRUE; 3030 PetscCall(PetscRegisterFinalize(PCBDDCFinalizePackage)); 3031 3032 /* general events */ 3033 PetscCall(PetscLogEventRegister("PCBDDCTopo",PC_CLASSID,&PC_BDDC_Topology[0])); 3034 PetscCall(PetscLogEventRegister("PCBDDCLKSP",PC_CLASSID,&PC_BDDC_LocalSolvers[0])); 3035 PetscCall(PetscLogEventRegister("PCBDDCLWor",PC_CLASSID,&PC_BDDC_LocalWork[0])); 3036 PetscCall(PetscLogEventRegister("PCBDDCCorr",PC_CLASSID,&PC_BDDC_CorrectionSetUp[0])); 3037 PetscCall(PetscLogEventRegister("PCBDDCASet",PC_CLASSID,&PC_BDDC_ApproxSetUp[0])); 3038 PetscCall(PetscLogEventRegister("PCBDDCAApp",PC_CLASSID,&PC_BDDC_ApproxApply[0])); 3039 PetscCall(PetscLogEventRegister("PCBDDCCSet",PC_CLASSID,&PC_BDDC_CoarseSetUp[0])); 3040 PetscCall(PetscLogEventRegister("PCBDDCCKSP",PC_CLASSID,&PC_BDDC_CoarseSolver[0])); 3041 PetscCall(PetscLogEventRegister("PCBDDCAdap",PC_CLASSID,&PC_BDDC_AdaptiveSetUp[0])); 3042 PetscCall(PetscLogEventRegister("PCBDDCScal",PC_CLASSID,&PC_BDDC_Scaling[0])); 3043 PetscCall(PetscLogEventRegister("PCBDDCSchr",PC_CLASSID,&PC_BDDC_Schurs[0])); 3044 PetscCall(PetscLogEventRegister("PCBDDCDirS",PC_CLASSID,&PC_BDDC_Solves[0][0])); 3045 PetscCall(PetscLogEventRegister("PCBDDCNeuS",PC_CLASSID,&PC_BDDC_Solves[0][1])); 3046 PetscCall(PetscLogEventRegister("PCBDDCCoaS",PC_CLASSID,&PC_BDDC_Solves[0][2])); 3047 for (i=1;i<PETSC_PCBDDC_MAXLEVELS;i++) { 3048 char ename[32]; 3049 3050 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCTopo l%02d",i)); 3051 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_Topology[i])); 3052 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCLKSP l%02d",i)); 3053 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_LocalSolvers[i])); 3054 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCLWor l%02d",i)); 3055 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_LocalWork[i])); 3056 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCCorr l%02d",i)); 3057 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_CorrectionSetUp[i])); 3058 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCASet l%02d",i)); 3059 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_ApproxSetUp[i])); 3060 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCAApp l%02d",i)); 3061 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_ApproxApply[i])); 3062 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCCSet l%02d",i)); 3063 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_CoarseSetUp[i])); 3064 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCCKSP l%02d",i)); 3065 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_CoarseSolver[i])); 3066 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCAdap l%02d",i)); 3067 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_AdaptiveSetUp[i])); 3068 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCScal l%02d",i)); 3069 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_Scaling[i])); 3070 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCSchr l%02d",i)); 3071 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_Schurs[i])); 3072 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCDirS l%02d",i)); 3073 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_Solves[i][0])); 3074 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCNeuS l%02d",i)); 3075 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_Solves[i][1])); 3076 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCBDDCCoaS l%02d",i)); 3077 PetscCall(PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_Solves[i][2])); 3078 } 3079 PetscFunctionReturn(0); 3080 } 3081 3082 /*@C 3083 PCBDDCFinalizePackage - This function frees everything from the PCBDDC package. It is 3084 called from PetscFinalize() automatically. 3085 3086 Level: developer 3087 3088 .seealso: PetscFinalize() 3089 @*/ 3090 PetscErrorCode PCBDDCFinalizePackage(void) 3091 { 3092 PetscFunctionBegin; 3093 PCBDDCPackageInitialized = PETSC_FALSE; 3094 PetscFunctionReturn(0); 3095 } 3096