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