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 flg; 399 PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat)); 400 PetscCall(MatGetOption(Cmat, MAT_SPD, &flg)); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ? 401 if (flg) { 402 PetscCall(MatSetOption(mat, MAT_SPD,PETSC_TRUE)); 403 } else { 404 PetscCall(MatGetOption(Cmat, MAT_HERMITIAN, &flg)); 405 if (flg) { 406 PetscCall(MatSetOption(mat, MAT_HERMITIAN,PETSC_TRUE)); 407 } else { 408 #if !defined(PETSC_USE_COMPLEX) 409 PetscCall(MatGetOption(Cmat, MAT_SYMMETRIC, &flg)); 410 if (flg) PetscCall(MatSetOption(mat, MAT_SYMMETRIC,PETSC_TRUE)); 411 #endif 412 } 413 } 414 *a_Amat_crs = mat; 415 } 416 PetscCall(MatDestroy(&Cmat)); 417 418 /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0)); */ 419 /* prolongator */ 420 { 421 IS findices; 422 PetscInt Istart,Iend; 423 Mat Pnew; 424 425 PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend)); 426 /* PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 427 PetscCall(ISCreateStride(comm,Iend-Istart,Istart,1,&findices)); 428 PetscCall(ISSetBlockSize(findices,f_bs)); 429 PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew)); 430 PetscCall(ISDestroy(&findices)); 431 PetscCall(MatSetOption(Pnew,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE)); 432 433 /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0)); */ 434 PetscCall(MatDestroy(a_P_inout)); 435 436 /* output - repartitioned */ 437 *a_P_inout = Pnew; 438 } 439 PetscCall(ISDestroy(&new_eq_indices)); 440 441 *a_nactive_proc = new_size; /* output */ 442 443 /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */ 444 if (pc_gamg->cpu_pin_coarse_grids) { 445 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 446 static PetscInt llev = 2; 447 PetscCall(PetscInfo(pc,"%s: Pinning level %" PetscInt_FMT " to the CPU\n",((PetscObject)pc)->prefix,llev++)); 448 #endif 449 PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE)); 450 PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE)); 451 if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */ 452 Mat A = *a_Amat_crs, P = *a_P_inout; 453 PetscMPIInt size; 454 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 455 if (size > 1) { 456 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data; 457 PetscCall(VecBindToCPU(a->lvec,PETSC_TRUE)); 458 PetscCall(VecBindToCPU(p->lvec,PETSC_TRUE)); 459 } 460 } 461 } 462 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE],0,0,0,0)); 463 } // processor reduce 464 PetscFunctionReturn(0); 465 } 466 467 // used in GEO 468 PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2) 469 { 470 const char *prefix; 471 char addp[32]; 472 PC_MG *mg = (PC_MG*)a_pc->data; 473 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 474 475 PetscFunctionBegin; 476 PetscCall(PCGetOptionsPrefix(a_pc,&prefix)); 477 PetscCall(PetscInfo(a_pc,"%s: Square Graph on level %" PetscInt_FMT "\n",((PetscObject)a_pc)->prefix,pc_gamg->current_level+1)); 478 PetscCall(MatProductCreate(Gmat1,Gmat1,NULL,Gmat2)); 479 PetscCall(MatSetOptionsPrefix(*Gmat2,prefix)); 480 PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%" PetscInt_FMT "_",pc_gamg->current_level)); 481 PetscCall(MatAppendOptionsPrefix(*Gmat2,addp)); 482 if ((*Gmat2)->structurally_symmetric) { 483 PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AB)); 484 } else { 485 PetscCall(MatSetOption(Gmat1,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE)); 486 PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AtB)); 487 } 488 PetscCall(MatProductSetFromOptions(*Gmat2)); 489 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0)); 490 PetscCall(MatProductSymbolic(*Gmat2)); 491 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0)); 492 PetscCall(MatProductClear(*Gmat2)); 493 /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */ 494 (*Gmat2)->assembled = PETSC_TRUE; 495 PetscFunctionReturn(0); 496 } 497 498 /* -------------------------------------------------------------------------- */ 499 /* 500 PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 501 by setting data structures and options. 502 503 Input Parameter: 504 . pc - the preconditioner context 505 506 */ 507 PetscErrorCode PCSetUp_GAMG(PC pc) 508 { 509 PC_MG *mg = (PC_MG*)pc->data; 510 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 511 Mat Pmat = pc->pmat; 512 PetscInt fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS]; 513 MPI_Comm comm; 514 PetscMPIInt rank,size,nactivepe; 515 Mat Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS]; 516 IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 517 PetscLogDouble nnz0=0.,nnztot=0.; 518 MatInfo info; 519 PetscBool is_last = PETSC_FALSE; 520 521 PetscFunctionBegin; 522 PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 523 PetscCallMPI(MPI_Comm_rank(comm,&rank)); 524 PetscCallMPI(MPI_Comm_size(comm,&size)); 525 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0)); 526 if (pc->setupcalled) { 527 if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) { 528 /* reset everything */ 529 PetscCall(PCReset_MG(pc)); 530 pc->setupcalled = 0; 531 } else { 532 PC_MG_Levels **mglevels = mg->levels; 533 /* just do Galerkin grids */ 534 Mat B,dA,dB; 535 if (pc_gamg->Nlevels > 1) { 536 PetscInt gl; 537 /* currently only handle case where mat and pmat are the same on coarser levels */ 538 PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB)); 539 /* (re)set to get dirty flag */ 540 PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB)); 541 542 for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) { 543 MatReuse reuse = MAT_INITIAL_MATRIX ; 544 #if defined(GAMG_STAGES) 545 PetscCall(PetscLogStagePush(gamg_stages[gl])); 546 #endif 547 /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 548 PetscCall(KSPGetOperators(mglevels[level]->smoothd,NULL,&B)); 549 if (B->product) { 550 if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) { 551 reuse = MAT_REUSE_MATRIX; 552 } 553 } 554 if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A)); 555 if (reuse == MAT_REUSE_MATRIX) { 556 PetscCall(PetscInfo(pc,"%s: RAP after first solve, reuse matrix level %" PetscInt_FMT "\n",((PetscObject)pc)->prefix,level)); 557 } else { 558 PetscCall(PetscInfo(pc,"%s: RAP after first solve, new matrix level %" PetscInt_FMT "\n",((PetscObject)pc)->prefix,level)); 559 } 560 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0)); 561 PetscCall(MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B)); 562 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0)); 563 if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B; 564 PetscCall(KSPSetOperators(mglevels[level]->smoothd,B,B)); 565 dB = B; 566 #if defined(GAMG_STAGES) 567 PetscCall(PetscLogStagePop()); 568 #endif 569 } 570 } 571 572 PetscCall(PCSetUp_MG(pc)); 573 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0)); 574 PetscFunctionReturn(0); 575 } 576 } 577 578 if (!pc_gamg->data) { 579 if (pc_gamg->orig_data) { 580 PetscCall(MatGetBlockSize(Pmat, &bs)); 581 PetscCall(MatGetLocalSize(Pmat, &qq, NULL)); 582 583 pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 584 pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 585 pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 586 587 PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data)); 588 for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 589 } else { 590 PetscCheck(pc_gamg->ops->createdefaultdata,comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 591 PetscCall(pc_gamg->ops->createdefaultdata(pc,Pmat)); 592 } 593 } 594 595 /* cache original data for reuse */ 596 if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 597 PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data)); 598 for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 599 pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 600 pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 601 } 602 603 /* get basic dims */ 604 PetscCall(MatGetBlockSize(Pmat, &bs)); 605 PetscCall(MatGetSize(Pmat, &M, &N)); 606 607 PetscCall(MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info)); /* global reduction */ 608 nnz0 = info.nz_used; 609 nnztot = info.nz_used; 610 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)); 611 612 /* Get A_i and R_i */ 613 for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) { 614 pc_gamg->current_level = level; 615 PetscCheck(level < PETSC_MG_MAXLEVELS,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %" PetscInt_FMT,level); 616 level1 = level + 1; 617 #if defined(GAMG_STAGES) 618 if (!gamg_stages[level]) { 619 char str[32]; 620 sprintf(str,"GAMG Level %d",(int)level); 621 PetscCall(PetscLogStageRegister(str, &gamg_stages[level])); 622 } 623 PetscCall(PetscLogStagePush(gamg_stages[level])); 624 #endif 625 { /* construct prolongator */ 626 Mat Gmat; 627 PetscCoarsenData *agg_lists; 628 Mat Prol11; 629 630 PetscCall(pc_gamg->ops->graph(pc,Aarr[level], &Gmat)); 631 PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists)); 632 PetscCall(pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11)); 633 634 /* could have failed to create new level */ 635 if (Prol11) { 636 const char *prefix; 637 char addp[32]; 638 639 /* get new block size of coarse matrices */ 640 PetscCall(MatGetBlockSizes(Prol11, NULL, &bs)); 641 642 if (pc_gamg->ops->optprolongator) { 643 /* smooth */ 644 PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11)); 645 } 646 647 if (pc_gamg->use_aggs_in_asm) { 648 PetscInt bs; 649 PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // not timed directly, ugly, could remove, but good ASM method 650 PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level])); 651 } 652 653 PetscCall(PCGetOptionsPrefix(pc,&prefix)); 654 PetscCall(MatSetOptionsPrefix(Prol11,prefix)); 655 PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level)); 656 PetscCall(MatAppendOptionsPrefix(Prol11,addp)); 657 /* Always generate the transpose with CUDA 658 Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */ 659 PetscCall(MatSetOption(Prol11,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE)); 660 PetscCall(MatSetFromOptions(Prol11)); 661 Parr[level1] = Prol11; 662 } else Parr[level1] = NULL; /* failed to coarsen */ 663 664 PetscCall(MatDestroy(&Gmat)); 665 PetscCall(PetscCDDestroy(agg_lists)); 666 } /* construct prolongator scope */ 667 if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 668 if (!Parr[level1]) { /* failed to coarsen */ 669 PetscCall(PetscInfo(pc,"%s: Stop gridding, level %" PetscInt_FMT "\n",((PetscObject)pc)->prefix,level)); 670 #if defined(GAMG_STAGES) 671 PetscCall(PetscLogStagePop()); 672 #endif 673 break; 674 } 675 PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */ 676 PetscCheck(!is_last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?"); 677 if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 678 if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE; 679 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL],0,0,0,0)); 680 PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last)); 681 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL],0,0,0,0)); 682 683 PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */ 684 PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info)); 685 nnztot += info.nz_used; 686 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)); 687 688 #if defined(GAMG_STAGES) 689 PetscCall(PetscLogStagePop()); 690 #endif 691 /* stop if one node or one proc -- could pull back for singular problems */ 692 if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) { 693 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)); 694 level++; 695 break; 696 } 697 } /* levels */ 698 PetscCall(PetscFree(pc_gamg->data)); 699 700 PetscCall(PetscInfo(pc,"%s: %" PetscInt_FMT " levels, grid complexity = %g\n",((PetscObject)pc)->prefix,level+1,nnztot/nnz0)); 701 pc_gamg->Nlevels = level + 1; 702 fine_level = level; 703 PetscCall(PCMGSetLevels(pc,pc_gamg->Nlevels,NULL)); 704 705 if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 706 707 /* set default smoothers & set operators */ 708 for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) { 709 KSP smoother; 710 PC subpc; 711 712 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 713 PetscCall(KSPGetPC(smoother, &subpc)); 714 715 PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 716 /* set ops */ 717 PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level])); 718 PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level+1])); 719 720 /* set defaults */ 721 PetscCall(KSPSetType(smoother, KSPCHEBYSHEV)); 722 723 /* set blocks for ASM smoother that uses the 'aggregates' */ 724 if (pc_gamg->use_aggs_in_asm) { 725 PetscInt sz; 726 IS *iss; 727 728 sz = nASMBlocksArr[level]; 729 iss = ASMLocalIDsArr[level]; 730 PetscCall(PCSetType(subpc, PCASM)); 731 PetscCall(PCASMSetOverlap(subpc, 0)); 732 PetscCall(PCASMSetType(subpc,PC_ASM_BASIC)); 733 if (!sz) { 734 IS is; 735 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is)); 736 PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is)); 737 PetscCall(ISDestroy(&is)); 738 } else { 739 PetscInt kk; 740 PetscCall(PCASMSetLocalSubdomains(subpc, sz, NULL, iss)); 741 for (kk=0; kk<sz; kk++) { 742 PetscCall(ISDestroy(&iss[kk])); 743 } 744 PetscCall(PetscFree(iss)); 745 } 746 ASMLocalIDsArr[level] = NULL; 747 nASMBlocksArr[level] = 0; 748 } else { 749 PetscCall(PCSetType(subpc, PCJACOBI)); 750 } 751 } 752 { 753 /* coarse grid */ 754 KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 755 Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 756 757 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 758 PetscCall(KSPSetOperators(smoother, Lmat, Lmat)); 759 if (!pc_gamg->use_parallel_coarse_grid_solver) { 760 PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE)); 761 PetscCall(KSPGetPC(smoother, &subpc)); 762 PetscCall(PCSetType(subpc, PCBJACOBI)); 763 PetscCall(PCSetUp(subpc)); 764 PetscCall(PCBJacobiGetSubKSP(subpc,&ii,&first,&k2)); 765 PetscCheck(ii == 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %" PetscInt_FMT " is not one",ii); 766 PetscCall(KSPGetPC(k2[0],&pc2)); 767 PetscCall(PCSetType(pc2, PCLU)); 768 PetscCall(PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS)); 769 PetscCall(KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1)); 770 PetscCall(KSPSetType(k2[0], KSPPREONLY)); 771 } 772 } 773 774 /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 775 PetscObjectOptionsBegin((PetscObject)pc); 776 PetscCall(PCSetFromOptions_MG(PetscOptionsObject,pc)); 777 PetscOptionsEnd(); 778 PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL)); 779 780 /* setup cheby eigen estimates from SA */ 781 if (pc_gamg->use_sa_esteig) { 782 for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) { 783 KSP smoother; 784 PetscBool ischeb; 785 786 PetscCall(PCMGGetSmoother(pc, lidx, &smoother)); 787 PetscCall(PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb)); 788 if (ischeb) { 789 KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data; 790 791 // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG 792 if (mg->max_eigen_DinvA[level] > 0) { 793 // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it. 794 // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.) 795 PetscReal emax,emin; 796 797 emin = mg->min_eigen_DinvA[level]; 798 emax = mg->max_eigen_DinvA[level]; 799 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)); 800 cheb->emin_provided = emin; 801 cheb->emax_provided = emax; 802 } 803 } 804 } 805 } 806 807 PetscCall(PCSetUp_MG(pc)); 808 809 /* clean up */ 810 for (level=1; level<pc_gamg->Nlevels; level++) { 811 PetscCall(MatDestroy(&Parr[level])); 812 PetscCall(MatDestroy(&Aarr[level])); 813 } 814 } else { 815 KSP smoother; 816 817 PetscCall(PetscInfo(pc,"%s: One level solver used (system is seen as DD). Using default solver.\n",((PetscObject)pc)->prefix)); 818 PetscCall(PCMGGetSmoother(pc, 0, &smoother)); 819 PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0])); 820 PetscCall(KSPSetType(smoother, KSPPREONLY)); 821 PetscCall(PCSetUp_MG(pc)); 822 } 823 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP],0,0,0,0)); 824 PetscFunctionReturn(0); 825 } 826 827 /* ------------------------------------------------------------------------- */ 828 /* 829 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 830 that was created with PCCreate_GAMG(). 831 832 Input Parameter: 833 . pc - the preconditioner context 834 835 Application Interface Routine: PCDestroy() 836 */ 837 PetscErrorCode PCDestroy_GAMG(PC pc) 838 { 839 PC_MG *mg = (PC_MG*)pc->data; 840 PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 841 842 PetscFunctionBegin; 843 PetscCall(PCReset_GAMG(pc)); 844 if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc)); 845 PetscCall(PetscFree(pc_gamg->ops)); 846 PetscCall(PetscFree(pc_gamg->gamg_type_name)); 847 PetscCall(PetscFree(pc_gamg)); 848 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",NULL)); 849 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",NULL)); 850 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",NULL)); 851 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",NULL)); 852 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",NULL)); 853 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",NULL)); 854 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",NULL)); 855 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",NULL)); 856 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",NULL)); 857 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",NULL)); 858 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",NULL)); 859 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",NULL)); 860 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",NULL)); 861 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",NULL)); 862 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",NULL)); 863 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",NULL)); 864 PetscCall(PCDestroy_MG(pc)); 865 PetscFunctionReturn(0); 866 } 867 868 /*@ 869 PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction. 870 871 Logically Collective on PC 872 873 Input Parameters: 874 + pc - the preconditioner context 875 - n - the number of equations 876 877 Options Database Key: 878 . -pc_gamg_process_eq_limit <limit> - set the limit 879 880 Notes: 881 GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process 882 that has degrees of freedom 883 884 Level: intermediate 885 886 .seealso: `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()` 887 @*/ 888 PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 889 { 890 PetscFunctionBegin; 891 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 892 PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n)); 893 PetscFunctionReturn(0); 894 } 895 896 static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 897 { 898 PC_MG *mg = (PC_MG*)pc->data; 899 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 900 901 PetscFunctionBegin; 902 if (n>0) pc_gamg->min_eq_proc = n; 903 PetscFunctionReturn(0); 904 } 905 906 /*@ 907 PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid. 908 909 Collective on PC 910 911 Input Parameters: 912 + pc - the preconditioner context 913 - n - maximum number of equations to aim for 914 915 Options Database Key: 916 . -pc_gamg_coarse_eq_limit <limit> - set the limit 917 918 Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 919 has less than 1000 unknowns. 920 921 Level: intermediate 922 923 .seealso: `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()` 924 @*/ 925 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 926 { 927 PetscFunctionBegin; 928 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 929 PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n)); 930 PetscFunctionReturn(0); 931 } 932 933 static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 934 { 935 PC_MG *mg = (PC_MG*)pc->data; 936 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 937 938 PetscFunctionBegin; 939 if (n>0) pc_gamg->coarse_eq_limit = n; 940 PetscFunctionReturn(0); 941 } 942 943 /*@ 944 PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids 945 946 Collective on PC 947 948 Input Parameters: 949 + pc - the preconditioner context 950 - n - PETSC_TRUE or PETSC_FALSE 951 952 Options Database Key: 953 . -pc_gamg_repartition <true,false> - turn on the repartitioning 954 955 Notes: 956 this will generally improve the loading balancing of the work on each level 957 958 Level: intermediate 959 960 @*/ 961 PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) 962 { 963 PetscFunctionBegin; 964 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 965 PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n)); 966 PetscFunctionReturn(0); 967 } 968 969 static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) 970 { 971 PC_MG *mg = (PC_MG*)pc->data; 972 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 973 974 PetscFunctionBegin; 975 pc_gamg->repart = n; 976 PetscFunctionReturn(0); 977 } 978 979 /*@ 980 PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Chebyshev smoother 981 982 Collective on PC 983 984 Input Parameters: 985 + pc - the preconditioner context 986 - n - number of its 987 988 Options Database Key: 989 . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate 990 991 Notes: 992 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$. 993 Eigenvalue estimates (based on a few CG or GMRES iterations) are computed to choose $\omega$ so that this is a stable smoothing operation. 994 If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused. 995 This option is only used when the smoother uses Jacobi, and should be turned off if a different PCJacobiType is used. 996 It became default in PETSc 3.17. 997 998 Level: advanced 999 1000 .seealso: `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()` 1001 @*/ 1002 PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n) 1003 { 1004 PetscFunctionBegin; 1005 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1006 PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n)); 1007 PetscFunctionReturn(0); 1008 } 1009 1010 static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n) 1011 { 1012 PC_MG *mg = (PC_MG*)pc->data; 1013 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1014 1015 PetscFunctionBegin; 1016 pc_gamg->use_sa_esteig = n; 1017 PetscFunctionReturn(0); 1018 } 1019 1020 /*@ 1021 PCGAMGSetEigenvalues - Set eigenvalues 1022 1023 Collective on PC 1024 1025 Input Parameters: 1026 + pc - the preconditioner context 1027 - emax - max eigenvalue 1028 . emin - min eigenvalue 1029 1030 Options Database Key: 1031 . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues 1032 1033 Level: intermediate 1034 1035 .seealso: `PCGAMGSetUseSAEstEig()` 1036 @*/ 1037 PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin) 1038 { 1039 PetscFunctionBegin; 1040 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1041 PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin)); 1042 PetscFunctionReturn(0); 1043 } 1044 1045 static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin) 1046 { 1047 PC_MG *mg = (PC_MG*)pc->data; 1048 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1049 1050 PetscFunctionBegin; 1051 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); 1052 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); 1053 pc_gamg->emax = emax; 1054 pc_gamg->emin = emin; 1055 PetscFunctionReturn(0); 1056 } 1057 1058 /*@ 1059 PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner 1060 1061 Collective on PC 1062 1063 Input Parameters: 1064 + pc - the preconditioner context 1065 - n - PETSC_TRUE or PETSC_FALSE 1066 1067 Options Database Key: 1068 . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation 1069 1070 Level: intermediate 1071 1072 Notes: 1073 May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1074 rebuilding the preconditioner quicker. 1075 1076 @*/ 1077 PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 1078 { 1079 PetscFunctionBegin; 1080 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1081 PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n)); 1082 PetscFunctionReturn(0); 1083 } 1084 1085 static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 1086 { 1087 PC_MG *mg = (PC_MG*)pc->data; 1088 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1089 1090 PetscFunctionBegin; 1091 pc_gamg->reuse_prol = n; 1092 PetscFunctionReturn(0); 1093 } 1094 1095 /*@ 1096 PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner. 1097 1098 Collective on PC 1099 1100 Input Parameters: 1101 + pc - the preconditioner context 1102 - flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not 1103 1104 Options Database Key: 1105 . -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains 1106 1107 Level: intermediate 1108 1109 @*/ 1110 PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) 1111 { 1112 PetscFunctionBegin; 1113 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1114 PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg)); 1115 PetscFunctionReturn(0); 1116 } 1117 1118 static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) 1119 { 1120 PC_MG *mg = (PC_MG*)pc->data; 1121 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1122 1123 PetscFunctionBegin; 1124 pc_gamg->use_aggs_in_asm = flg; 1125 PetscFunctionReturn(0); 1126 } 1127 1128 /*@ 1129 PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver 1130 1131 Collective on PC 1132 1133 Input Parameters: 1134 + pc - the preconditioner context 1135 - flg - PETSC_TRUE to not force coarse grid onto one processor 1136 1137 Options Database Key: 1138 . -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver 1139 1140 Level: intermediate 1141 1142 .seealso: `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()` 1143 @*/ 1144 PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg) 1145 { 1146 PetscFunctionBegin; 1147 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1148 PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg)); 1149 PetscFunctionReturn(0); 1150 } 1151 1152 static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) 1153 { 1154 PC_MG *mg = (PC_MG*)pc->data; 1155 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1156 1157 PetscFunctionBegin; 1158 pc_gamg->use_parallel_coarse_grid_solver = flg; 1159 PetscFunctionReturn(0); 1160 } 1161 1162 /*@ 1163 PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU 1164 1165 Collective on PC 1166 1167 Input Parameters: 1168 + pc - the preconditioner context 1169 - flg - PETSC_TRUE to pin coarse grids to CPU 1170 1171 Options Database Key: 1172 . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU 1173 1174 Level: intermediate 1175 1176 .seealso: `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetUseParallelCoarseGridSolve()` 1177 @*/ 1178 PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) 1179 { 1180 PetscFunctionBegin; 1181 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1182 PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg)); 1183 PetscFunctionReturn(0); 1184 } 1185 1186 static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) 1187 { 1188 PC_MG *mg = (PC_MG*)pc->data; 1189 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1190 1191 PetscFunctionBegin; 1192 pc_gamg->cpu_pin_coarse_grids = flg; 1193 PetscFunctionReturn(0); 1194 } 1195 1196 /*@ 1197 PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type) 1198 1199 Collective on PC 1200 1201 Input Parameters: 1202 + pc - the preconditioner context 1203 - flg - Layout type 1204 1205 Options Database Key: 1206 . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering 1207 1208 Level: intermediate 1209 1210 .seealso: `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()` 1211 @*/ 1212 PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) 1213 { 1214 PetscFunctionBegin; 1215 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1216 PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg)); 1217 PetscFunctionReturn(0); 1218 } 1219 1220 static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) 1221 { 1222 PC_MG *mg = (PC_MG*)pc->data; 1223 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1224 1225 PetscFunctionBegin; 1226 pc_gamg->layout_type = flg; 1227 PetscFunctionReturn(0); 1228 } 1229 1230 /*@ 1231 PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use 1232 1233 Not collective on PC 1234 1235 Input Parameters: 1236 + pc - the preconditioner 1237 - n - the maximum number of levels to use 1238 1239 Options Database Key: 1240 . -pc_mg_levels <n> - set the maximum number of levels to allow 1241 1242 Level: intermediate 1243 1244 @*/ 1245 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 1246 { 1247 PetscFunctionBegin; 1248 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1249 PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n)); 1250 PetscFunctionReturn(0); 1251 } 1252 1253 static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 1254 { 1255 PC_MG *mg = (PC_MG*)pc->data; 1256 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1257 1258 PetscFunctionBegin; 1259 pc_gamg->Nlevels = n; 1260 PetscFunctionReturn(0); 1261 } 1262 1263 /*@ 1264 PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 1265 1266 Not collective on PC 1267 1268 Input Parameters: 1269 + pc - the preconditioner context 1270 . 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 1271 - n - number of threshold values provided in array 1272 1273 Options Database Key: 1274 . -pc_gamg_threshold <threshold> - the threshold to drop edges 1275 1276 Notes: 1277 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. 1278 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. 1279 1280 If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening. 1281 In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold(). 1282 If n is greater than the total number of levels, the excess entries in threshold will not be used. 1283 1284 Level: intermediate 1285 1286 .seealso: `PCGAMGFilterGraph()`, `PCGAMGSetAggressiveLevels()` 1287 @*/ 1288 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) 1289 { 1290 PetscFunctionBegin; 1291 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1292 if (n) PetscValidRealPointer(v,2); 1293 PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n)); 1294 PetscFunctionReturn(0); 1295 } 1296 1297 static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) 1298 { 1299 PC_MG *mg = (PC_MG*)pc->data; 1300 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1301 PetscInt i; 1302 PetscFunctionBegin; 1303 for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i]; 1304 for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale; 1305 PetscFunctionReturn(0); 1306 } 1307 1308 /*@ 1309 PCGAMGSetRankReductionFactors - Set manual schedule for process reduction on coarse grids 1310 1311 Collective on PC 1312 1313 Input Parameters: 1314 + pc - the preconditioner context 1315 . v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda 1316 - n - number of values provided in array 1317 1318 Options Database Key: 1319 . -pc_gamg_rank_reduction_factors <factors> - provide the schedule 1320 1321 Level: intermediate 1322 1323 .seealso: `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()` 1324 @*/ 1325 PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n) 1326 { 1327 PetscFunctionBegin; 1328 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1329 if (n) PetscValidIntPointer(v,2); 1330 PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n)); 1331 PetscFunctionReturn(0); 1332 } 1333 1334 static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n) 1335 { 1336 PC_MG *mg = (PC_MG*)pc->data; 1337 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1338 PetscInt i; 1339 PetscFunctionBegin; 1340 for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i]; 1341 for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */ 1342 PetscFunctionReturn(0); 1343 } 1344 1345 /*@ 1346 PCGAMGSetThresholdScale - Relative threshold reduction at each level 1347 1348 Not collective on PC 1349 1350 Input Parameters: 1351 + pc - the preconditioner context 1352 - scale - the threshold value reduction, usually < 1.0 1353 1354 Options Database Key: 1355 . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level 1356 1357 Notes: 1358 The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold(). 1359 This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold(). 1360 1361 Level: advanced 1362 1363 .seealso: `PCGAMGSetThreshold()` 1364 @*/ 1365 PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) 1366 { 1367 PetscFunctionBegin; 1368 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1369 PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v)); 1370 PetscFunctionReturn(0); 1371 } 1372 1373 static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) 1374 { 1375 PC_MG *mg = (PC_MG*)pc->data; 1376 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1377 PetscFunctionBegin; 1378 pc_gamg->threshold_scale = v; 1379 PetscFunctionReturn(0); 1380 } 1381 1382 /*@C 1383 PCGAMGSetType - Set solution method 1384 1385 Collective on PC 1386 1387 Input Parameters: 1388 + pc - the preconditioner context 1389 - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1390 1391 Options Database Key: 1392 . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply 1393 1394 Level: intermediate 1395 1396 .seealso: `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType` 1397 @*/ 1398 PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1399 { 1400 PetscFunctionBegin; 1401 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1402 PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type)); 1403 PetscFunctionReturn(0); 1404 } 1405 1406 /*@C 1407 PCGAMGGetType - Get solution method 1408 1409 Collective on PC 1410 1411 Input Parameter: 1412 . pc - the preconditioner context 1413 1414 Output Parameter: 1415 . type - the type of algorithm used 1416 1417 Level: intermediate 1418 1419 .seealso: `PCGAMGSetType()`, `PCGAMGType` 1420 @*/ 1421 PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1422 { 1423 PetscFunctionBegin; 1424 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1425 PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type)); 1426 PetscFunctionReturn(0); 1427 } 1428 1429 static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1430 { 1431 PC_MG *mg = (PC_MG*)pc->data; 1432 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1433 1434 PetscFunctionBegin; 1435 *type = pc_gamg->type; 1436 PetscFunctionReturn(0); 1437 } 1438 1439 static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1440 { 1441 PC_MG *mg = (PC_MG*)pc->data; 1442 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1443 PetscErrorCode (*r)(PC); 1444 1445 PetscFunctionBegin; 1446 pc_gamg->type = type; 1447 PetscCall(PetscFunctionListFind(GAMGList,type,&r)); 1448 PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 1449 if (pc_gamg->ops->destroy) { 1450 PetscCall((*pc_gamg->ops->destroy)(pc)); 1451 PetscCall(PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps))); 1452 pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1453 /* cleaning up common data in pc_gamg - this should disapear someday */ 1454 pc_gamg->data_cell_cols = 0; 1455 pc_gamg->data_cell_rows = 0; 1456 pc_gamg->orig_data_cell_cols = 0; 1457 pc_gamg->orig_data_cell_rows = 0; 1458 PetscCall(PetscFree(pc_gamg->data)); 1459 pc_gamg->data_sz = 0; 1460 } 1461 PetscCall(PetscFree(pc_gamg->gamg_type_name)); 1462 PetscCall(PetscStrallocpy(type,&pc_gamg->gamg_type_name)); 1463 PetscCall((*r)(pc)); 1464 PetscFunctionReturn(0); 1465 } 1466 1467 static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer) 1468 { 1469 PC_MG *mg = (PC_MG*)pc->data; 1470 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1471 PetscReal gc=0, oc=0; 1472 1473 PetscFunctionBegin; 1474 PetscCall(PetscViewerASCIIPrintf(viewer," GAMG specific options\n")); 1475 PetscCall(PetscViewerASCIIPrintf(viewer," Threshold for dropping small values in graph on each level =")); 1476 for (PetscInt i=0;i<mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i])); 1477 PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); 1478 PetscCall(PetscViewerASCIIPrintf(viewer," Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale)); 1479 if (pc_gamg->use_aggs_in_asm) { 1480 PetscCall(PetscViewerASCIIPrintf(viewer," Using aggregates from coarsening process to define subdomains for PCASM\n")); 1481 } 1482 if (pc_gamg->use_parallel_coarse_grid_solver) { 1483 PetscCall(PetscViewerASCIIPrintf(viewer," Using parallel coarse grid solver (all coarse grid equations not put on one process)\n")); 1484 } 1485 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 1486 if (pc_gamg->cpu_pin_coarse_grids) { 1487 /* PetscCall(PetscViewerASCIIPrintf(viewer," Pinning coarse grids to the CPU)\n")); */ 1488 } 1489 #endif 1490 /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */ 1491 /* PetscCall(PetscViewerASCIIPrintf(viewer," Put reduced grids on processes in natural order (ie, 0,1,2...)\n")); */ 1492 /* } else { */ 1493 /* PetscCall(PetscViewerASCIIPrintf(viewer," Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n")); */ 1494 /* } */ 1495 if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc,viewer)); 1496 PetscCall(PCMGGetGridComplexity(pc,&gc,&oc)); 1497 PetscCall(PetscViewerASCIIPrintf(viewer," Complexity: grid = %g operator = %g\n",(double)gc,(double)oc)); 1498 PetscFunctionReturn(0); 1499 } 1500 1501 PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc) 1502 { 1503 PC_MG *mg = (PC_MG*)pc->data; 1504 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1505 PetscBool flag; 1506 MPI_Comm comm; 1507 char prefix[256],tname[32]; 1508 PetscInt i,n; 1509 const char *pcpre; 1510 static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL}; 1511 PetscFunctionBegin; 1512 PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 1513 PetscOptionsHeadBegin(PetscOptionsObject,"GAMG options"); 1514 PetscCall(PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag)); 1515 if (flag) PetscCall(PCGAMGSetType(pc,tname)); 1516 PetscCall(PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL)); 1517 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)); 1518 PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL)); 1519 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)); 1520 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)); 1521 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)); 1522 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)); 1523 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)); 1524 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)); 1525 PetscCall(PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL)); 1526 n = PETSC_MG_MAXLEVELS; 1527 PetscCall(PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag)); 1528 if (!flag || n < PETSC_MG_MAXLEVELS) { 1529 if (!flag) n = 1; 1530 i = n; 1531 do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS); 1532 } 1533 n = PETSC_MG_MAXLEVELS; 1534 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)); 1535 if (!flag) i = 0; 1536 else i = n; 1537 do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS); 1538 PetscCall(PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL)); 1539 { 1540 PetscReal eminmax[2] = {0., 0.}; 1541 n = 2; 1542 PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag)); 1543 if (flag) { 1544 PetscCheck(n == 2,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 1545 PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0])); 1546 } 1547 } 1548 /* set options for subtype */ 1549 if (pc_gamg->ops->setfromoptions) PetscCall((*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc)); 1550 1551 PetscCall(PCGetOptionsPrefix(pc, &pcpre)); 1552 PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "")); 1553 PetscOptionsHeadEnd(); 1554 PetscFunctionReturn(0); 1555 } 1556 1557 /* -------------------------------------------------------------------------- */ 1558 /*MC 1559 PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 1560 1561 Options Database Keys: 1562 + -pc_gamg_type <type> - one of agg, geo, or classical 1563 . -pc_gamg_repartition <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined 1564 - -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations 1565 + -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 1566 . -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> 1567 equations on each process that has degrees of freedom 1568 - -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 1569 + -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level 1570 . -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1571 1572 Options Database Keys for default Aggregation: 1573 + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 1574 . -pc_gamg_symmetrize_graph <true,default=false> - symmetrize the graph before computing the aggregation 1575 - -pc_gamg_square_graph <n,default=1> - alias for pc_gamg_aggressive_coarsening (deprecated) 1576 - -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest. 1577 1578 Multigrid options: 1579 + -pc_mg_cycles <v> - v or w, see PCMGSetCycleType() 1580 . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1581 . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1582 - -pc_mg_levels <levels> - Number of levels of multigrid to use. 1583 1584 Notes: 1585 In order to obtain good performance for PCGAMG for vector valued problems you must 1586 Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point 1587 Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator 1588 See the Users Manual Chapter 4 for more details 1589 1590 Level: intermediate 1591 1592 .seealso: `PCCreate()`, `PCSetType()`, `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, 1593 `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGSetUseSAEstEig()` 1594 M*/ 1595 1596 PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 1597 { 1598 PC_GAMG *pc_gamg; 1599 PC_MG *mg; 1600 1601 PetscFunctionBegin; 1602 /* register AMG type */ 1603 PetscCall(PCGAMGInitializePackage()); 1604 1605 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 1606 PetscCall(PCSetType(pc, PCMG)); 1607 PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG)); 1608 1609 /* create a supporting struct and attach it to pc */ 1610 PetscCall(PetscNewLog(pc,&pc_gamg)); 1611 PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL)); 1612 mg = (PC_MG*)pc->data; 1613 mg->innerctx = pc_gamg; 1614 1615 PetscCall(PetscNewLog(pc,&pc_gamg->ops)); 1616 1617 /* these should be in subctx but repartitioning needs simple arrays */ 1618 pc_gamg->data_sz = 0; 1619 pc_gamg->data = NULL; 1620 1621 /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 1622 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 1623 pc->ops->setup = PCSetUp_GAMG; 1624 pc->ops->reset = PCReset_GAMG; 1625 pc->ops->destroy = PCDestroy_GAMG; 1626 mg->view = PCView_GAMG; 1627 1628 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG)); 1629 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG)); 1630 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG)); 1631 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG)); 1632 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG)); 1633 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG)); 1634 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG)); 1635 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG)); 1636 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG)); 1637 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG)); 1638 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG)); 1639 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG)); 1640 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG)); 1641 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG)); 1642 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG)); 1643 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG)); 1644 pc_gamg->repart = PETSC_FALSE; 1645 pc_gamg->reuse_prol = PETSC_FALSE; 1646 pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1647 pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1648 pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1649 pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1650 pc_gamg->min_eq_proc = 50; 1651 pc_gamg->coarse_eq_limit = 50; 1652 PetscCall(PetscArrayzero(pc_gamg->threshold,PETSC_MG_MAXLEVELS)); 1653 pc_gamg->threshold_scale = 1.; 1654 pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 1655 pc_gamg->current_level = 0; /* don't need to init really */ 1656 pc_gamg->use_sa_esteig = PETSC_TRUE; 1657 pc_gamg->emin = 0; 1658 pc_gamg->emax = 0; 1659 1660 pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1661 1662 /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1663 PetscCall(PCGAMGSetType(pc,PCGAMGAGG)); 1664 PetscFunctionReturn(0); 1665 } 1666 1667 /*@C 1668 PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 1669 from PCInitializePackage(). 1670 1671 Level: developer 1672 1673 .seealso: `PetscInitialize()` 1674 @*/ 1675 PetscErrorCode PCGAMGInitializePackage(void) 1676 { 1677 PetscInt l; 1678 1679 PetscFunctionBegin; 1680 if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 1681 PCGAMGPackageInitialized = PETSC_TRUE; 1682 PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO)); 1683 PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG)); 1684 PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical)); 1685 PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage)); 1686 1687 /* general events */ 1688 PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP])); 1689 PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH])); 1690 PetscCall(PetscLogEventRegister(" PCGAMGFilter", PC_CLASSID, &petsc_gamg_setup_events[GAMG_FILTER])); 1691 PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN])); 1692 PetscCall(PetscLogEventRegister(" GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS])); 1693 PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL])); 1694 PetscCall(PetscLogEventRegister(" GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA])); 1695 PetscCall(PetscLogEventRegister(" GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB])); 1696 PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT])); 1697 PetscCall(PetscLogEventRegister(" GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM])); 1698 PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL])); 1699 PetscCall(PetscLogEventRegister(" GAMG PtAP",PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP])); 1700 PetscCall(PetscLogEventRegister(" GAMG Reduce",PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE])); 1701 PetscCall(PetscLogEventRegister(" GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART])); 1702 /* PetscCall(PetscLogEventRegister(" GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */ 1703 /* PetscCall(PetscLogEventRegister(" GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */ 1704 /* PetscCall(PetscLogEventRegister(" GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */ 1705 for (l=0;l<PETSC_MG_MAXLEVELS;l++) { 1706 char ename[32]; 1707 1708 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02" PetscInt_FMT,l)); 1709 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0])); 1710 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02" PetscInt_FMT,l)); 1711 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1])); 1712 PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02" PetscInt_FMT,l)); 1713 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2])); 1714 } 1715 #if defined(GAMG_STAGES) 1716 { /* create timer stages */ 1717 char str[32]; 1718 sprintf(str,"GAMG Level %d",0); 1719 PetscCall(PetscLogStageRegister(str, &gamg_stages[0])); 1720 } 1721 #endif 1722 PetscFunctionReturn(0); 1723 } 1724 1725 /*@C 1726 PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is 1727 called from PetscFinalize() automatically. 1728 1729 Level: developer 1730 1731 .seealso: `PetscFinalize()` 1732 @*/ 1733 PetscErrorCode PCGAMGFinalizePackage(void) 1734 { 1735 PetscFunctionBegin; 1736 PCGAMGPackageInitialized = PETSC_FALSE; 1737 PetscCall(PetscFunctionListDestroy(&GAMGList)); 1738 PetscFunctionReturn(0); 1739 } 1740 1741 /*@C 1742 PCGAMGRegister - Register a PCGAMG implementation. 1743 1744 Input Parameters: 1745 + type - string that will be used as the name of the GAMG type. 1746 - create - function for creating the gamg context. 1747 1748 Level: advanced 1749 1750 .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()` 1751 @*/ 1752 PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1753 { 1754 PetscFunctionBegin; 1755 PetscCall(PCGAMGInitializePackage()); 1756 PetscCall(PetscFunctionListAdd(&GAMGList,type,create)); 1757 PetscFunctionReturn(0); 1758 } 1759