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