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