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 <petsc-private/kspimpl.h> 7 #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */ 8 9 #if defined PETSC_GAMG_USE_LOG 10 PetscLogEvent petsc_gamg_setup_events[NUM_SET]; 11 #endif 12 13 #if defined PETSC_USE_LOG 14 PetscLogEvent PC_GAMGGgraph_AGG; 15 PetscLogEvent PC_GAMGGgraph_GEO; 16 PetscLogEvent PC_GAMGCoarsen_AGG; 17 PetscLogEvent PC_GAMGCoarsen_GEO; 18 PetscLogEvent PC_GAMGProlongator_AGG; 19 PetscLogEvent PC_GAMGProlongator_GEO; 20 PetscLogEvent PC_GAMGOptprol_AGG; 21 #endif 22 23 #define GAMG_MAXLEVELS 30 24 25 /* #define GAMG_STAGES */ 26 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 27 static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; 28 #endif 29 30 static PetscFunctionList GAMGList = 0; 31 static PetscBool PCGAMGPackageInitialized; 32 33 /* ----------------------------------------------------------------------------- */ 34 #undef __FUNCT__ 35 #define __FUNCT__ "PCReset_GAMG" 36 PetscErrorCode PCReset_GAMG(PC pc) 37 { 38 PetscErrorCode ierr; 39 PC_MG *mg = (PC_MG*)pc->data; 40 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 41 42 PetscFunctionBegin; 43 if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */ 44 PetscPrintf(PetscObjectComm((PetscObject)pc),"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__); 45 ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 46 } 47 pc_gamg->data = NULL; pc_gamg->data_sz = 0; 48 49 if (pc_gamg->orig_data) { 50 ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr); 51 } 52 PetscFunctionReturn(0); 53 } 54 55 /* -------------------------------------------------------------------------- */ 56 /* 57 PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 58 of active processors. 59 60 Input Parameter: 61 . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 62 'pc_gamg->data_sz' are changed via repartitioning/reduction. 63 . Amat_fine - matrix on this fine (k) level 64 . cr_bs - coarse block size 65 In/Output Parameter: 66 . a_P_inout - prolongation operator to the next level (k-->k-1) 67 . a_nactive_proc - number of active procs 68 Output Parameter: 69 . a_Amat_crs - coarse matrix that is created (k-1) 70 */ 71 72 #undef __FUNCT__ 73 #define __FUNCT__ "PCGAMGCreateLevel_GAMG" 74 static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs, 75 Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc, 76 IS * Pcolumnperm) 77 { 78 PetscErrorCode ierr; 79 PC_MG *mg = (PC_MG*)pc->data; 80 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 81 Mat Cmat,Pold=*a_P_inout; 82 MPI_Comm comm; 83 PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc; 84 PetscInt ncrs_eq,ncrs,f_bs; 85 86 PetscFunctionBegin; 87 ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr); 88 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 89 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 90 ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr); 91 /* RAP */ 92 ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr); 93 94 /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 95 ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr); 96 if (pc_gamg->data_cell_rows>0) { 97 ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 98 } 99 else { 100 PetscInt bs; 101 ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr); 102 ncrs = ncrs_eq/bs; 103 } 104 105 /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 106 { 107 PetscInt ncrs_eq_glob; 108 ierr = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr); 109 new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 110 if (new_size == 0) new_size = 1; /* not likely, posible? */ 111 else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 112 } 113 114 if (Pcolumnperm) *Pcolumnperm = NULL; 115 116 if (!pc_gamg->repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 117 else { 118 PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old; 119 IS is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices; 120 121 nloc_old = ncrs_eq/cr_bs; 122 if (ncrs_eq % cr_bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs); 123 #if defined PETSC_GAMG_USE_LOG 124 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 125 #endif 126 /* make 'is_eq_newproc' */ 127 ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr); 128 if (pc_gamg->repart) { 129 /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */ 130 Mat adj; 131 132 if (pc_gamg->verbose>0) { 133 if (pc_gamg->verbose==1) PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,ncrs_eq); 134 else { 135 PetscInt n; 136 ierr = MPI_Allreduce(&ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr); 137 PetscPrintf(comm,"\t[%d]%s repartition: size (active): %d --> %d, neq = %d\n",rank,__FUNCT__,*a_nactive_proc,new_size,n); 138 } 139 } 140 141 /* get 'adj' */ 142 if (cr_bs == 1) { 143 ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 144 } else { 145 /* make a scalar matrix to partition (no Stokes here) */ 146 Mat tMat; 147 PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 148 const PetscScalar *vals; 149 const PetscInt *idx; 150 PetscInt *d_nnz, *o_nnz, M, N; 151 static PetscInt llev = 0; 152 153 ierr = PetscMalloc1(ncrs, &d_nnz);CHKERRQ(ierr); 154 ierr = PetscMalloc1(ncrs, &o_nnz);CHKERRQ(ierr); 155 ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr); 156 ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr); 157 for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 158 ierr = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 159 d_nnz[jj] = ncols/cr_bs; 160 o_nnz[jj] = ncols/cr_bs; 161 ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 162 if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 163 if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs; 164 } 165 166 ierr = MatCreate(comm, &tMat);CHKERRQ(ierr); 167 ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 168 ierr = MatSetType(tMat,MATAIJ);CHKERRQ(ierr); 169 ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 170 ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 171 ierr = PetscFree(d_nnz);CHKERRQ(ierr); 172 ierr = PetscFree(o_nnz);CHKERRQ(ierr); 173 174 for (ii = Istart_crs; ii < Iend_crs; ii++) { 175 PetscInt dest_row = ii/cr_bs; 176 ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 177 for (jj = 0; jj < ncols; jj++) { 178 PetscInt dest_col = idx[jj]/cr_bs; 179 PetscScalar v = 1.0; 180 ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr); 181 } 182 ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 183 } 184 ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 185 ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 186 187 if (llev++ == -1) { 188 PetscViewer viewer; char fname[32]; 189 ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr); 190 PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer); 191 ierr = MatView(tMat, viewer);CHKERRQ(ierr); 192 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 193 } 194 195 ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 196 197 ierr = MatDestroy(&tMat);CHKERRQ(ierr); 198 } /* create 'adj' */ 199 200 { /* partition: get newproc_idx */ 201 char prefix[256]; 202 const char *pcpre; 203 const PetscInt *is_idx; 204 MatPartitioning mpart; 205 IS proc_is; 206 PetscInt targetPE; 207 208 ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr); 209 ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr); 210 ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 211 ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 212 ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 213 ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 214 ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr); 215 ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr); 216 ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 217 218 /* collect IS info */ 219 ierr = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr); 220 ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr); 221 targetPE = 1; /* bring to "front" of machine */ 222 /*targetPE = size/new_size;*/ /* spread partitioning across machine */ 223 for (kk = jj = 0 ; kk < nloc_old ; kk++) { 224 for (ii = 0 ; ii < cr_bs ; ii++, jj++) { 225 newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */ 226 } 227 } 228 ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr); 229 ierr = ISDestroy(&proc_is);CHKERRQ(ierr); 230 } 231 ierr = MatDestroy(&adj);CHKERRQ(ierr); 232 233 ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr); 234 if (newproc_idx != 0) { 235 ierr = PetscFree(newproc_idx);CHKERRQ(ierr); 236 } 237 } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */ 238 239 PetscInt rfactor,targetPE; 240 /* find factor */ 241 if (new_size == 1) rfactor = size; /* easy */ 242 else { 243 PetscReal best_fact = 0.; 244 jj = -1; 245 for (kk = 1 ; kk <= size ; kk++) { 246 if (size%kk==0) { /* a candidate */ 247 PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size; 248 if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 249 if (fact > best_fact) { 250 best_fact = fact; jj = kk; 251 } 252 } 253 } 254 if (jj != -1) rfactor = jj; 255 else rfactor = 1; /* does this happen .. a prime */ 256 } 257 new_size = size/rfactor; 258 259 if (new_size==nactive) { 260 *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */ 261 ierr = PetscFree(counts);CHKERRQ(ierr); 262 if (pc_gamg->verbose>0) { 263 PetscPrintf(comm,"\t[%d]%s aggregate processors noop: new_size=%d, neq(loc)=%d\n",rank,__FUNCT__,new_size,ncrs_eq); 264 } 265 PetscFunctionReturn(0); 266 } 267 268 if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",rank,__FUNCT__,ncrs_eq); 269 targetPE = rank/rfactor; 270 ierr = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr); 271 } /* end simple 'is_eq_newproc' */ 272 273 /* 274 Create an index set from the is_eq_newproc index set to indicate the mapping TO 275 */ 276 ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr); 277 is_eq_num_prim = is_eq_num; 278 /* 279 Determine how many equations/vertices are assigned to each processor 280 */ 281 ierr = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr); 282 ncrs_eq_new = counts[rank]; 283 ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr); 284 ncrs_new = ncrs_eq_new/cr_bs; /* eqs */ 285 286 ierr = PetscFree(counts);CHKERRQ(ierr); 287 #if defined PETSC_GAMG_USE_LOG 288 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 289 #endif 290 /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxilary data into some complex abstracted thing */ 291 { 292 Vec src_crd, dest_crd; 293 const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols; 294 VecScatter vecscat; 295 PetscScalar *array; 296 IS isscat; 297 298 /* move data (for primal equations only) */ 299 /* Create a vector to contain the newly ordered element information */ 300 ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr); 301 ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr); 302 ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */ 303 /* 304 There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 305 a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 306 */ 307 ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr); 308 ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 309 for (ii=0,jj=0; ii<ncrs; ii++) { 310 PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 311 for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 312 } 313 ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 314 ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr); 315 ierr = PetscFree(tidx);CHKERRQ(ierr); 316 /* 317 Create a vector to contain the original vertex information for each element 318 */ 319 ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr); 320 for (jj=0; jj<ndata_cols; jj++) { 321 const PetscInt stride0=ncrs*pc_gamg->data_cell_rows; 322 for (ii=0; ii<ncrs; ii++) { 323 for (kk=0; kk<ndata_rows; kk++) { 324 PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 325 PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 326 ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr); 327 } 328 } 329 } 330 ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr); 331 ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr); 332 /* 333 Scatter the element vertex information (still in the original vertex ordering) 334 to the correct processor 335 */ 336 ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr); 337 ierr = ISDestroy(&isscat);CHKERRQ(ierr); 338 ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 339 ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 340 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 341 ierr = VecDestroy(&src_crd);CHKERRQ(ierr); 342 /* 343 Put the element vertex data into a new allocation of the gdata->ele 344 */ 345 ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 346 ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr); 347 348 pc_gamg->data_sz = node_data_sz*ncrs_new; 349 strideNew = ncrs_new*ndata_rows; 350 351 ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr); 352 for (jj=0; jj<ndata_cols; jj++) { 353 for (ii=0; ii<ncrs_new; ii++) { 354 for (kk=0; kk<ndata_rows; kk++) { 355 PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 356 pc_gamg->data[ix] = PetscRealPart(array[jx]); 357 } 358 } 359 } 360 ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr); 361 ierr = VecDestroy(&dest_crd);CHKERRQ(ierr); 362 } 363 /* move A and P (columns) with new layout */ 364 #if defined PETSC_GAMG_USE_LOG 365 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 366 #endif 367 368 /* 369 Invert for MatGetSubMatrix 370 */ 371 ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr); 372 ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */ 373 ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr); 374 if (is_eq_num != is_eq_num_prim) { 375 ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 376 } 377 if (Pcolumnperm) { 378 ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr); 379 *Pcolumnperm = new_eq_indices; 380 } 381 ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr); 382 #if defined PETSC_GAMG_USE_LOG 383 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 384 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 385 #endif 386 /* 'a_Amat_crs' output */ 387 { 388 Mat mat; 389 ierr = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr); 390 *a_Amat_crs = mat; 391 392 if (!PETSC_TRUE) { 393 PetscInt cbs, rbs; 394 ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr); 395 ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 396 ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr); 397 ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",rank,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr); 398 } 399 } 400 ierr = MatDestroy(&Cmat);CHKERRQ(ierr); 401 402 #if defined PETSC_GAMG_USE_LOG 403 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 404 #endif 405 /* prolongator */ 406 { 407 IS findices; 408 PetscInt Istart,Iend; 409 Mat Pnew; 410 ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr); 411 #if defined PETSC_GAMG_USE_LOG 412 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 413 #endif 414 ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr); 415 ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr); 416 ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr); 417 ierr = ISDestroy(&findices);CHKERRQ(ierr); 418 419 if (!PETSC_TRUE) { 420 PetscInt cbs, rbs; 421 ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr); 422 ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 423 ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr); 424 ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 425 } 426 #if defined PETSC_GAMG_USE_LOG 427 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 428 #endif 429 ierr = MatDestroy(a_P_inout);CHKERRQ(ierr); 430 431 /* output - repartitioned */ 432 *a_P_inout = Pnew; 433 } 434 ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr); 435 436 *a_nactive_proc = new_size; /* output */ 437 } 438 439 /* outout matrix data */ 440 if (!PETSC_TRUE) { 441 PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs; 442 if (llev==0) { 443 sprintf(fname,"Cmat_%d.m",llev++); 444 PetscViewerASCIIOpen(comm,fname,&viewer); 445 ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 446 ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr); 447 ierr = PetscViewerDestroy(&viewer); 448 } 449 sprintf(fname,"Cmat_%d.m",llev++); 450 PetscViewerASCIIOpen(comm,fname,&viewer); 451 ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 452 ierr = MatView(Cmat, viewer);CHKERRQ(ierr); 453 ierr = PetscViewerDestroy(&viewer); 454 } 455 PetscFunctionReturn(0); 456 } 457 458 /* -------------------------------------------------------------------------- */ 459 /* 460 PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 461 by setting data structures and options. 462 463 Input Parameter: 464 . pc - the preconditioner context 465 466 Application Interface Routine: PCSetUp() 467 468 Notes: 469 The interface routine PCSetUp() is not usually called directly by 470 the user, but instead is called by PCApply() if necessary. 471 */ 472 #undef __FUNCT__ 473 #define __FUNCT__ "PCSetUp_GAMG" 474 PetscErrorCode PCSetUp_GAMG(PC pc) 475 { 476 PetscErrorCode ierr; 477 PC_MG *mg = (PC_MG*)pc->data; 478 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 479 Mat Pmat = pc->pmat; 480 PetscInt fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS]; 481 MPI_Comm comm; 482 PetscMPIInt rank,size,nactivepe; 483 Mat Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS]; 484 PetscReal emaxs[GAMG_MAXLEVELS]; 485 IS *ASMLocalIDsArr[GAMG_MAXLEVELS]; 486 PetscLogDouble nnz0=0.,nnztot=0.; 487 MatInfo info; 488 PetscBool redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol); 489 490 PetscFunctionBegin; 491 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 492 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 493 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 494 495 if (pc_gamg->verbose>2) PetscPrintf(comm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",rank,__FUNCT__,pc_gamg->setup_count,pc->setupcalled); 496 497 if (pc_gamg->setup_count++ > 0) { 498 if (redo_mesh_setup) { 499 /* reset everything */ 500 ierr = PCReset_MG(pc);CHKERRQ(ierr); 501 pc->setupcalled = 0; 502 } else { 503 PC_MG_Levels **mglevels = mg->levels; 504 /* just do Galerkin grids */ 505 Mat B,dA,dB; 506 507 if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet"); 508 if (pc_gamg->Nlevels > 1) { 509 /* currently only handle case where mat and pmat are the same on coarser levels */ 510 ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 511 /* (re)set to get dirty flag */ 512 ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr); 513 514 for (level=pc_gamg->Nlevels-2; level>=0; level--) { 515 /* the first time through the matrix structure has changed from repartitioning */ 516 if (pc_gamg->setup_count==2) { 517 ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 518 ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); 519 520 mglevels[level]->A = B; 521 } else { 522 ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr); 523 ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 524 } 525 ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr); 526 dB = B; 527 } 528 } 529 530 ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 531 532 PetscFunctionReturn(0); 533 } 534 } 535 536 if (!pc_gamg->data) { 537 if (pc_gamg->orig_data) { 538 ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 539 ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr); 540 541 pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 542 pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 543 pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 544 545 ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr); 546 for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 547 } else { 548 if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 549 ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr); 550 } 551 } 552 553 /* cache original data for reuse */ 554 if (!pc_gamg->orig_data && redo_mesh_setup) { 555 ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr); 556 for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 557 pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 558 pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 559 } 560 561 /* get basic dims */ 562 ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 563 564 ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr); 565 if (pc_gamg->verbose) { 566 PetscInt NN = M; 567 if (pc_gamg->verbose==1) { 568 ierr = MatGetInfo(Pmat,MAT_LOCAL,&info);CHKERRQ(ierr); 569 ierr = MatGetLocalSize(Pmat, &NN, &qq);CHKERRQ(ierr); 570 if (!NN) NN=1; 571 } else { 572 ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); 573 } 574 nnz0 = info.nz_used; 575 nnztot = info.nz_used; 576 ierr = PetscPrintf(comm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 577 rank,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols, 578 (int)(nnz0/(PetscReal)NN),size);CHKERRQ(ierr); 579 } 580 581 /* Get A_i and R_i */ 582 for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */ 583 level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); 584 level++) { 585 pc_gamg->firstCoarsen = (level ? PETSC_FALSE : PETSC_TRUE); 586 level1 = level + 1; 587 #if defined PETSC_GAMG_USE_LOG 588 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 589 #if (defined GAMG_STAGES) 590 ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr); 591 #endif 592 #endif 593 { /* construct prolongator */ 594 Mat Gmat; 595 PetscCoarsenData *agg_lists; 596 Mat Prol11; 597 598 ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr); 599 ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr); 600 ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr); 601 602 /* could have failed to create new level */ 603 if (Prol11) { 604 /* get new block size of coarse matrices */ 605 ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr); 606 607 if (pc_gamg->ops->optprol) { 608 /* smooth */ 609 ierr = pc_gamg->ops->optprol(pc, Aarr[level], &Prol11);CHKERRQ(ierr); 610 } 611 612 Parr[level1] = Prol11; 613 } else Parr[level1] = NULL; 614 615 if (pc_gamg->use_aggs_in_gasm) { 616 PetscInt bs; 617 ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr); 618 ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr); 619 } 620 621 ierr = MatDestroy(&Gmat);CHKERRQ(ierr); 622 ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr); 623 } /* construct prolongator scope */ 624 #if defined PETSC_GAMG_USE_LOG 625 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 626 #endif 627 /* cache eigen estimate */ 628 if (pc_gamg->emax_id != -1) { 629 PetscBool flag; 630 ierr = PetscObjectComposedDataGetReal((PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag);CHKERRQ(ierr); 631 if (!flag) emaxs[level] = -1.; 632 } else emaxs[level] = -1.; 633 if (level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 634 if (!Parr[level1]) { 635 if (pc_gamg->verbose) { 636 ierr = PetscPrintf(comm,"\t[%d]%s stop gridding, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr); 637 } 638 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 639 ierr = PetscLogStagePop();CHKERRQ(ierr); 640 #endif 641 break; 642 } 643 #if defined PETSC_GAMG_USE_LOG 644 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 645 #endif 646 647 ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, 648 &Parr[level1], &Aarr[level1], &nactivepe, NULL);CHKERRQ(ierr); 649 650 #if defined PETSC_GAMG_USE_LOG 651 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 652 #endif 653 ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr); 654 655 if (pc_gamg->verbose > 0) { 656 PetscInt NN = M; 657 if (pc_gamg->verbose==1) { 658 ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info);CHKERRQ(ierr); 659 ierr = MatGetLocalSize(Aarr[level1], &NN, &qq);CHKERRQ(ierr); 660 if (!NN) NN=1; 661 } else { 662 ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr); 663 } 664 665 nnztot += info.nz_used; 666 ierr = PetscPrintf(comm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 667 rank,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols, 668 (int)(info.nz_used/(PetscReal)NN), nactivepe);CHKERRQ(ierr); 669 } 670 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 671 ierr = PetscLogStagePop();CHKERRQ(ierr); 672 #endif 673 /* stop if one node or one proc -- could pull back for singular problems */ 674 if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M < 2) ) { 675 if (pc_gamg->verbose) { 676 ierr = PetscPrintf(comm,"\t[%d]%s HARD stop of coarsening ?????????, level %d\n",rank,__FUNCT__,level);CHKERRQ(ierr); 677 } 678 level++; 679 break; 680 } 681 } /* levels */ 682 pc_gamg->firstCoarsen = PETSC_FALSE; 683 684 if (pc_gamg->data) { 685 ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 686 pc_gamg->data = NULL; 687 } 688 689 if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0); 690 pc_gamg->Nlevels = level + 1; 691 fine_level = level; 692 ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 693 694 /* simple setup */ 695 if (!PETSC_TRUE) { 696 PC_MG_Levels **mglevels = mg->levels; 697 for (lidx=0,level=pc_gamg->Nlevels-1; 698 lidx<fine_level; 699 lidx++, level--) { 700 ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr); 701 ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level]);CHKERRQ(ierr); 702 ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 703 ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 704 } 705 ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0]);CHKERRQ(ierr); 706 707 ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 708 } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 709 /* set default smoothers & set operators */ 710 for (lidx = 1, level = pc_gamg->Nlevels-2; 711 lidx <= fine_level; 712 lidx++, level--) { 713 KSP smoother; 714 PC subpc; 715 716 ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 717 ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 718 719 ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 720 /* set ops */ 721 ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr); 722 ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr); 723 724 /* set defaults */ 725 ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr); 726 727 /* set blocks for GASM smoother that uses the 'aggregates' */ 728 if (pc_gamg->use_aggs_in_gasm) { 729 PetscInt sz; 730 IS *is; 731 732 sz = nASMBlocksArr[level]; 733 is = ASMLocalIDsArr[level]; 734 ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr); 735 ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr); 736 if (sz==0) { 737 IS is; 738 PetscInt my0,kk; 739 ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr); 740 ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr); 741 ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr); 742 ierr = ISDestroy(&is);CHKERRQ(ierr); 743 } else { 744 PetscInt kk; 745 ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr); 746 for (kk=0; kk<sz; kk++) { 747 ierr = ISDestroy(&is[kk]);CHKERRQ(ierr); 748 } 749 ierr = PetscFree(is);CHKERRQ(ierr); 750 } 751 ASMLocalIDsArr[level] = NULL; 752 nASMBlocksArr[level] = 0; 753 ierr = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr); 754 } else { 755 ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr); 756 } 757 } 758 { 759 /* coarse grid */ 760 KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 761 Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 762 ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 763 ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr); 764 ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 765 ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 766 ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr); 767 ierr = PCSetUp(subpc);CHKERRQ(ierr); 768 ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 769 if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii); 770 ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 771 ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr); 772 ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 773 ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 774 /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in 775 * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to 776 * view every subdomain as though they were different. */ 777 ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE; 778 } 779 780 /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 781 ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 782 ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); 783 ierr = PetscOptionsEnd();CHKERRQ(ierr); 784 if (mg->galerkin != 2) SETERRQ(comm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used."); 785 786 /* create cheby smoothers */ 787 for (lidx = 1, level = pc_gamg->Nlevels-2; 788 lidx <= fine_level; 789 lidx++, level--) { 790 KSP smoother; 791 PetscBool flag,flag2; 792 PC subpc; 793 794 ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 795 ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 796 797 /* do my own cheby */ 798 ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr); 799 if (flag) { 800 PetscReal emax, emin; 801 ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr); 802 ierr = PetscObjectTypeCompare((PetscObject)subpc, PCSOR, &flag2);CHKERRQ(ierr); 803 if ((flag||flag2) && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC but lets acccept SOR because it is close and safe (always lower) */ 804 else { /* eigen estimate 'emax' -- this is done in cheby */ 805 KSP eksp; 806 Mat Lmat = Aarr[level]; 807 Vec bb, xx; 808 809 ierr = MatCreateVecs(Lmat, &bb, 0);CHKERRQ(ierr); 810 ierr = MatCreateVecs(Lmat, &xx, 0);CHKERRQ(ierr); 811 { 812 PetscRandom rctx; 813 ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr); 814 ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 815 ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 816 ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 817 } 818 819 /* zeroing out BC rows -- needed for crazy matrices */ 820 { 821 PetscInt Istart,Iend,ncols,jj,Ii; 822 PetscScalar zero = 0.0; 823 ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr); 824 for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) { 825 ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr); 826 if (ncols <= 1) { 827 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr); 828 } 829 ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr); 830 } 831 ierr = VecAssemblyBegin(bb);CHKERRQ(ierr); 832 ierr = VecAssemblyEnd(bb);CHKERRQ(ierr); 833 } 834 835 ierr = KSPCreate(comm, &eksp);CHKERRQ(ierr); 836 ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr); 837 ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 838 ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 839 ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 840 ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 841 842 ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 843 ierr = KSPSetOperators(eksp, Lmat, Lmat);CHKERRQ(ierr); 844 ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 845 846 /* set PC type to be same as smoother */ 847 ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr); 848 849 /* solve - keep stuff out of logging */ 850 ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 851 ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 852 ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 853 ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 854 ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 855 856 ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 857 858 ierr = VecDestroy(&xx);CHKERRQ(ierr); 859 ierr = VecDestroy(&bb);CHKERRQ(ierr); 860 ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 861 862 if (pc_gamg->verbose > 0) { 863 PetscInt N1, tt; 864 ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr); 865 PetscPrintf(comm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1); 866 } 867 } 868 { 869 PetscInt N1, N0; 870 ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr); 871 ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr); 872 /* heuristic - is this crap? */ 873 /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */ 874 emin = emax * pc_gamg->eigtarget[0]; 875 emax *= pc_gamg->eigtarget[1]; 876 } 877 ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); 878 } /* setup checby flag */ 879 } /* non-coarse levels */ 880 881 /* clean up */ 882 for (level=1; level<pc_gamg->Nlevels; level++) { 883 ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 884 ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 885 } 886 887 ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 888 } else { 889 KSP smoother; 890 if (pc_gamg->verbose) PetscPrintf(comm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",rank,__FUNCT__); 891 ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 892 ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr); 893 ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 894 ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 895 } 896 PetscFunctionReturn(0); 897 } 898 899 /* ------------------------------------------------------------------------- */ 900 /* 901 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 902 that was created with PCCreate_GAMG(). 903 904 Input Parameter: 905 . pc - the preconditioner context 906 907 Application Interface Routine: PCDestroy() 908 */ 909 #undef __FUNCT__ 910 #define __FUNCT__ "PCDestroy_GAMG" 911 PetscErrorCode PCDestroy_GAMG(PC pc) 912 { 913 PetscErrorCode ierr; 914 PC_MG *mg = (PC_MG*)pc->data; 915 PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 916 917 PetscFunctionBegin; 918 ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 919 if (pc_gamg->ops->destroy) { 920 ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 921 } 922 ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr); 923 ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 924 ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 925 ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 926 PetscFunctionReturn(0); 927 } 928 929 930 #undef __FUNCT__ 931 #define __FUNCT__ "PCGAMGSetProcEqLim" 932 /*@ 933 PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via processor reduction. 934 935 Logically Collective on PC 936 937 Input Parameters: 938 + pc - the preconditioner context 939 - n - the number of equations 940 941 942 Options Database Key: 943 . -pc_gamg_process_eq_limit <limit> 944 945 Level: intermediate 946 947 Concepts: Unstructured multrigrid preconditioner 948 949 .seealso: () 950 @*/ 951 PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 952 { 953 PetscErrorCode ierr; 954 955 PetscFunctionBegin; 956 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 957 ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 958 PetscFunctionReturn(0); 959 } 960 961 #undef __FUNCT__ 962 #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 963 static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 964 { 965 PC_MG *mg = (PC_MG*)pc->data; 966 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 967 968 PetscFunctionBegin; 969 if (n>0) pc_gamg->min_eq_proc = n; 970 PetscFunctionReturn(0); 971 } 972 973 #undef __FUNCT__ 974 #define __FUNCT__ "PCGAMGSetCoarseEqLim" 975 /*@ 976 PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 977 978 Collective on PC 979 980 Input Parameters: 981 + pc - the preconditioner context 982 - n - maximum number of equations to aim for 983 984 Options Database Key: 985 . -pc_gamg_coarse_eq_limit <limit> 986 987 Level: intermediate 988 989 Concepts: Unstructured multrigrid preconditioner 990 991 @*/ 992 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 993 { 994 PetscErrorCode ierr; 995 996 PetscFunctionBegin; 997 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 998 ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 999 PetscFunctionReturn(0); 1000 } 1001 1002 #undef __FUNCT__ 1003 #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 1004 static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 1005 { 1006 PC_MG *mg = (PC_MG*)pc->data; 1007 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1008 1009 PetscFunctionBegin; 1010 if (n>0) pc_gamg->coarse_eq_limit = n; 1011 PetscFunctionReturn(0); 1012 } 1013 1014 #undef __FUNCT__ 1015 #define __FUNCT__ "PCGAMGSetRepartitioning" 1016 /*@ 1017 PCGAMGSetRepartitioning - Repartition the coarse grids 1018 1019 Collective on PC 1020 1021 Input Parameters: 1022 + pc - the preconditioner context 1023 - n - PETSC_TRUE or PETSC_FALSE 1024 1025 Options Database Key: 1026 . -pc_gamg_repartition <true,false> 1027 1028 Level: intermediate 1029 1030 Concepts: Unstructured multrigrid preconditioner 1031 1032 .seealso: () 1033 @*/ 1034 PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 1035 { 1036 PetscErrorCode ierr; 1037 1038 PetscFunctionBegin; 1039 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1040 ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1041 PetscFunctionReturn(0); 1042 } 1043 1044 #undef __FUNCT__ 1045 #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 1046 static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 1047 { 1048 PC_MG *mg = (PC_MG*)pc->data; 1049 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1050 1051 PetscFunctionBegin; 1052 pc_gamg->repart = n; 1053 PetscFunctionReturn(0); 1054 } 1055 1056 #undef __FUNCT__ 1057 #define __FUNCT__ "PCGAMGSetReuseInterpolation" 1058 /*@ 1059 PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding 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> 1069 1070 Level: intermediate 1071 1072 Concepts: Unstructured multrigrid preconditioner 1073 1074 .seealso: () 1075 @*/ 1076 PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 1077 { 1078 PetscErrorCode ierr; 1079 1080 PetscFunctionBegin; 1081 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1082 ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1083 PetscFunctionReturn(0); 1084 } 1085 1086 #undef __FUNCT__ 1087 #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG" 1088 static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 1089 { 1090 PC_MG *mg = (PC_MG*)pc->data; 1091 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1092 1093 PetscFunctionBegin; 1094 pc_gamg->reuse_prol = n; 1095 PetscFunctionReturn(0); 1096 } 1097 1098 #undef __FUNCT__ 1099 #define __FUNCT__ "PCGAMGSetUseASMAggs" 1100 /*@ 1101 PCGAMGSetUseASMAggs - 1102 1103 Collective on PC 1104 1105 Input Parameters: 1106 . pc - the preconditioner context 1107 1108 1109 Options Database Key: 1110 . -pc_gamg_use_agg_gasm 1111 1112 Level: intermediate 1113 1114 Concepts: Unstructured multrigrid preconditioner 1115 1116 .seealso: () 1117 @*/ 1118 PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n) 1119 { 1120 PetscErrorCode ierr; 1121 1122 PetscFunctionBegin; 1123 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1124 ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1125 PetscFunctionReturn(0); 1126 } 1127 1128 #undef __FUNCT__ 1129 #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG" 1130 static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n) 1131 { 1132 PC_MG *mg = (PC_MG*)pc->data; 1133 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1134 1135 PetscFunctionBegin; 1136 pc_gamg->use_aggs_in_gasm = n; 1137 PetscFunctionReturn(0); 1138 } 1139 1140 #undef __FUNCT__ 1141 #define __FUNCT__ "PCGAMGSetNlevels" 1142 /*@ 1143 PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use 1144 1145 Not collective on PC 1146 1147 Input Parameters: 1148 + pc - the preconditioner 1149 - n - the maximum number of levels to use 1150 1151 Options Database Key: 1152 . -pc_mg_levels 1153 1154 Level: intermediate 1155 1156 Concepts: Unstructured multrigrid preconditioner 1157 1158 .seealso: () 1159 @*/ 1160 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 1161 { 1162 PetscErrorCode ierr; 1163 1164 PetscFunctionBegin; 1165 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1166 ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1167 PetscFunctionReturn(0); 1168 } 1169 1170 #undef __FUNCT__ 1171 #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 1172 static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 1173 { 1174 PC_MG *mg = (PC_MG*)pc->data; 1175 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1176 1177 PetscFunctionBegin; 1178 pc_gamg->Nlevels = n; 1179 PetscFunctionReturn(0); 1180 } 1181 1182 #undef __FUNCT__ 1183 #define __FUNCT__ "PCGAMGSetThreshold" 1184 /*@ 1185 PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 1186 1187 Not collective on PC 1188 1189 Input Parameters: 1190 + pc - the preconditioner context 1191 - threshold - the threshold value, 0.0 is the lowest possible 1192 1193 Options Database Key: 1194 . -pc_gamg_threshold <threshold> 1195 1196 Level: intermediate 1197 1198 Concepts: Unstructured multrigrid preconditioner 1199 1200 .seealso: () 1201 @*/ 1202 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 1203 { 1204 PetscErrorCode ierr; 1205 1206 PetscFunctionBegin; 1207 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1208 ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 1209 PetscFunctionReturn(0); 1210 } 1211 1212 #undef __FUNCT__ 1213 #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 1214 static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 1215 { 1216 PC_MG *mg = (PC_MG*)pc->data; 1217 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1218 1219 PetscFunctionBegin; 1220 pc_gamg->threshold = n; 1221 PetscFunctionReturn(0); 1222 } 1223 1224 #undef __FUNCT__ 1225 #define __FUNCT__ "PCGAMGSetType" 1226 /*@ 1227 PCGAMGSetType - Set solution method 1228 1229 Collective on PC 1230 1231 Input Parameters: 1232 + pc - the preconditioner context 1233 - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1234 1235 Options Database Key: 1236 . -pc_gamg_type <agg,geo,classical> 1237 1238 Level: intermediate 1239 1240 Concepts: Unstructured multrigrid preconditioner 1241 1242 .seealso: PCGAMGGetType(), PCGAMG 1243 @*/ 1244 PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1245 { 1246 PetscErrorCode ierr; 1247 1248 PetscFunctionBegin; 1249 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1250 ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 1251 PetscFunctionReturn(0); 1252 } 1253 1254 #undef __FUNCT__ 1255 #define __FUNCT__ "PCGAMGGetType" 1256 /*@ 1257 PCGAMGGetType - Get solution method 1258 1259 Collective on PC 1260 1261 Input Parameter: 1262 . pc - the preconditioner context 1263 1264 Output Parameter: 1265 . type - the type of algorithm used 1266 1267 Level: intermediate 1268 1269 Concepts: Unstructured multrigrid preconditioner 1270 1271 .seealso: PCGAMGSetType() 1272 @*/ 1273 PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1274 { 1275 PetscErrorCode ierr; 1276 1277 PetscFunctionBegin; 1278 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1279 ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 1280 PetscFunctionReturn(0); 1281 } 1282 1283 #undef __FUNCT__ 1284 #define __FUNCT__ "PCGAMGGetType_GAMG" 1285 static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1286 { 1287 PC_MG *mg = (PC_MG*)pc->data; 1288 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1289 1290 PetscFunctionBegin; 1291 *type = pc_gamg->type; 1292 PetscFunctionReturn(0); 1293 } 1294 1295 #undef __FUNCT__ 1296 #define __FUNCT__ "PCGAMGSetType_GAMG" 1297 static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1298 { 1299 PetscErrorCode ierr,(*r)(PC); 1300 PC_MG *mg = (PC_MG*)pc->data; 1301 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1302 1303 PetscFunctionBegin; 1304 pc_gamg->type = type; 1305 ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr); 1306 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 1307 if (pc_gamg->ops->destroy) { 1308 /* there was something here - kill it */ 1309 ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 1310 ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr); 1311 pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1312 /* cleaning up common data in pc_gamg - this should disapear someday */ 1313 pc_gamg->data_cell_cols = 0; 1314 pc_gamg->data_cell_rows = 0; 1315 pc_gamg->orig_data_cell_cols = 0; 1316 pc_gamg->orig_data_cell_rows = 0; 1317 if (pc_gamg->data_sz) { 1318 ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 1319 pc_gamg->data_sz = 0; 1320 } else if (pc_gamg->data) { 1321 ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); /* can this happen ? */ 1322 } 1323 } 1324 ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 1325 ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr); 1326 ierr = (*r)(pc);CHKERRQ(ierr); 1327 PetscFunctionReturn(0); 1328 } 1329 1330 #undef __FUNCT__ 1331 #define __FUNCT__ "PCSetFromOptions_GAMG" 1332 PetscErrorCode PCSetFromOptions_GAMG(PetscOptions *PetscOptionsObject,PC pc) 1333 { 1334 PetscErrorCode ierr; 1335 PC_MG *mg = (PC_MG*)pc->data; 1336 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1337 PetscBool flag; 1338 PetscInt two = 2; 1339 MPI_Comm comm; 1340 1341 PetscFunctionBegin; 1342 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1343 ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr); 1344 { 1345 /* -pc_gamg_type */ 1346 { 1347 char tname[256]; 1348 ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr); 1349 if (flag) { 1350 ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr); 1351 } 1352 } 1353 /* -pc_gamg_verbose */ 1354 ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG","none", pc_gamg->verbose,&pc_gamg->verbose, NULL);CHKERRQ(ierr); 1355 /* -pc_gamg_repartition */ 1356 ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr); 1357 /* -pc_gamg_reuse_interpolation */ 1358 ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr); 1359 /* -pc_gamg_use_agg_gasm */ 1360 ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm","Use aggregation agragates for GASM smoother","PCGAMGUseASMAggs",pc_gamg->use_aggs_in_gasm,&pc_gamg->use_aggs_in_gasm,NULL);CHKERRQ(ierr); 1361 /* -pc_gamg_process_eq_limit */ 1362 ierr = 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);CHKERRQ(ierr); 1363 /* -pc_gamg_coarse_eq_limit */ 1364 ierr = 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);CHKERRQ(ierr); 1365 /* -pc_gamg_threshold */ 1366 ierr = PetscOptionsReal("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&pc_gamg->threshold,&flag);CHKERRQ(ierr); 1367 if (flag && pc_gamg->verbose) { 1368 ierr = PetscPrintf(comm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr); 1369 } 1370 /* -pc_gamg_eigtarget */ 1371 ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr); 1372 ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 1373 1374 /* set options for subtype */ 1375 if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);} 1376 } 1377 ierr = PetscOptionsTail();CHKERRQ(ierr); 1378 PetscFunctionReturn(0); 1379 } 1380 1381 /* -------------------------------------------------------------------------- */ 1382 /*MC 1383 PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 1384 1385 Options Database Keys: 1386 Multigrid options(inherited) 1387 + -pc_mg_cycles <v>: v or w (PCMGSetCycleType()) 1388 . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1389 . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 1390 - -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade 1391 1392 1393 Notes: In order to obtain good performance for PCGAMG for vector valued problems you must 1394 $ Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point 1395 $ Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator 1396 $ See the Users Manual Chapter 4 for more details 1397 1398 Level: intermediate 1399 1400 Concepts: algebraic multigrid 1401 1402 .seealso: PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(), 1403 PCGAMGSetCoarseEqLim(), PCGAMGSetRepartitioning(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseASMAggs(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType() 1404 M*/ 1405 1406 #undef __FUNCT__ 1407 #define __FUNCT__ "PCCreate_GAMG" 1408 PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 1409 { 1410 PetscErrorCode ierr; 1411 PC_GAMG *pc_gamg; 1412 PC_MG *mg; 1413 #if defined PETSC_GAMG_USE_LOG 1414 static long count = 0; 1415 #endif 1416 1417 PetscFunctionBegin; 1418 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 1419 ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 1420 ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr); 1421 1422 /* create a supporting struct and attach it to pc */ 1423 ierr = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr); 1424 mg = (PC_MG*)pc->data; 1425 mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 1426 mg->innerctx = pc_gamg; 1427 1428 ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr); 1429 1430 pc_gamg->setup_count = 0; 1431 /* these should be in subctx but repartitioning needs simple arrays */ 1432 pc_gamg->data_sz = 0; 1433 pc_gamg->data = 0; 1434 1435 /* register AMG type */ 1436 ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 1437 1438 /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 1439 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 1440 pc->ops->setup = PCSetUp_GAMG; 1441 pc->ops->reset = PCReset_GAMG; 1442 pc->ops->destroy = PCDestroy_GAMG; 1443 1444 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr); 1445 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr); 1446 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr); 1447 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr); 1448 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr); 1449 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr); 1450 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr); 1451 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr); 1452 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr); 1453 pc_gamg->repart = PETSC_FALSE; 1454 pc_gamg->reuse_prol = PETSC_FALSE; 1455 pc_gamg->use_aggs_in_gasm = PETSC_FALSE; 1456 pc_gamg->min_eq_proc = 50; 1457 pc_gamg->coarse_eq_limit = 50; 1458 pc_gamg->threshold = 0.; 1459 pc_gamg->Nlevels = GAMG_MAXLEVELS; 1460 pc_gamg->verbose = 0; 1461 pc_gamg->emax_id = -1; 1462 pc_gamg->firstCoarsen = PETSC_FALSE; 1463 pc_gamg->eigtarget[0] = 0.05; 1464 pc_gamg->eigtarget[1] = 1.05; 1465 pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 1466 1467 /* private events */ 1468 #if defined PETSC_GAMG_USE_LOG 1469 if (count++ == 0) { 1470 ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr); 1471 ierr = PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr); 1472 /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 1473 /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 1474 /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 1475 ierr = PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr); 1476 ierr = PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr); 1477 ierr = PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr); 1478 ierr = PetscLogEventRegister(" search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr); 1479 ierr = PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr); 1480 ierr = PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr); 1481 ierr = PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr); 1482 ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr); 1483 ierr = PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr); 1484 ierr = PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr); 1485 ierr = PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr); 1486 ierr = PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr); 1487 1488 /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 1489 /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 1490 /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 1491 /* create timer stages */ 1492 #if defined GAMG_STAGES 1493 { 1494 char str[32]; 1495 PetscInt lidx; 1496 sprintf(str,"MG Level %d (finest)",0); 1497 ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr); 1498 for (lidx=1; lidx<9; lidx++) { 1499 sprintf(str,"MG Level %d",lidx); 1500 ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr); 1501 } 1502 } 1503 #endif 1504 } 1505 #endif 1506 /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1507 ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr); 1508 PetscFunctionReturn(0); 1509 } 1510 1511 #undef __FUNCT__ 1512 #define __FUNCT__ "PCGAMGInitializePackage" 1513 /*@C 1514 PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 1515 from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG() 1516 when using static libraries. 1517 1518 Level: developer 1519 1520 .keywords: PC, PCGAMG, initialize, package 1521 .seealso: PetscInitialize() 1522 @*/ 1523 PetscErrorCode PCGAMGInitializePackage(void) 1524 { 1525 PetscErrorCode ierr; 1526 1527 PetscFunctionBegin; 1528 if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 1529 PCGAMGPackageInitialized = PETSC_TRUE; 1530 ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr); 1531 ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr); 1532 ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr); 1533 ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr); 1534 1535 /* general events */ 1536 #if defined PETSC_USE_LOG 1537 ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr); 1538 ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr); 1539 ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr); 1540 ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr); 1541 ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr); 1542 ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr); 1543 ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr); 1544 #endif 1545 1546 PetscFunctionReturn(0); 1547 } 1548 1549 #undef __FUNCT__ 1550 #define __FUNCT__ "PCGAMGFinalizePackage" 1551 /*@C 1552 PCGAMGFinalizePackage - This function destroys everything in the PCGAMG package. It is 1553 called from PetscFinalize(). 1554 1555 Level: developer 1556 1557 .keywords: Petsc, destroy, package 1558 .seealso: PetscFinalize() 1559 @*/ 1560 PetscErrorCode PCGAMGFinalizePackage(void) 1561 { 1562 PetscErrorCode ierr; 1563 1564 PetscFunctionBegin; 1565 PCGAMGPackageInitialized = PETSC_FALSE; 1566 ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr); 1567 PetscFunctionReturn(0); 1568 } 1569 1570 #undef __FUNCT__ 1571 #define __FUNCT__ "PCGAMGRegister" 1572 /*@C 1573 PCGAMGRegister - Register a PCGAMG implementation. 1574 1575 Input Parameters: 1576 + type - string that will be used as the name of the GAMG type. 1577 - create - function for creating the gamg context. 1578 1579 Level: advanced 1580 1581 .seealso: PCGAMGGetContext(), PCGAMGSetContext() 1582 @*/ 1583 PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1584 { 1585 PetscErrorCode ierr; 1586 1587 PetscFunctionBegin; 1588 ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 1589 ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr); 1590 PetscFunctionReturn(0); 1591 } 1592 1593