1 /* 2 GAMG geometric-algebric multigrid PC - Mark Adams 2011 3 */ 4 #include <petsc/private/matimpl.h> 5 #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 6 #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h> /*I "petscksp.h" I*/ 7 8 #if defined(PETSC_HAVE_CUDA) 9 #include <cuda_runtime.h> 10 #endif 11 12 #if defined(PETSC_HAVE_HIP) 13 #include <hip/hip_runtime.h> 14 #endif 15 16 PetscLogEvent petsc_gamg_setup_events[GAMG_NUM_SET]; 17 PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3]; 18 19 // #define GAMG_STAGES 20 #if defined(GAMG_STAGES) 21 static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS]; 22 #endif 23 24 static PetscFunctionList GAMGList = NULL; 25 static PetscBool PCGAMGPackageInitialized; 26 27 /* ----------------------------------------------------------------------------- */ 28 PetscErrorCode PCReset_GAMG(PC pc) 29 { 30 PC_MG *mg = (PC_MG*)pc->data; 31 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 32 33 PetscFunctionBegin; 34 PetscCall(PetscFree(pc_gamg->data)); 35 pc_gamg->data_sz = 0; 36 PetscCall(PetscFree(pc_gamg->orig_data)); 37 for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS ; level++) { 38 mg->min_eigen_DinvA[level] = 0; 39 mg->max_eigen_DinvA[level] = 0; 40 } 41 pc_gamg->emin = 0; 42 pc_gamg->emax = 0; 43 PetscFunctionReturn(0); 44 } 45 46 /* -------------------------------------------------------------------------- */ 47 /* 48 PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 49 of active processors. 50 51 Input Parameter: 52 . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 53 'pc_gamg->data_sz' are changed via repartitioning/reduction. 54 . Amat_fine - matrix on this fine (k) level 55 . cr_bs - coarse block size 56 In/Output Parameter: 57 . a_P_inout - prolongation operator to the next level (k-->k-1) 58 . a_nactive_proc - number of active procs 59 Output Parameter: 60 . a_Amat_crs - coarse matrix that is created (k-1) 61 */ 62 static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,IS * Pcolumnperm, PetscBool is_last) 63 { 64 PC_MG *mg = (PC_MG*)pc->data; 65 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 66 Mat Cmat,Pold=*a_P_inout; 67 MPI_Comm comm; 68 PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc; 69 PetscInt ncrs_eq,ncrs,f_bs; 70 71 PetscFunctionBegin; 72 PetscCall(PetscObjectGetComm((PetscObject)Amat_fine,&comm)); 73 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 74 PetscCallMPI(MPI_Comm_size(comm, &size)); 75 PetscCall(MatGetBlockSize(Amat_fine, &f_bs)); 76 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP],0,0,0,0)); 77 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0)); 78 PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat)); 79 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0)); 80 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP],0,0,0,0)); 81 82 if (Pcolumnperm) *Pcolumnperm = NULL; 83 84 /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 85 PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL)); 86 if (pc_gamg->data_cell_rows>0) { 87 ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 88 } else { 89 PetscInt bs; 90 PetscCall(MatGetBlockSize(Cmat, &bs)); 91 ncrs = ncrs_eq/bs; 92 } 93 /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 94 if (pc_gamg->level_reduction_factors[pc_gamg->current_level] == 0 && PetscDefined(HAVE_CUDA) && pc_gamg->current_level==0) { /* 0 turns reducing to 1 process/device on; do for HIP, etc. */ 95 #if defined(PETSC_HAVE_CUDA) 96 PetscShmComm pshmcomm; 97 PetscMPIInt locrank; 98 MPI_Comm loccomm; 99 PetscInt s_nnodes,r_nnodes, new_new_size; 100 cudaError_t cerr; 101 int devCount; 102 PetscCall(PetscShmCommGet(comm,&pshmcomm)); 103 PetscCall(PetscShmCommGetMpiShmComm(pshmcomm,&loccomm)); 104 PetscCallMPI(MPI_Comm_rank(loccomm, &locrank)); 105 s_nnodes = !locrank; 106 PetscCallMPI(MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm)); 107 PetscCheck((size%r_nnodes) == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of nodes np=%d nnodes%" PetscInt_FMT,size,r_nnodes); 108 devCount = 0; 109 cerr = cudaGetDeviceCount(&devCount); 110 cudaGetLastError(); /* Reset the last error */ 111 if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */ 112 new_new_size = r_nnodes * devCount; 113 new_size = new_new_size; 114 PetscCall(PetscInfo(pc,"%s: Fine grid with Cuda. %" PetscInt_FMT " nodes. Change new active set size %d --> %d (devCount=%d #nodes=%" PetscInt_FMT ")\n",((PetscObject)pc)->prefix,r_nnodes,nactive,new_size,devCount,r_nnodes)); 115 } else { 116 PetscCall(PetscInfo(pc,"%s: With Cuda but no device. Use heuristics.",((PetscObject)pc)->prefix)); 117 goto HEURISTIC; 118 } 119 #else 120 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here"); 121 #endif 122 } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) { 123 PetscCheck(nactive%pc_gamg->level_reduction_factors[pc_gamg->current_level] == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of active process %d wrt reduction factor %" PetscInt_FMT,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]); 124 new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level]; 125 PetscCall(PetscInfo(pc,"%s: Manually setting reduction to %d active processes (%d/%" PetscInt_FMT ")\n",((PetscObject)pc)->prefix,new_size,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level])); 126 } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) { 127 new_size = 1; 128 PetscCall(PetscInfo(pc,"%s: Force coarsest grid reduction to %d active processes\n",((PetscObject)pc)->prefix,new_size)); 129 } else { 130 PetscInt ncrs_eq_glob; 131 #if defined(PETSC_HAVE_CUDA) 132 HEURISTIC: 133 #endif 134 PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL)); 135 new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 136 if (!new_size) new_size = 1; /* not likely, posible? */ 137 else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 138 PetscCall(PetscInfo(pc,"%s: Coarse grid reduction from %d to %d active processes\n",((PetscObject)pc)->prefix,nactive,new_size)); 139 } 140 if (new_size==nactive) { 141 *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 142 if (new_size < size) { 143 /* odd case where multiple coarse grids are on one processor or no coarsening ... */ 144 PetscCall(PetscInfo(pc,"%s: reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",((PetscObject)pc)->prefix,nactive)); 145 if (pc_gamg->cpu_pin_coarse_grids) { 146 PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE)); 147 PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE)); 148 } 149 } 150 /* we know that the grid structure can be reused in MatPtAP */ 151 } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */ 152 PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1; 153 IS is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices; 154 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REDUCE],0,0,0,0)); 155 nloc_old = ncrs_eq/cr_bs; 156 PetscCheck(ncrs_eq % cr_bs == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %" PetscInt_FMT " not divisible by cr_bs %" PetscInt_FMT,ncrs_eq,cr_bs); 157 /* get new_size and rfactor */ 158 if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) { 159 /* find factor */ 160 if (new_size == 1) rfactor = size; /* don't modify */ 161 else { 162 PetscReal best_fact = 0.; 163 jj = -1; 164 for (kk = 1 ; kk <= size ; kk++) { 165 if (!(size%kk)) { /* a candidate */ 166 PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size; 167 if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 168 if (fact > best_fact) { 169 best_fact = fact; jj = kk; 170 } 171 } 172 } 173 if (jj != -1) rfactor = jj; 174 else rfactor = 1; /* a prime */ 175 if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1; 176 else expand_factor = rfactor; 177 } 178 new_size = size/rfactor; /* make new size one that is factor */ 179 if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */ 180 *a_Amat_crs = Cmat; 181 PetscCall(PetscInfo(pc,"%s: Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%" PetscInt_FMT "\n",((PetscObject)pc)->prefix,new_size,ncrs_eq)); 182 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE],0,0,0,0)); 183 PetscFunctionReturn(0); 184 } 185 } 186 /* make 'is_eq_newproc' */ 187 PetscCall(PetscMalloc1(size, &counts)); 188 if (pc_gamg->repart) { /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */ 189 Mat adj; 190 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REPART],0,0,0,0)); 191 PetscCall(PetscInfo(pc,"%s: Repartition: size (active): %d --> %d, %" PetscInt_FMT " local equations, using %s process layout\n",((PetscObject)pc)->prefix,*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread")); 192 /* get 'adj' */ 193 if (cr_bs == 1) { 194 PetscCall(MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj)); 195 } else { 196 /* make a scalar matrix to partition (no Stokes here) */ 197 Mat tMat; 198 PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 199 const PetscScalar *vals; 200 const PetscInt *idx; 201 PetscInt *d_nnz, *o_nnz, M, N; 202 static PetscInt llev = 0; /* ugly but just used for debugging */ 203 MatType mtype; 204 205 PetscCall(PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz)); 206 PetscCall(MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs)); 207 PetscCall(MatGetSize(Cmat, &M, &N)); 208 for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 209 PetscCall(MatGetRow(Cmat,Ii,&ncols,NULL,NULL)); 210 d_nnz[jj] = ncols/cr_bs; 211 o_nnz[jj] = ncols/cr_bs; 212 PetscCall(MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL)); 213 if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 214 if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs; 215 } 216 217 PetscCall(MatGetType(Amat_fine,&mtype)); 218 PetscCall(MatCreate(comm, &tMat)); 219 PetscCall(MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE)); 220 PetscCall(MatSetType(tMat,mtype)); 221 PetscCall(MatSeqAIJSetPreallocation(tMat,0,d_nnz)); 222 PetscCall(MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz)); 223 PetscCall(PetscFree2(d_nnz,o_nnz)); 224 225 for (ii = Istart_crs; ii < Iend_crs; ii++) { 226 PetscInt dest_row = ii/cr_bs; 227 PetscCall(MatGetRow(Cmat,ii,&ncols,&idx,&vals)); 228 for (jj = 0; jj < ncols; jj++) { 229 PetscInt dest_col = idx[jj]/cr_bs; 230 PetscScalar v = 1.0; 231 PetscCall(MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES)); 232 } 233 PetscCall(MatRestoreRow(Cmat,ii,&ncols,&idx,&vals)); 234 } 235 PetscCall(MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY)); 236 PetscCall(MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY)); 237 238 if (llev++ == -1) { 239 PetscViewer viewer; char fname[32]; 240 PetscCall(PetscSNPrintf(fname,sizeof(fname),"part_mat_%" PetscInt_FMT ".mat",llev)); 241 PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer); 242 PetscCall(MatView(tMat, viewer)); 243 PetscCall(PetscViewerDestroy(&viewer)); 244 } 245 PetscCall(MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj)); 246 PetscCall(MatDestroy(&tMat)); 247 } /* create 'adj' */ 248 249 { /* partition: get newproc_idx */ 250 char prefix[256]; 251 const char *pcpre; 252 const PetscInt *is_idx; 253 MatPartitioning mpart; 254 IS proc_is; 255 256 PetscCall(MatPartitioningCreate(comm, &mpart)); 257 PetscCall(MatPartitioningSetAdjacency(mpart, adj)); 258 PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 259 PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "")); 260 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix)); 261 PetscCall(MatPartitioningSetFromOptions(mpart)); 262 PetscCall(MatPartitioningSetNParts(mpart, new_size)); 263 PetscCall(MatPartitioningApply(mpart, &proc_is)); 264 PetscCall(MatPartitioningDestroy(&mpart)); 265 266 /* collect IS info */ 267 PetscCall(PetscMalloc1(ncrs_eq, &newproc_idx)); 268 PetscCall(ISGetIndices(proc_is, &is_idx)); 269 for (kk = jj = 0 ; kk < nloc_old ; kk++) { 270 for (ii = 0 ; ii < cr_bs ; ii++, jj++) { 271 newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ 272 } 273 } 274 PetscCall(ISRestoreIndices(proc_is, &is_idx)); 275 PetscCall(ISDestroy(&proc_is)); 276 } 277 PetscCall(MatDestroy(&adj)); 278 279 PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc)); 280 PetscCall(PetscFree(newproc_idx)); 281 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART],0,0,0,0)); 282 } else { /* simple aggregation of parts -- 'is_eq_newproc' */ 283 PetscInt targetPE; 284 PetscCheck(new_size!=nactive,PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen"); 285 PetscCall(PetscInfo(pc,"%s: Number of equations (loc) %" PetscInt_FMT " with simple aggregation\n",((PetscObject)pc)->prefix,ncrs_eq)); 286 targetPE = (rank/rfactor)*expand_factor; 287 PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc)); 288 } /* end simple 'is_eq_newproc' */ 289 290 /* 291 Create an index set from the is_eq_newproc index set to indicate the mapping TO 292 */ 293 PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num)); 294 is_eq_num_prim = is_eq_num; 295 /* 296 Determine how many equations/vertices are assigned to each processor 297 */ 298 PetscCall(ISPartitioningCount(is_eq_newproc, size, counts)); 299 ncrs_eq_new = counts[rank]; 300 PetscCall(ISDestroy(&is_eq_newproc)); 301 ncrs_new = ncrs_eq_new/cr_bs; 302 303 PetscCall(PetscFree(counts)); 304 /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxiliary data into some complex abstracted thing */ 305 { 306 Vec src_crd, dest_crd; 307 const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols; 308 VecScatter vecscat; 309 PetscScalar *array; 310 IS isscat; 311 /* move data (for primal equations only) */ 312 /* Create a vector to contain the newly ordered element information */ 313 PetscCall(VecCreate(comm, &dest_crd)); 314 PetscCall(VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE)); 315 PetscCall(VecSetType(dest_crd,VECSTANDARD)); /* this is needed! */ 316 /* 317 There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 318 a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 319 */ 320 PetscCall(PetscMalloc1(ncrs*node_data_sz, &tidx)); 321 PetscCall(ISGetIndices(is_eq_num_prim, &idx)); 322 for (ii=0,jj=0; ii<ncrs; ii++) { 323 PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 324 for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 325 } 326 PetscCall(ISRestoreIndices(is_eq_num_prim, &idx)); 327 PetscCall(ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat)); 328 PetscCall(PetscFree(tidx)); 329 /* 330 Create a vector to contain the original vertex information for each element 331 */ 332 PetscCall(VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd)); 333 for (jj=0; jj<ndata_cols; jj++) { 334 const PetscInt stride0=ncrs*pc_gamg->data_cell_rows; 335 for (ii=0; ii<ncrs; ii++) { 336 for (kk=0; kk<ndata_rows; kk++) { 337 PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 338 PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 339 PetscCall(VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES)); 340 } 341 } 342 } 343 PetscCall(VecAssemblyBegin(src_crd)); 344 PetscCall(VecAssemblyEnd(src_crd)); 345 /* 346 Scatter the element vertex information (still in the original vertex ordering) 347 to the correct processor 348 */ 349 PetscCall(VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat)); 350 PetscCall(ISDestroy(&isscat)); 351 PetscCall(VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD)); 352 PetscCall(VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD)); 353 PetscCall(VecScatterDestroy(&vecscat)); 354 PetscCall(VecDestroy(&src_crd)); 355 /* 356 Put the element vertex data into a new allocation of the gdata->ele 357 */ 358 PetscCall(PetscFree(pc_gamg->data)); 359 PetscCall(PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data)); 360 361 pc_gamg->data_sz = node_data_sz*ncrs_new; 362 strideNew = ncrs_new*ndata_rows; 363 364 PetscCall(VecGetArray(dest_crd, &array)); 365 for (jj=0; jj<ndata_cols; jj++) { 366 for (ii=0; ii<ncrs_new; ii++) { 367 for (kk=0; kk<ndata_rows; kk++) { 368 PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 369 pc_gamg->data[ix] = PetscRealPart(array[jx]); 370 } 371 } 372 } 373 PetscCall(VecRestoreArray(dest_crd, &array)); 374 PetscCall(VecDestroy(&dest_crd)); 375 } 376 /* move A and P (columns) with new layout */ 377 /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0)); */ 378 /* 379 Invert for MatCreateSubMatrix 380 */ 381 PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices)); 382 PetscCall(ISSort(new_eq_indices)); /* is this needed? */ 383 PetscCall(ISSetBlockSize(new_eq_indices, cr_bs)); 384 if (is_eq_num != is_eq_num_prim) { 385 PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */ 386 } 387 if (Pcolumnperm) { 388 PetscCall(PetscObjectReference((PetscObject)new_eq_indices)); 389 *Pcolumnperm = new_eq_indices; 390 } 391 PetscCall(ISDestroy(&is_eq_num)); 392 /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0)); */ 393 /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 394 395 /* 'a_Amat_crs' output */ 396 { 397 Mat mat; 398 PetscBool isset,isspd,isher; 399 #if !defined(PETSC_USE_COMPLEX) 400 PetscBool issym; 401 #endif 402 403 PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat)); 404 PetscCall(MatIsSPDKnown(Cmat, &isset, &isspd)); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ? 405 if (isset) PetscCall(MatSetOption(mat, MAT_SPD,isspd)); 406 else { 407 PetscCall(MatIsHermitianKnown(Cmat,&isset,&isher)); 408 if (isset) PetscCall(MatSetOption(mat, MAT_HERMITIAN,isher)); 409 else { 410 #if !defined(PETSC_USE_COMPLEX) 411 PetscCall(MatIsSymmetricKnown(Cmat,&isset, &issym)); 412 if (isset) PetscCall(MatSetOption(mat, MAT_SYMMETRIC,issym)); 413 #endif 414 } 415 } 416 *a_Amat_crs = mat; 417 } 418 PetscCall(MatDestroy(&Cmat)); 419 420 /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 421 /* prolongator */ 422 { 423 IS findices; 424 PetscInt Istart,Iend; 425 Mat Pnew; 426 427 PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend)); 428 /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 429 PetscCall(ISCreateStride(comm,Iend-Istart,Istart,1,&findices)); 430 PetscCall(ISSetBlockSize(findices,f_bs)); 431 PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew)); 432 PetscCall(ISDestroy(&findices)); 433 PetscCall(MatSetOption(Pnew,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE)); 434 435 /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 436 PetscCall(MatDestroy(a_P_inout)); 437 438 /* output - repartitioned */ 439 *a_P_inout = Pnew; 440 } 441 PetscCall(ISDestroy(&new_eq_indices)); 442 443 *a_nactive_proc = new_size; /* output */ 444 445 /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */ 446 if (pc_gamg->cpu_pin_coarse_grids) { 447 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 448 static PetscInt llev = 2; 449 PetscCall(PetscInfo(pc,"%s: Pinning level %" PetscInt_FMT " to the CPU\n",((PetscObject)pc)->prefix,llev++)); 450 #endif 451 PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE)); 452 PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE)); 453 if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */ 454 Mat A = *a_Amat_crs, P = *a_P_inout; 455 PetscMPIInt size; 456 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 457 if (size > 1) { 458 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data; 459 PetscCall(VecBindToCPU(a->lvec,PETSC_TRUE)); 460 PetscCall(VecBindToCPU(p->lvec,PETSC_TRUE)); 461 } 462 } 463 } 464 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE],0,0,0,0)); 465 } // processor reduce 466 PetscFunctionReturn(0); 467 } 468 469 // used in GEO 470 PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2) 471 { 472 const char *prefix; 473 char addp[32]; 474 PC_MG *mg = (PC_MG*)a_pc->data; 475 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 476 477 PetscFunctionBegin; 478 PetscCall(PCGetOptionsPrefix(a_pc,&prefix)); 479 PetscCall(PetscInfo(a_pc,"%s: Square Graph on level %" PetscInt_FMT "\n",((PetscObject)a_pc)->prefix,pc_gamg->current_level+1)); 480 PetscCall(MatProductCreate(Gmat1,Gmat1,NULL,Gmat2)); 481 PetscCall(MatSetOptionsPrefix(*Gmat2,prefix)); 482 PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%" PetscInt_FMT "_",pc_gamg->current_level)); 483 PetscCall(MatAppendOptionsPrefix(*Gmat2,addp)); 484 if ((*Gmat2)->structurally_symmetric == PETSC_BOOL3_TRUE) { 485 PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AB)); 486 } else { 487 PetscCall(MatSetOption(Gmat1,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE)); 488 PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AtB)); 489 } 490 PetscCall(MatProductSetFromOptions(*Gmat2)); 491 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0)); 492 PetscCall(MatProductSymbolic(*Gmat2)); 493 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0)); 494 PetscCall(MatProductClear(*Gmat2)); 495 /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */ 496 (*Gmat2)->assembled = PETSC_TRUE; 497 PetscFunctionReturn(0); 498 } 499 500 /* -------------------------------------------------------------------------- */ 501 /* 502 PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 503 by setting data structures and options. 504 505 Input Parameter: 506 . pc - the preconditioner context 507 508 */ 509 PetscErrorCode PCSetUp_GAMG(PC pc) 510 { 511 PC_MG *mg = (PC_MG*)pc->data; 512 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 513 Mat Pmat = pc->pmat; 514 PetscInt fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS]; 515 MPI_Comm comm; 516 PetscMPIInt rank,size,nactivepe; 517 Mat Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS]; 518 IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 519 PetscLogDouble nnz0=0.,nnztot=0.; 520 MatInfo info; 521 PetscBool is_last = PETSC_FALSE; 522 523 PetscFunctionBegin; 524 PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 525 PetscCallMPI(MPI_Comm_rank(comm,&rank)); 526 PetscCallMPI(MPI_Comm_size(comm,&size)); 527 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0)); 528 if (pc->setupcalled) { 529 if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) { 530 /* reset everything */ 531 PetscCall(PCReset_MG(pc)); 532 pc->setupcalled = 0; 533 } else { 534 PC_MG_Levels **mglevels = mg->levels; 535 /* just do Galerkin grids */ 536 Mat B,dA,dB; 537 if (pc_gamg->Nlevels > 1) { 538 PetscInt gl; 539 /* currently only handle case where mat and pmat are the same on coarser levels */ 540 PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB)); 541 /* (re)set to get dirty flag */ 542 PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB)); 543 544 for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) { 545 MatReuse reuse = MAT_INITIAL_MATRIX ; 546 #if defined(GAMG_STAGES) 547 PetscCall(PetscLogStagePush(gamg_stages[gl])); 548 #endif 549 /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 550 PetscCall(KSPGetOperators(mglevels[level]->smoothd,NULL,&B)); 551 if (B->product) { 552 if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) { 553 reuse = MAT_REUSE_MATRIX; 554 } 555 } 556 if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A)); 557 if (reuse == MAT_REUSE_MATRIX) { 558 PetscCall(PetscInfo(pc,"%s: RAP after first solve, reuse matrix level %" PetscInt_FMT "\n",((PetscObject)pc)->prefix,level)); 559 } else { 560 PetscCall(PetscInfo(pc,"%s: RAP after first solve, new matrix level %" PetscInt_FMT "\n",((PetscObject)pc)->prefix,level)); 561 } 562 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0)); 563 PetscCall(MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B)); 564 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0)); 565 if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B; 566 PetscCall(KSPSetOperators(mglevels[level]->smoothd,B,B)); 567 dB = B; 568 #if defined(GAMG_STAGES) 569 PetscCall(PetscLogStagePop()); 570 #endif 571 } 572 } 573 574 PetscCall(PCSetUp_MG(pc)); 575 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0)); 576 PetscFunctionReturn(0); 577 } 578 } 579 580 if (!pc_gamg->data) { 581 if (pc_gamg->orig_data) { 582 PetscCall(MatGetBlockSize(Pmat, &bs)); 583 PetscCall(MatGetLocalSize(Pmat, &qq, NULL)); 584 585 pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 586 pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 587 pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 588 589 PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data)); 590 for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 591 } else { 592 PetscCheck(pc_gamg->ops->createdefaultdata,comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 593 PetscCall(pc_gamg->ops->createdefaultdata(pc,Pmat)); 594 } 595 } 596 597 /* cache original data for reuse */ 598 if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 599 PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data)); 600 for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 601 pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 602 pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 603 } 604 605 /* get basic dims */ 606 PetscCall(MatGetBlockSize(Pmat, &bs)); 607 PetscCall(MatGetSize(Pmat, &M, &N)); 608 609 PetscCall(MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info)); /* global reduction */ 610 nnz0 = info.nz_used; 611 nnztot = info.nz_used; 612 PetscCall(PetscInfo(pc,"%s: level %d) N=%" PetscInt_FMT ", n data rows=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", np=%d\n",((PetscObject)pc)->prefix,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(PetscInt)(nnz0/(PetscReal)M+0.5),size)); 613 614 /* Get A_i and R_i */ 615 for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) { 616 pc_gamg->current_level = level; 617 PetscCheck(level < PETSC_MG_MAXLEVELS,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %" PetscInt_FMT,level); 618 level1 = level + 1; 619 #if defined(GAMG_STAGES) 620 if (!gamg_stages[level]) { 621 char str[32]; 622 sprintf(str,"GAMG Level %d",(int)level); 623 PetscCall(PetscLogStageRegister(str, &gamg_stages[level])); 624 } 625 PetscCall(PetscLogStagePush(gamg_stages[level])); 626 #endif 627 { /* construct prolongator */ 628 Mat Gmat; 629 PetscCoarsenData *agg_lists; 630 Mat Prol11; 631 632 PetscCall(pc_gamg->ops->graph(pc,Aarr[level], &Gmat)); 633 PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists)); 634 PetscCall(pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11)); 635 636 /* could have failed to create new level */ 637 if (Prol11) { 638 const char *prefix; 639 char addp[32]; 640 641 /* get new block size of coarse matrices */ 642 PetscCall(MatGetBlockSizes(Prol11, NULL, &bs)); 643 644 if (pc_gamg->ops->optprolongator) { 645 /* smooth */ 646 PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11)); 647 } 648 649 if (pc_gamg->use_aggs_in_asm) { 650 PetscInt bs; 651 PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // not timed directly, ugly, could remove, but good ASM method 652 PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level])); 653 } 654 655 PetscCall(PCGetOptionsPrefix(pc,&prefix)); 656 PetscCall(MatSetOptionsPrefix(Prol11,prefix)); 657 PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level)); 658 PetscCall(MatAppendOptionsPrefix(Prol11,addp)); 659 /* Always generate the transpose with CUDA 660 Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */ 661 PetscCall(MatSetOption(Prol11,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE)); 662 PetscCall(MatSetFromOptions(Prol11)); 663 Parr[level1] = Prol11; 664 } else Parr[level1] = NULL; /* failed to coarsen */ 665 666 PetscCall(MatDestroy(&Gmat)); 667 PetscCall(PetscCDDestroy(agg_lists)); 668 } /* construct prolongator scope */ 669 if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 670 if (!Parr[level1]) { /* failed to coarsen */ 671 PetscCall(PetscInfo(pc,"%s: Stop gridding, level %" PetscInt_FMT "\n",((PetscObject)pc)->prefix,level)); 672 #if defined(GAMG_STAGES) 673 PetscCall(PetscLogStagePop()); 674 #endif 675 break; 676 } 677 PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */ 678 PetscCheck(!is_last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?"); 679 if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 680 if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE; 681 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL],0,0,0,0)); 682 PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last)); 683 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL],0,0,0,0)); 684 685 PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */ 686 PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info)); 687 nnztot += info.nz_used; 688 PetscCall(PetscInfo(pc,"%s: %d) N=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", %d active pes\n",((PetscObject)pc)->prefix,(int)level1,M,pc_gamg->data_cell_cols,(PetscInt)(info.nz_used/(PetscReal)M),nactivepe)); 689 690 #if defined(GAMG_STAGES) 691 PetscCall(PetscLogStagePop()); 692 #endif 693 /* stop if one node or one proc -- could pull back for singular problems */ 694 if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) { 695 PetscCall(PetscInfo(pc,"%s: HARD stop of coarsening on level %" PetscInt_FMT ". Grid too small: %" PetscInt_FMT " block nodes\n",((PetscObject)pc)->prefix,level,M/bs)); 696 level++; 697 break; 698 } 699 } /* levels */ 700 PetscCall(PetscFree(pc_gamg->data)); 701 702 PetscCall(PetscInfo(pc,"%s: %" PetscInt_FMT " levels, grid complexity = %g\n",((PetscObject)pc)->prefix,level+1,nnztot/nnz0)); 703 pc_gamg->Nlevels = level + 1; 704 fine_level = level; 705 PetscCall(PCMGSetLevels(pc,pc_gamg->Nlevels,NULL)); 706 707 if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 708 709 /* set default smoothers & set operators */ 710 for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) { 711 KSP smoother; 712 PC subpc; 713 714 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 715 PetscCall(KSPGetPC(smoother, &subpc)); 716 717 PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 718 /* set ops */ 719 PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level])); 720 PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level+1])); 721 722 /* set defaults */ 723 PetscCall(KSPSetType(smoother, KSPCHEBYSHEV)); 724 725 /* set blocks for ASM smoother that uses the 'aggregates' */ 726 if (pc_gamg->use_aggs_in_asm) { 727 PetscInt sz; 728 IS *iss; 729 730 sz = nASMBlocksArr[level]; 731 iss = ASMLocalIDsArr[level]; 732 PetscCall(PCSetType(subpc, PCASM)); 733 PetscCall(PCASMSetOverlap(subpc, 0)); 734 PetscCall(PCASMSetType(subpc,PC_ASM_BASIC)); 735 if (!sz) { 736 IS is; 737 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is)); 738 PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is)); 739 PetscCall(ISDestroy(&is)); 740 } else { 741 PetscInt kk; 742 PetscCall(PCASMSetLocalSubdomains(subpc, sz, NULL, iss)); 743 for (kk=0; kk<sz; kk++) { 744 PetscCall(ISDestroy(&iss[kk])); 745 } 746 PetscCall(PetscFree(iss)); 747 } 748 ASMLocalIDsArr[level] = NULL; 749 nASMBlocksArr[level] = 0; 750 } else { 751 PetscCall(PCSetType(subpc, PCJACOBI)); 752 } 753 } 754 { 755 /* coarse grid */ 756 KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 757 Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 758 759 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 760 PetscCall(KSPSetOperators(smoother, Lmat, Lmat)); 761 if (!pc_gamg->use_parallel_coarse_grid_solver) { 762 PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 763 PetscCall(KSPGetPC(smoother, &subpc)); 764 PetscCall(PCSetType(subpc, PCBJACOBI)); 765 PetscCall(PCSetUp(subpc)); 766 PetscCall(PCBJacobiGetSubKSP(subpc,&ii,&first,&k2)); 767 PetscCheck(ii == 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %" PetscInt_FMT " is not one",ii); 768 PetscCall(KSPGetPC(k2[0],&pc2)); 769 PetscCall(PCSetType(pc2, PCLU)); 770 PetscCall(PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS)); 771 PetscCall(KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1)); 772 PetscCall(KSPSetType(k2[0], KSPPREONLY)); 773 } 774 } 775 776 /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 777 PetscObjectOptionsBegin((PetscObject)pc); 778 PetscCall(PCSetFromOptions_MG(PetscOptionsObject,pc)); 779 PetscOptionsEnd(); 780 PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL)); 781 782 /* setup cheby eigen estimates from SA */ 783 if (pc_gamg->use_sa_esteig) { 784 for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) { 785 KSP smoother; 786 PetscBool ischeb; 787 788 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 789 PetscCall(PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb)); 790 if (ischeb) { 791 KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data; 792 793 // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG 794 if (mg->max_eigen_DinvA[level] > 0) { 795 // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it. 796 // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.) 797 PetscReal emax,emin; 798 799 emin = mg->min_eigen_DinvA[level]; 800 emax = mg->max_eigen_DinvA[level]; 801 PetscCall(PetscInfo(pc,"%s: PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %" PetscInt_FMT " (N=%" PetscInt_FMT ") with emax = %g emin = %g\n",((PetscObject)pc)->prefix,level,Aarr[level]->rmap->N,(double)emax,(double)emin)); 802 cheb->emin_provided = emin; 803 cheb->emax_provided = emax; 804 } 805 } 806 } 807 } 808 809 PetscCall(PCSetUp_MG(pc)); 810 811 /* clean up */ 812 for (level=1; level<pc_gamg->Nlevels; level++) { 813 PetscCall(MatDestroy(&Parr[level])); 814 PetscCall(MatDestroy(&Aarr[level])); 815 } 816 } else { 817 KSP smoother; 818 819 PetscCall(PetscInfo(pc,"%s: One level solver used (system is seen as DD). Using default solver.\n",((PetscObject)pc)->prefix)); 820 PetscCall(PCMGGetSmoother(pc, 0, &smoother)); 821 PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0])); 822 PetscCall(KSPSetType(smoother, KSPPREONLY)); 823 PetscCall(PCSetUp_MG(pc)); 824 } 825 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0)); 826 PetscFunctionReturn(0); 827 } 828 829 /* ------------------------------------------------------------------------- */ 830 /* 831 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 832 that was created with PCCreate_GAMG(). 833 834 Input Parameter: 835 . pc - the preconditioner context 836 837 Application Interface Routine: PCDestroy() 838 */ 839 PetscErrorCode PCDestroy_GAMG(PC pc) 840 { 841 PC_MG *mg = (PC_MG*)pc->data; 842 PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 843 844 PetscFunctionBegin; 845 PetscCall(PCReset_GAMG(pc)); 846 if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc)); 847 PetscCall(PetscFree(pc_gamg->ops)); 848 PetscCall(PetscFree(pc_gamg->gamg_type_name)); 849 PetscCall(PetscFree(pc_gamg)); 850 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",NULL)); 851 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",NULL)); 852 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",NULL)); 853 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",NULL)); 854 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",NULL)); 855 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",NULL)); 856 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",NULL)); 857 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",NULL)); 858 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",NULL)); 859 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",NULL)); 860 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",NULL)); 861 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",NULL)); 862 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",NULL)); 863 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",NULL)); 864 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",NULL)); 865 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",NULL)); 866 PetscCall(PCDestroy_MG(pc)); 867 PetscFunctionReturn(0); 868 } 869 870 /*@ 871 PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction. 872 873 Logically Collective on PC 874 875 Input Parameters: 876 + pc - the preconditioner context 877 - n - the number of equations 878 879 Options Database Key: 880 . -pc_gamg_process_eq_limit <limit> - set the limit 881 882 Notes: 883 GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process 884 that has degrees of freedom 885 886 Level: intermediate 887 888 .seealso: `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()` 889 @*/ 890 PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 891 { 892 PetscFunctionBegin; 893 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 894 PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n)); 895 PetscFunctionReturn(0); 896 } 897 898 static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 899 { 900 PC_MG *mg = (PC_MG*)pc->data; 901 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 902 903 PetscFunctionBegin; 904 if (n>0) pc_gamg->min_eq_proc = n; 905 PetscFunctionReturn(0); 906 } 907 908 /*@ 909 PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid. 910 911 Collective on PC 912 913 Input Parameters: 914 + pc - the preconditioner context 915 - n - maximum number of equations to aim for 916 917 Options Database Key: 918 . -pc_gamg_coarse_eq_limit <limit> - set the limit 919 920 Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 921 has less than 1000 unknowns. 922 923 Level: intermediate 924 925 .seealso: `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()` 926 @*/ 927 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 928 { 929 PetscFunctionBegin; 930 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 931 PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n)); 932 PetscFunctionReturn(0); 933 } 934 935 static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 936 { 937 PC_MG *mg = (PC_MG*)pc->data; 938 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 939 940 PetscFunctionBegin; 941 if (n>0) pc_gamg->coarse_eq_limit = n; 942 PetscFunctionReturn(0); 943 } 944 945 /*@ 946 PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids 947 948 Collective on PC 949 950 Input Parameters: 951 + pc - the preconditioner context 952 - n - PETSC_TRUE or PETSC_FALSE 953 954 Options Database Key: 955 . -pc_gamg_repartition <true,false> - turn on the repartitioning 956 957 Notes: 958 this will generally improve the loading balancing of the work on each level 959 960 Level: intermediate 961 962 @*/ 963 PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) 964 { 965 PetscFunctionBegin; 966 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 967 PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n)); 968 PetscFunctionReturn(0); 969 } 970 971 static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) 972 { 973 PC_MG *mg = (PC_MG*)pc->data; 974 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 975 976 PetscFunctionBegin; 977 pc_gamg->repart = n; 978 PetscFunctionReturn(0); 979 } 980 981 /*@ 982 PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Chebyshev smoother 983 984 Collective on PC 985 986 Input Parameters: 987 + pc - the preconditioner context 988 - n - number of its 989 990 Options Database Key: 991 . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate 992 993 Notes: 994 Smoothed aggregation constructs the smoothed prolongator $P = (I - \omega D^{-1} A) T$ where $T$ is the tentative prolongator and $D$ is the diagonal of $A$. 995 Eigenvalue estimates (based on a few CG or GMRES iterations) are computed to choose $\omega$ so that this is a stable smoothing operation. 996 If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused. 997 This option is only used when the smoother uses Jacobi, and should be turned off if a different PCJacobiType is used. 998 It became default in PETSc 3.17. 999 1000 Level: advanced 1001 1002 .seealso: `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()` 1003 @*/ 1004 PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n) 1005 { 1006 PetscFunctionBegin; 1007 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1008 PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n)); 1009 PetscFunctionReturn(0); 1010 } 1011 1012 static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n) 1013 { 1014 PC_MG *mg = (PC_MG*)pc->data; 1015 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1016 1017 PetscFunctionBegin; 1018 pc_gamg->use_sa_esteig = n; 1019 PetscFunctionReturn(0); 1020 } 1021 1022 /*@ 1023 PCGAMGSetEigenvalues - Set eigenvalues 1024 1025 Collective on PC 1026 1027 Input Parameters: 1028 + pc - the preconditioner context 1029 - emax - max eigenvalue 1030 . emin - min eigenvalue 1031 1032 Options Database Key: 1033 . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues 1034 1035 Level: intermediate 1036 1037 .seealso: `PCGAMGSetUseSAEstEig()` 1038 @*/ 1039 PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin) 1040 { 1041 PetscFunctionBegin; 1042 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1043 PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin)); 1044 PetscFunctionReturn(0); 1045 } 1046 1047 static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin) 1048 { 1049 PC_MG *mg = (PC_MG*)pc->data; 1050 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1051 1052 PetscFunctionBegin; 1053 PetscCheck(emax > emin,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin); 1054 PetscCheck(emax*emin > 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin); 1055 pc_gamg->emax = emax; 1056 pc_gamg->emin = emin; 1057 PetscFunctionReturn(0); 1058 } 1059 1060 /*@ 1061 PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner 1062 1063 Collective on PC 1064 1065 Input Parameters: 1066 + pc - the preconditioner context 1067 - n - PETSC_TRUE or PETSC_FALSE 1068 1069 Options Database Key: 1070 . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation 1071 1072 Level: intermediate 1073 1074 Notes: 1075 May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1076 rebuilding the preconditioner quicker. 1077 1078 @*/ 1079 PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 1080 { 1081 PetscFunctionBegin; 1082 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1083 PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n)); 1084 PetscFunctionReturn(0); 1085 } 1086 1087 static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 1088 { 1089 PC_MG *mg = (PC_MG*)pc->data; 1090 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1091 1092 PetscFunctionBegin; 1093 pc_gamg->reuse_prol = n; 1094 PetscFunctionReturn(0); 1095 } 1096 1097 /*@ 1098 PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner. 1099 1100 Collective on PC 1101 1102 Input Parameters: 1103 + pc - the preconditioner context 1104 - flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not 1105 1106 Options Database Key: 1107 . -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains 1108 1109 Level: intermediate 1110 1111 @*/ 1112 PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) 1113 { 1114 PetscFunctionBegin; 1115 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1116 PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg)); 1117 PetscFunctionReturn(0); 1118 } 1119 1120 static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) 1121 { 1122 PC_MG *mg = (PC_MG*)pc->data; 1123 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1124 1125 PetscFunctionBegin; 1126 pc_gamg->use_aggs_in_asm = flg; 1127 PetscFunctionReturn(0); 1128 } 1129 1130 /*@ 1131 PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver 1132 1133 Collective on PC 1134 1135 Input Parameters: 1136 + pc - the preconditioner context 1137 - flg - PETSC_TRUE to not force coarse grid onto one processor 1138 1139 Options Database Key: 1140 . -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver 1141 1142 Level: intermediate 1143 1144 .seealso: `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()` 1145 @*/ 1146 PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg) 1147 { 1148 PetscFunctionBegin; 1149 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1150 PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg)); 1151 PetscFunctionReturn(0); 1152 } 1153 1154 static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) 1155 { 1156 PC_MG *mg = (PC_MG*)pc->data; 1157 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1158 1159 PetscFunctionBegin; 1160 pc_gamg->use_parallel_coarse_grid_solver = flg; 1161 PetscFunctionReturn(0); 1162 } 1163 1164 /*@ 1165 PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU 1166 1167 Collective on PC 1168 1169 Input Parameters: 1170 + pc - the preconditioner context 1171 - flg - PETSC_TRUE to pin coarse grids to CPU 1172 1173 Options Database Key: 1174 . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU 1175 1176 Level: intermediate 1177 1178 .seealso: `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetUseParallelCoarseGridSolve()` 1179 @*/ 1180 PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) 1181 { 1182 PetscFunctionBegin; 1183 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1184 PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg)); 1185 PetscFunctionReturn(0); 1186 } 1187 1188 static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) 1189 { 1190 PC_MG *mg = (PC_MG*)pc->data; 1191 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1192 1193 PetscFunctionBegin; 1194 pc_gamg->cpu_pin_coarse_grids = flg; 1195 PetscFunctionReturn(0); 1196 } 1197 1198 /*@ 1199 PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type) 1200 1201 Collective on PC 1202 1203 Input Parameters: 1204 + pc - the preconditioner context 1205 - flg - Layout type 1206 1207 Options Database Key: 1208 . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering 1209 1210 Level: intermediate 1211 1212 .seealso: `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()` 1213 @*/ 1214 PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) 1215 { 1216 PetscFunctionBegin; 1217 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1218 PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg)); 1219 PetscFunctionReturn(0); 1220 } 1221 1222 static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) 1223 { 1224 PC_MG *mg = (PC_MG*)pc->data; 1225 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1226 1227 PetscFunctionBegin; 1228 pc_gamg->layout_type = flg; 1229 PetscFunctionReturn(0); 1230 } 1231 1232 /*@ 1233 PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use 1234 1235 Not collective on PC 1236 1237 Input Parameters: 1238 + pc - the preconditioner 1239 - n - the maximum number of levels to use 1240 1241 Options Database Key: 1242 . -pc_mg_levels <n> - set the maximum number of levels to allow 1243 1244 Level: intermediate 1245 1246 @*/ 1247 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 1248 { 1249 PetscFunctionBegin; 1250 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1251 PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n)); 1252 PetscFunctionReturn(0); 1253 } 1254 1255 static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 1256 { 1257 PC_MG *mg = (PC_MG*)pc->data; 1258 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1259 1260 PetscFunctionBegin; 1261 pc_gamg->Nlevels = n; 1262 PetscFunctionReturn(0); 1263 } 1264 1265 /*@ 1266 PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 1267 1268 Not collective on PC 1269 1270 Input Parameters: 1271 + pc - the preconditioner context 1272 . v - array of threshold values for finest n levels; 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph 1273 - n - number of threshold values provided in array 1274 1275 Options Database Key: 1276 . -pc_gamg_threshold <threshold> - the threshold to drop edges 1277 1278 Notes: 1279 Increasing the threshold decreases the rate of coarsening. Conversely reducing the threshold increases the rate of coarsening (aggressive coarsening) and thereby reduces the complexity of the coarse grids, and generally results in slower solver converge rates. Reducing coarse grid complexity reduced the complexity of Galerkin coarse grid construction considerably. 1280 Before coarsening or aggregating the graph, GAMG removes small values from the graph with this threshold, and thus reducing the coupling in the graph and a different (perhaps better) coarser set of points. 1281 1282 If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening. 1283 In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold(). 1284 If n is greater than the total number of levels, the excess entries in threshold will not be used. 1285 1286 Level: intermediate 1287 1288 .seealso: `PCGAMGFilterGraph()`, `PCGAMGSetAggressiveLevels()` 1289 @*/ 1290 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) 1291 { 1292 PetscFunctionBegin; 1293 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1294 if (n) PetscValidRealPointer(v,2); 1295 PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n)); 1296 PetscFunctionReturn(0); 1297 } 1298 1299 static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) 1300 { 1301 PC_MG *mg = (PC_MG*)pc->data; 1302 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1303 PetscInt i; 1304 PetscFunctionBegin; 1305 for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i]; 1306 for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale; 1307 PetscFunctionReturn(0); 1308 } 1309 1310 /*@ 1311 PCGAMGSetRankReductionFactors - Set manual schedule for process reduction on coarse grids 1312 1313 Collective on PC 1314 1315 Input Parameters: 1316 + pc - the preconditioner context 1317 . v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda 1318 - n - number of values provided in array 1319 1320 Options Database Key: 1321 . -pc_gamg_rank_reduction_factors <factors> - provide the schedule 1322 1323 Level: intermediate 1324 1325 .seealso: `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()` 1326 @*/ 1327 PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n) 1328 { 1329 PetscFunctionBegin; 1330 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1331 if (n) PetscValidIntPointer(v,2); 1332 PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n)); 1333 PetscFunctionReturn(0); 1334 } 1335 1336 static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n) 1337 { 1338 PC_MG *mg = (PC_MG*)pc->data; 1339 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1340 PetscInt i; 1341 PetscFunctionBegin; 1342 for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i]; 1343 for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */ 1344 PetscFunctionReturn(0); 1345 } 1346 1347 /*@ 1348 PCGAMGSetThresholdScale - Relative threshold reduction at each level 1349 1350 Not collective on PC 1351 1352 Input Parameters: 1353 + pc - the preconditioner context 1354 - scale - the threshold value reduction, usually < 1.0 1355 1356 Options Database Key: 1357 . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level 1358 1359 Notes: 1360 The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold(). 1361 This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold(). 1362 1363 Level: advanced 1364 1365 .seealso: `PCGAMGSetThreshold()` 1366 @*/ 1367 PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) 1368 { 1369 PetscFunctionBegin; 1370 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1371 PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v)); 1372 PetscFunctionReturn(0); 1373 } 1374 1375 static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) 1376 { 1377 PC_MG *mg = (PC_MG*)pc->data; 1378 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1379 PetscFunctionBegin; 1380 pc_gamg->threshold_scale = v; 1381 PetscFunctionReturn(0); 1382 } 1383 1384 /*@C 1385 PCGAMGSetType - Set solution method 1386 1387 Collective on PC 1388 1389 Input Parameters: 1390 + pc - the preconditioner context 1391 - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1392 1393 Options Database Key: 1394 . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply 1395 1396 Level: intermediate 1397 1398 .seealso: `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType` 1399 @*/ 1400 PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1401 { 1402 PetscFunctionBegin; 1403 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1404 PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type)); 1405 PetscFunctionReturn(0); 1406 } 1407 1408 /*@C 1409 PCGAMGGetType - Get solution method 1410 1411 Collective on PC 1412 1413 Input Parameter: 1414 . pc - the preconditioner context 1415 1416 Output Parameter: 1417 . type - the type of algorithm used 1418 1419 Level: intermediate 1420 1421 .seealso: `PCGAMGSetType()`, `PCGAMGType` 1422 @*/ 1423 PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1424 { 1425 PetscFunctionBegin; 1426 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1427 PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type)); 1428 PetscFunctionReturn(0); 1429 } 1430 1431 static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1432 { 1433 PC_MG *mg = (PC_MG*)pc->data; 1434 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1435 1436 PetscFunctionBegin; 1437 *type = pc_gamg->type; 1438 PetscFunctionReturn(0); 1439 } 1440 1441 static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1442 { 1443 PC_MG *mg = (PC_MG*)pc->data; 1444 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1445 PetscErrorCode (*r)(PC); 1446 1447 PetscFunctionBegin; 1448 pc_gamg->type = type; 1449 PetscCall(PetscFunctionListFind(GAMGList,type,&r)); 1450 PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 1451 if (pc_gamg->ops->destroy) { 1452 PetscCall((*pc_gamg->ops->destroy)(pc)); 1453 PetscCall(PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps))); 1454 pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1455 /* cleaning up common data in pc_gamg - this should disapear someday */ 1456 pc_gamg->data_cell_cols = 0; 1457 pc_gamg->data_cell_rows = 0; 1458 pc_gamg->orig_data_cell_cols = 0; 1459 pc_gamg->orig_data_cell_rows = 0; 1460 PetscCall(PetscFree(pc_gamg->data)); 1461 pc_gamg->data_sz = 0; 1462 } 1463 PetscCall(PetscFree(pc_gamg->gamg_type_name)); 1464 PetscCall(PetscStrallocpy(type,&pc_gamg->gamg_type_name)); 1465 PetscCall((*r)(pc)); 1466 PetscFunctionReturn(0); 1467 } 1468 1469 static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer) 1470 { 1471 PC_MG *mg = (PC_MG*)pc->data; 1472 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1473 PetscReal gc=0, oc=0; 1474 1475 PetscFunctionBegin; 1476 PetscCall(PetscViewerASCIIPrintf(viewer," GAMG specific options\n")); 1477 PetscCall(PetscViewerASCIIPrintf(viewer," Threshold for dropping small values in graph on each level =")); 1478 for (PetscInt i=0;i<mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i])); 1479 PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); 1480 PetscCall(PetscViewerASCIIPrintf(viewer," Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale)); 1481 if (pc_gamg->use_aggs_in_asm) { 1482 PetscCall(PetscViewerASCIIPrintf(viewer," Using aggregates from coarsening process to define subdomains for PCASM\n")); 1483 } 1484 if (pc_gamg->use_parallel_coarse_grid_solver) { 1485 PetscCall(PetscViewerASCIIPrintf(viewer," Using parallel coarse grid solver (all coarse grid equations not put on one process)\n")); 1486 } 1487 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 1488 if (pc_gamg->cpu_pin_coarse_grids) { 1489 /* PetscCall(PetscViewerASCIIPrintf(viewer," Pinning coarse grids to the CPU)\n")); */ 1490 } 1491 #endif 1492 /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */ 1493 /* PetscCall(PetscViewerASCIIPrintf(viewer," Put reduced grids on processes in natural order (ie, 0,1,2...)\n")); */ 1494 /* } else { */ 1495 /* PetscCall(PetscViewerASCIIPrintf(viewer," Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n")); */ 1496 /* } */ 1497 if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc,viewer)); 1498 PetscCall(PCMGGetGridComplexity(pc,&gc,&oc)); 1499 PetscCall(PetscViewerASCIIPrintf(viewer," Complexity: grid = %g operator = %g\n",(double)gc,(double)oc)); 1500 PetscFunctionReturn(0); 1501 } 1502 1503 PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc) 1504 { 1505 PC_MG *mg = (PC_MG*)pc->data; 1506 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1507 PetscBool flag; 1508 MPI_Comm comm; 1509 char prefix[256],tname[32]; 1510 PetscInt i,n; 1511 const char *pcpre; 1512 static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL}; 1513 PetscFunctionBegin; 1514 PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 1515 PetscOptionsHeadBegin(PetscOptionsObject,"GAMG options"); 1516 PetscCall(PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag)); 1517 if (flag) PetscCall(PCGAMGSetType(pc,tname)); 1518 PetscCall(PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL)); 1519 PetscCall(PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",pc_gamg->use_sa_esteig,&pc_gamg->use_sa_esteig,NULL)); 1520 PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL)); 1521 PetscCall(PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation aggregates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL)); 1522 PetscCall(PetscOptionsBool("-pc_gamg_use_parallel_coarse_grid_solver","Use parallel coarse grid solver (otherwise put last grid on one process)","PCGAMGSetUseParallelCoarseGridSolve",pc_gamg->use_parallel_coarse_grid_solver,&pc_gamg->use_parallel_coarse_grid_solver,NULL)); 1523 PetscCall(PetscOptionsBool("-pc_gamg_cpu_pin_coarse_grids","Pin coarse grids to the CPU","PCGAMGSetCpuPinCoarseGrids",pc_gamg->cpu_pin_coarse_grids,&pc_gamg->cpu_pin_coarse_grids,NULL)); 1524 PetscCall(PetscOptionsEnum("-pc_gamg_coarse_grid_layout_type","compact: place reduced grids on processes in natural order; spread: distribute to whole machine for more memory bandwidth","PCGAMGSetCoarseGridLayoutType",LayoutTypes,(PetscEnum)pc_gamg->layout_type,(PetscEnum*)&pc_gamg->layout_type,NULL)); 1525 PetscCall(PetscOptionsInt("-pc_gamg_process_eq_limit","Limit (goal) on number of equations per process on coarse grids","PCGAMGSetProcEqLim",pc_gamg->min_eq_proc,&pc_gamg->min_eq_proc,NULL)); 1526 PetscCall(PetscOptionsInt("-pc_gamg_coarse_eq_limit","Limit on number of equations for the coarse grid","PCGAMGSetCoarseEqLim",pc_gamg->coarse_eq_limit,&pc_gamg->coarse_eq_limit,NULL)); 1527 PetscCall(PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL)); 1528 n = PETSC_MG_MAXLEVELS; 1529 PetscCall(PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag)); 1530 if (!flag || n < PETSC_MG_MAXLEVELS) { 1531 if (!flag) n = 1; 1532 i = n; 1533 do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS); 1534 } 1535 n = PETSC_MG_MAXLEVELS; 1536 PetscCall(PetscOptionsIntArray("-pc_gamg_rank_reduction_factors","Manual schedule of coarse grid reduction factors that overrides internal heuristics (0 for first reduction puts one process/device)","PCGAMGSetRankReductionFactors",pc_gamg->level_reduction_factors,&n,&flag)); 1537 if (!flag) i = 0; 1538 else i = n; 1539 do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS); 1540 PetscCall(PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL)); 1541 { 1542 PetscReal eminmax[2] = {0., 0.}; 1543 n = 2; 1544 PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag)); 1545 if (flag) { 1546 PetscCheck(n == 2,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 1547 PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0])); 1548 } 1549 } 1550 /* set options for subtype */ 1551 if (pc_gamg->ops->setfromoptions) PetscCall((*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc)); 1552 1553 PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 1554 PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "")); 1555 PetscOptionsHeadEnd(); 1556 PetscFunctionReturn(0); 1557 } 1558 1559 /* -------------------------------------------------------------------------- */ 1560 /*MC 1561 PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 1562 1563 Options Database Keys: 1564 + -pc_gamg_type <type> - one of agg, geo, or classical 1565 . -pc_gamg_repartition <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined 1566 - -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations 1567 + -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother 1568 . -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> 1569 equations on each process that has degrees of freedom 1570 - -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 1571 + -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level 1572 . -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1573 1574 Options Database Keys for default Aggregation: 1575 + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 1576 . -pc_gamg_symmetrize_graph <true,default=false> - symmetrize the graph before computing the aggregation 1577 - -pc_gamg_square_graph <n,default=1> - alias for pc_gamg_aggressive_coarsening (deprecated) 1578 - -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest. 1579 1580 Multigrid options: 1581 + -pc_mg_cycles <v> - v or w, see PCMGSetCycleType() 1582 . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1583 . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1584 - -pc_mg_levels <levels> - Number of levels of multigrid to use. 1585 1586 Notes: 1587 In order to obtain good performance for PCGAMG for vector valued problems you must 1588 Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point 1589 Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator 1590 See the Users Manual Chapter 4 for more details 1591 1592 Level: intermediate 1593 1594 .seealso: `PCCreate()`, `PCSetType()`, `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, 1595 `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGSetUseSAEstEig()` 1596 M*/ 1597 1598 PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 1599 { 1600 PC_GAMG *pc_gamg; 1601 PC_MG *mg; 1602 1603 PetscFunctionBegin; 1604 /* register AMG type */ 1605 PetscCall(PCGAMGInitializePackage()); 1606 1607 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 1608 PetscCall(PCSetType(pc, PCMG)); 1609 PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG)); 1610 1611 /* create a supporting struct and attach it to pc */ 1612 PetscCall(PetscNewLog(pc,&pc_gamg)); 1613 PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL)); 1614 mg = (PC_MG*)pc->data; 1615 mg->innerctx = pc_gamg; 1616 1617 PetscCall(PetscNewLog(pc,&pc_gamg->ops)); 1618 1619 /* these should be in subctx but repartitioning needs simple arrays */ 1620 pc_gamg->data_sz = 0; 1621 pc_gamg->data = NULL; 1622 1623 /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 1624 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 1625 pc->ops->setup = PCSetUp_GAMG; 1626 pc->ops->reset = PCReset_GAMG; 1627 pc->ops->destroy = PCDestroy_GAMG; 1628 mg->view = PCView_GAMG; 1629 1630 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG)); 1631 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG)); 1632 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG)); 1633 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG)); 1634 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG)); 1635 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG)); 1636 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG)); 1637 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG)); 1638 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG)); 1639 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG)); 1640 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG)); 1641 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG)); 1642 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG)); 1643 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG)); 1644 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG)); 1645 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG)); 1646 pc_gamg->repart = PETSC_FALSE; 1647 pc_gamg->reuse_prol = PETSC_FALSE; 1648 pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1649 pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1650 pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1651 pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1652 pc_gamg->min_eq_proc = 50; 1653 pc_gamg->coarse_eq_limit = 50; 1654 PetscCall(PetscArrayzero(pc_gamg->threshold,PETSC_MG_MAXLEVELS)); 1655 pc_gamg->threshold_scale = 1.; 1656 pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 1657 pc_gamg->current_level = 0; /* don't need to init really */ 1658 pc_gamg->use_sa_esteig = PETSC_TRUE; 1659 pc_gamg->emin = 0; 1660 pc_gamg->emax = 0; 1661 1662 pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1663 1664 /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1665 PetscCall(PCGAMGSetType(pc,PCGAMGAGG)); 1666 PetscFunctionReturn(0); 1667 } 1668 1669 /*@C 1670 PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 1671 from PCInitializePackage(). 1672 1673 Level: developer 1674 1675 .seealso: `PetscInitialize()` 1676 @*/ 1677 PetscErrorCode PCGAMGInitializePackage(void) 1678 { 1679 PetscInt l; 1680 1681 PetscFunctionBegin; 1682 if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 1683 PCGAMGPackageInitialized = PETSC_TRUE; 1684 PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO)); 1685 PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG)); 1686 PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical)); 1687 PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage)); 1688 1689 /* general events */ 1690 PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP])); 1691 PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH])); 1692 PetscCall(PetscLogEventRegister(" PCGAMGFilter", PC_CLASSID, &petsc_gamg_setup_events[GAMG_FILTER])); 1693 PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN])); 1694 PetscCall(PetscLogEventRegister(" GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS])); 1695 PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL])); 1696 PetscCall(PetscLogEventRegister(" GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA])); 1697 PetscCall(PetscLogEventRegister(" GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB])); 1698 PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT])); 1699 PetscCall(PetscLogEventRegister(" GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM])); 1700 PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL])); 1701 PetscCall(PetscLogEventRegister(" GAMG PtAP",PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP])); 1702 PetscCall(PetscLogEventRegister(" GAMG Reduce",PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE])); 1703 PetscCall(PetscLogEventRegister(" GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART])); 1704 /* PetscCall(PetscLogEventRegister(" GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */ 1705 /* PetscCall(PetscLogEventRegister(" GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */ 1706 /* PetscCall(PetscLogEventRegister(" GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */ 1707 for (l=0;l<PETSC_MG_MAXLEVELS;l++) { 1708 char ename[32]; 1709 1710 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02" PetscInt_FMT,l)); 1711 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0])); 1712 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02" PetscInt_FMT,l)); 1713 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1])); 1714 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02" PetscInt_FMT,l)); 1715 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2])); 1716 } 1717 #if defined(GAMG_STAGES) 1718 { /* create timer stages */ 1719 char str[32]; 1720 sprintf(str,"GAMG Level %d",0); 1721 PetscCall(PetscLogStageRegister(str, &gamg_stages[0])); 1722 } 1723 #endif 1724 PetscFunctionReturn(0); 1725 } 1726 1727 /*@C 1728 PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is 1729 called from PetscFinalize() automatically. 1730 1731 Level: developer 1732 1733 .seealso: `PetscFinalize()` 1734 @*/ 1735 PetscErrorCode PCGAMGFinalizePackage(void) 1736 { 1737 PetscFunctionBegin; 1738 PCGAMGPackageInitialized = PETSC_FALSE; 1739 PetscCall(PetscFunctionListDestroy(&GAMGList)); 1740 PetscFunctionReturn(0); 1741 } 1742 1743 /*@C 1744 PCGAMGRegister - Register a PCGAMG implementation. 1745 1746 Input Parameters: 1747 + type - string that will be used as the name of the GAMG type. 1748 - create - function for creating the gamg context. 1749 1750 Level: advanced 1751 1752 .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()` 1753 @*/ 1754 PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1755 { 1756 PetscFunctionBegin; 1757 PetscCall(PCGAMGInitializePackage()); 1758 PetscCall(PetscFunctionListAdd(&GAMGList,type,create)); 1759 PetscFunctionReturn(0); 1760 } 1761