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