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