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