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