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 8 #if defined PETSC_GAMG_USE_LOG 9 PetscLogEvent petsc_gamg_setup_events[NUM_SET]; 10 #endif 11 12 #if defined PETSC_USE_LOG 13 PetscLogEvent PC_GAMGGgraph_AGG; 14 PetscLogEvent PC_GAMGGgraph_GEO; 15 PetscLogEvent PC_GAMGCoarsen_AGG; 16 PetscLogEvent PC_GAMGCoarsen_GEO; 17 PetscLogEvent PC_GAMGProlongator_AGG; 18 PetscLogEvent PC_GAMGProlongator_GEO; 19 PetscLogEvent PC_GAMGOptprol_AGG; 20 PetscLogEvent PC_GAMGKKTProl_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 32 /* ----------------------------------------------------------------------------- */ 33 #undef __FUNCT__ 34 #define __FUNCT__ "PCReset_GAMG" 35 PetscErrorCode PCReset_GAMG(PC pc) 36 { 37 PetscErrorCode ierr; 38 PC_MG *mg = (PC_MG*)pc->data; 39 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 40 41 PetscFunctionBegin; 42 if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */ 43 PetscPrintf(((PetscObject)pc)->comm,"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__); 44 ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 45 } 46 pc_gamg->data = NULL; pc_gamg->data_sz = 0; 47 48 if (pc_gamg->orig_data) { 49 ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr); 50 } 51 PetscFunctionReturn(0); 52 } 53 54 /* private 2x2 Mat Nest for Stokes */ 55 typedef struct { 56 Mat A11,A21,A12,Amat; 57 IS prim_is,constr_is; 58 } GAMGKKTMat; 59 60 #undef __FUNCT__ 61 #define __FUNCT__ "GAMGKKTMatCreate" 62 static PetscErrorCode GAMGKKTMatCreate(Mat A, PetscBool iskkt, GAMGKKTMat *out) 63 { 64 PetscFunctionBegin; 65 out->Amat = A; 66 if (iskkt) { 67 PetscErrorCode ierr; 68 IS is_constraint, is_prime; 69 PetscInt nmin,nmax; 70 71 ierr = MatGetOwnershipRange(A, &nmin, &nmax);CHKERRQ(ierr); 72 ierr = MatFindZeroDiagonals(A, &is_constraint);CHKERRQ(ierr); 73 ierr = ISComplement(is_constraint, nmin, nmax, &is_prime);CHKERRQ(ierr); 74 75 out->prim_is = is_prime; 76 out->constr_is = is_constraint; 77 78 ierr = MatGetSubMatrix(A, is_prime, is_prime, MAT_INITIAL_MATRIX, &out->A11);CHKERRQ(ierr); 79 ierr = MatGetSubMatrix(A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12);CHKERRQ(ierr); 80 ierr = MatGetSubMatrix(A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21);CHKERRQ(ierr); 81 } else { 82 out->A11 = A; 83 out->A21 = NULL; 84 out->A12 = NULL; 85 out->prim_is = NULL; 86 out->constr_is = NULL; 87 } 88 PetscFunctionReturn(0); 89 } 90 91 #undef __FUNCT__ 92 #define __FUNCT__ "GAMGKKTMatDestroy" 93 static PetscErrorCode GAMGKKTMatDestroy(GAMGKKTMat *mat) 94 { 95 PetscErrorCode ierr; 96 97 PetscFunctionBegin; 98 if (mat->A11 && mat->A11 != mat->Amat) { 99 ierr = MatDestroy(&mat->A11);CHKERRQ(ierr); 100 } 101 ierr = MatDestroy(&mat->A21);CHKERRQ(ierr); 102 ierr = MatDestroy(&mat->A12);CHKERRQ(ierr); 103 104 ierr = ISDestroy(&mat->prim_is);CHKERRQ(ierr); 105 ierr = ISDestroy(&mat->constr_is);CHKERRQ(ierr); 106 PetscFunctionReturn(0); 107 } 108 109 /* -------------------------------------------------------------------------- */ 110 /* 111 createLevel: create coarse op with RAP. repartition and/or reduce number 112 of active processors. 113 114 Input Parameter: 115 . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 116 'pc_gamg->data_sz' are changed via repartitioning/reduction. 117 . Amat_fine - matrix on this fine (k) level 118 . cr_bs - coarse block size 119 . isLast - 120 . stokes - 121 In/Output Parameter: 122 . a_P_inout - prolongation operator to the next level (k-->k-1) 123 . a_nactive_proc - number of active procs 124 Output Parameter: 125 . a_Amat_crs - coarse matrix that is created (k-1) 126 */ 127 128 #undef __FUNCT__ 129 #define __FUNCT__ "createLevel" 130 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) 131 { 132 PetscErrorCode ierr; 133 PC_MG *mg = (PC_MG*)pc->data; 134 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 135 const PetscBool repart = pc_gamg->repart; 136 const PetscInt min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit; 137 Mat Cmat,Pold=*a_P_inout; 138 MPI_Comm wcomm = ((PetscObject)Amat_fine)->comm; 139 PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc; 140 PetscInt ncrs_eq,ncrs_prim,f_bs; 141 142 PetscFunctionBegin; 143 ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr); 144 ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr); 145 ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr); 146 /* RAP */ 147 ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr); 148 149 /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/ 150 ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 151 if (pc_gamg->data_sz % (pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data_sz %D not divisible by (pc_gamg->data_cell_cols %D *pc_gamg->data_cell_rows %D)",pc_gamg->data_sz,pc_gamg->data_cell_cols,pc_gamg->data_cell_rows); 152 ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr); 153 154 /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 155 { 156 PetscInt ncrs_eq_glob; 157 ierr = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr); 158 new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 159 if (new_size == 0 || ncrs_eq_glob < coarse_max) new_size = 1; 160 else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 161 if (isLast) new_size = 1; 162 } 163 164 if (!repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 165 else { 166 const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols; 167 PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old; 168 IS is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices; 169 VecScatter vecscat; 170 PetscScalar *array; 171 Vec src_crd, dest_crd; 172 173 nloc_old = ncrs_eq/cr_bs; 174 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); 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 558 if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet"); 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); 850 if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii); 851 ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 852 ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr); 853 ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 854 } 855 856 /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 857 ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 858 ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); 859 ierr = PetscOptionsEnd();CHKERRQ(ierr); 860 if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used."); 861 862 /* create cheby smoothers */ 863 for (lidx = 1, level = pc_gamg->Nlevels-2; 864 lidx <= fine_level; 865 lidx++, level--) { 866 KSP smoother; 867 PetscBool flag; 868 PC subpc; 869 870 ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 871 ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 872 873 /* create field split PC, get subsmoother */ 874 if (stokes) { 875 KSP *ksps; 876 PetscInt nn; 877 ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps);CHKERRQ(ierr); 878 smoother = ksps[0]; 879 ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 880 ierr = PetscFree(ksps);CHKERRQ(ierr); 881 } 882 883 /* do my own cheby */ 884 ierr = PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &flag);CHKERRQ(ierr); 885 if (flag) { 886 PetscReal emax, emin; 887 ierr = PetscObjectTypeCompare((PetscObject)subpc, PCJACOBI, &flag);CHKERRQ(ierr); 888 if (flag && emaxs[level] > 0.0) emax=emaxs[level]; /* eigen estimate only for diagnal PC */ 889 else { /* eigen estimate 'emax' */ 890 KSP eksp; 891 Mat Lmat = Aarr[level]; 892 Vec bb, xx; 893 894 ierr = MatGetVecs(Lmat, &bb, 0);CHKERRQ(ierr); 895 ierr = MatGetVecs(Lmat, &xx, 0);CHKERRQ(ierr); 896 { 897 PetscRandom rctx; 898 ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 899 ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 900 ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 901 ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 902 } 903 904 /* zeroing out BC rows -- needed for crazy matrices */ 905 { 906 PetscInt Istart,Iend,ncols,jj,Ii; 907 PetscScalar zero = 0.0; 908 ierr = MatGetOwnershipRange(Lmat, &Istart, &Iend);CHKERRQ(ierr); 909 for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) { 910 ierr = MatGetRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr); 911 if (ncols <= 1) { 912 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr); 913 } 914 ierr = MatRestoreRow(Lmat,Ii,&ncols,0,0);CHKERRQ(ierr); 915 } 916 ierr = VecAssemblyBegin(bb);CHKERRQ(ierr); 917 ierr = VecAssemblyEnd(bb);CHKERRQ(ierr); 918 } 919 920 ierr = KSPCreate(wcomm, &eksp);CHKERRQ(ierr); 921 ierr = KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10);CHKERRQ(ierr); 922 ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 923 ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 924 ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 925 ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 926 927 ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 928 ierr = KSPSetOperators(eksp, Lmat, Lmat, SAME_NONZERO_PATTERN);CHKERRQ(ierr); 929 ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 930 931 /* set PC type to be same as smoother */ 932 ierr = KSPSetPC(eksp, subpc);CHKERRQ(ierr); 933 934 /* solve - keep stuff out of logging */ 935 ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 936 ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 937 ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 938 ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 939 ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 940 941 ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 942 943 ierr = VecDestroy(&xx);CHKERRQ(ierr); 944 ierr = VecDestroy(&bb);CHKERRQ(ierr); 945 ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 946 947 if (pc_gamg->verbose > 0) { 948 PetscInt N1, tt; 949 ierr = MatGetSize(Aarr[level], &N1, &tt);CHKERRQ(ierr); 950 PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1); 951 } 952 } 953 { 954 PetscInt N1, N0; 955 ierr = MatGetSize(Aarr[level], &N1, NULL);CHKERRQ(ierr); 956 ierr = MatGetSize(Aarr[level+1], &N0, NULL);CHKERRQ(ierr); 957 /* heuristic - is this crap? */ 958 /* emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); */ 959 emin = emax * pc_gamg->eigtarget[0]; 960 emax *= pc_gamg->eigtarget[1]; 961 } 962 ierr = KSPChebyshevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); 963 } /* setup checby flag */ 964 } /* non-coarse levels */ 965 966 /* clean up */ 967 for (level=1; level<pc_gamg->Nlevels; level++) { 968 ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 969 ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 970 } 971 972 ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 973 974 if (PETSC_TRUE) { 975 KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 976 ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 977 ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 978 } 979 } else { 980 KSP smoother; 981 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",rank,__FUNCT__); 982 ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 983 ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN);CHKERRQ(ierr); 984 ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 985 ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 986 } 987 PetscFunctionReturn(0); 988 } 989 990 /* ------------------------------------------------------------------------- */ 991 /* 992 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 993 that was created with PCCreate_GAMG(). 994 995 Input Parameter: 996 . pc - the preconditioner context 997 998 Application Interface Routine: PCDestroy() 999 */ 1000 #undef __FUNCT__ 1001 #define __FUNCT__ "PCDestroy_GAMG" 1002 PetscErrorCode PCDestroy_GAMG(PC pc) 1003 { 1004 PetscErrorCode ierr; 1005 PC_MG *mg = (PC_MG*)pc->data; 1006 PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 1007 1008 PetscFunctionBegin; 1009 ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 1010 ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 1011 ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 1012 PetscFunctionReturn(0); 1013 } 1014 1015 1016 #undef __FUNCT__ 1017 #define __FUNCT__ "PCGAMGSetProcEqLim" 1018 /*@ 1019 PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 1020 processor reduction. 1021 1022 Not Collective on PC 1023 1024 Input Parameters: 1025 . pc - the preconditioner context 1026 1027 1028 Options Database Key: 1029 . -pc_gamg_process_eq_limit 1030 1031 Level: intermediate 1032 1033 Concepts: Unstructured multrigrid preconditioner 1034 1035 .seealso: () 1036 @*/ 1037 PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 1038 { 1039 PetscErrorCode ierr; 1040 1041 PetscFunctionBegin; 1042 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1043 ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1044 PetscFunctionReturn(0); 1045 } 1046 1047 EXTERN_C_BEGIN 1048 #undef __FUNCT__ 1049 #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 1050 PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 1051 { 1052 PC_MG *mg = (PC_MG*)pc->data; 1053 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1054 1055 PetscFunctionBegin; 1056 if (n>0) pc_gamg->min_eq_proc = n; 1057 PetscFunctionReturn(0); 1058 } 1059 EXTERN_C_END 1060 1061 #undef __FUNCT__ 1062 #define __FUNCT__ "PCGAMGSetCoarseEqLim" 1063 /*@ 1064 PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 1065 1066 Collective on PC 1067 1068 Input Parameters: 1069 . pc - the preconditioner context 1070 1071 1072 Options Database Key: 1073 . -pc_gamg_coarse_eq_limit 1074 1075 Level: intermediate 1076 1077 Concepts: Unstructured multrigrid preconditioner 1078 1079 .seealso: () 1080 @*/ 1081 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 1082 { 1083 PetscErrorCode ierr; 1084 1085 PetscFunctionBegin; 1086 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1087 ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1088 PetscFunctionReturn(0); 1089 } 1090 1091 EXTERN_C_BEGIN 1092 #undef __FUNCT__ 1093 #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 1094 PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 1095 { 1096 PC_MG *mg = (PC_MG*)pc->data; 1097 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1098 1099 PetscFunctionBegin; 1100 if (n>0) pc_gamg->coarse_eq_limit = n; 1101 PetscFunctionReturn(0); 1102 } 1103 EXTERN_C_END 1104 1105 #undef __FUNCT__ 1106 #define __FUNCT__ "PCGAMGSetRepartitioning" 1107 /*@ 1108 PCGAMGSetRepartitioning - Repartition the coarse grids 1109 1110 Collective on PC 1111 1112 Input Parameters: 1113 . pc - the preconditioner context 1114 1115 1116 Options Database Key: 1117 . -pc_gamg_repartition 1118 1119 Level: intermediate 1120 1121 Concepts: Unstructured multrigrid preconditioner 1122 1123 .seealso: () 1124 @*/ 1125 PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 1126 { 1127 PetscErrorCode ierr; 1128 1129 PetscFunctionBegin; 1130 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1131 ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1132 PetscFunctionReturn(0); 1133 } 1134 1135 EXTERN_C_BEGIN 1136 #undef __FUNCT__ 1137 #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 1138 PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 1139 { 1140 PC_MG *mg = (PC_MG*)pc->data; 1141 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1142 1143 PetscFunctionBegin; 1144 pc_gamg->repart = n; 1145 PetscFunctionReturn(0); 1146 } 1147 EXTERN_C_END 1148 1149 #undef __FUNCT__ 1150 #define __FUNCT__ "PCGAMGSetReuseProl" 1151 /*@ 1152 PCGAMGSetReuseProl - Reuse prlongation 1153 1154 Collective on PC 1155 1156 Input Parameters: 1157 . pc - the preconditioner context 1158 1159 1160 Options Database Key: 1161 . -pc_gamg_reuse_interpolation 1162 1163 Level: intermediate 1164 1165 Concepts: Unstructured multrigrid preconditioner 1166 1167 .seealso: () 1168 @*/ 1169 PetscErrorCode PCGAMGSetReuseProl(PC pc, PetscBool n) 1170 { 1171 PetscErrorCode ierr; 1172 1173 PetscFunctionBegin; 1174 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1175 ierr = PetscTryMethod(pc,"PCGAMGSetReuseProl_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1176 PetscFunctionReturn(0); 1177 } 1178 1179 EXTERN_C_BEGIN 1180 #undef __FUNCT__ 1181 #define __FUNCT__ "PCGAMGSetReuseProl_GAMG" 1182 PetscErrorCode PCGAMGSetReuseProl_GAMG(PC pc, PetscBool n) 1183 { 1184 PC_MG *mg = (PC_MG*)pc->data; 1185 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1186 1187 PetscFunctionBegin; 1188 pc_gamg->reuse_prol = n; 1189 PetscFunctionReturn(0); 1190 } 1191 EXTERN_C_END 1192 1193 #undef __FUNCT__ 1194 #define __FUNCT__ "PCGAMGSetUseASMAggs" 1195 /*@ 1196 PCGAMGSetUseASMAggs - 1197 1198 Collective on PC 1199 1200 Input Parameters: 1201 . pc - the preconditioner context 1202 1203 1204 Options Database Key: 1205 . -pc_gamg_use_agg_gasm 1206 1207 Level: intermediate 1208 1209 Concepts: Unstructured multrigrid preconditioner 1210 1211 .seealso: () 1212 @*/ 1213 PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n) 1214 { 1215 PetscErrorCode ierr; 1216 1217 PetscFunctionBegin; 1218 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1219 ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1220 PetscFunctionReturn(0); 1221 } 1222 1223 EXTERN_C_BEGIN 1224 #undef __FUNCT__ 1225 #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG" 1226 PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n) 1227 { 1228 PC_MG *mg = (PC_MG*)pc->data; 1229 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1230 1231 PetscFunctionBegin; 1232 pc_gamg->use_aggs_in_gasm = n; 1233 PetscFunctionReturn(0); 1234 } 1235 EXTERN_C_END 1236 1237 #undef __FUNCT__ 1238 #define __FUNCT__ "PCGAMGSetNlevels" 1239 /*@ 1240 PCGAMGSetNlevels - 1241 1242 Not collective on PC 1243 1244 Input Parameters: 1245 . pc - the preconditioner context 1246 1247 1248 Options Database Key: 1249 . -pc_mg_levels 1250 1251 Level: intermediate 1252 1253 Concepts: Unstructured multrigrid preconditioner 1254 1255 .seealso: () 1256 @*/ 1257 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 1258 { 1259 PetscErrorCode ierr; 1260 1261 PetscFunctionBegin; 1262 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1263 ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1264 PetscFunctionReturn(0); 1265 } 1266 1267 EXTERN_C_BEGIN 1268 #undef __FUNCT__ 1269 #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 1270 PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 1271 { 1272 PC_MG *mg = (PC_MG*)pc->data; 1273 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1274 1275 PetscFunctionBegin; 1276 pc_gamg->Nlevels = n; 1277 PetscFunctionReturn(0); 1278 } 1279 EXTERN_C_END 1280 1281 #undef __FUNCT__ 1282 #define __FUNCT__ "PCGAMGSetThreshold" 1283 /*@ 1284 PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 1285 1286 Not collective on PC 1287 1288 Input Parameters: 1289 . pc - the preconditioner context 1290 1291 1292 Options Database Key: 1293 . -pc_gamg_threshold 1294 1295 Level: intermediate 1296 1297 Concepts: Unstructured multrigrid preconditioner 1298 1299 .seealso: () 1300 @*/ 1301 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 1302 { 1303 PetscErrorCode ierr; 1304 1305 PetscFunctionBegin; 1306 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1307 ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 1308 PetscFunctionReturn(0); 1309 } 1310 1311 EXTERN_C_BEGIN 1312 #undef __FUNCT__ 1313 #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 1314 PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 1315 { 1316 PC_MG *mg = (PC_MG*)pc->data; 1317 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1318 1319 PetscFunctionBegin; 1320 pc_gamg->threshold = n; 1321 PetscFunctionReturn(0); 1322 } 1323 EXTERN_C_END 1324 1325 #undef __FUNCT__ 1326 #define __FUNCT__ "PCGAMGSetType" 1327 /*@ 1328 PCGAMGSetType - Set solution method - calls sub create method 1329 1330 Collective on PC 1331 1332 Input Parameters: 1333 . pc - the preconditioner context 1334 1335 1336 Options Database Key: 1337 . -pc_gamg_type 1338 1339 Level: intermediate 1340 1341 Concepts: Unstructured multrigrid preconditioner 1342 1343 .seealso: () 1344 @*/ 1345 PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1346 { 1347 PetscErrorCode ierr; 1348 1349 PetscFunctionBegin; 1350 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1351 ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 1352 PetscFunctionReturn(0); 1353 } 1354 1355 EXTERN_C_BEGIN 1356 #undef __FUNCT__ 1357 #define __FUNCT__ "PCGAMGSetType_GAMG" 1358 PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1359 { 1360 PetscErrorCode ierr,(*r)(PC); 1361 1362 PetscFunctionBegin; 1363 ierr = PetscFunctionListFind(((PetscObject)pc)->comm,GAMGList,type,PETSC_FALSE,(PetscVoidStarFunction)&r);CHKERRQ(ierr); 1364 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 1365 ierr = (*r)(pc);CHKERRQ(ierr); 1366 PetscFunctionReturn(0); 1367 } 1368 EXTERN_C_END 1369 1370 #undef __FUNCT__ 1371 #define __FUNCT__ "PCSetFromOptions_GAMG" 1372 PetscErrorCode PCSetFromOptions_GAMG(PC pc) 1373 { 1374 PetscErrorCode ierr; 1375 PC_MG *mg = (PC_MG*)pc->data; 1376 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1377 PetscBool flag; 1378 PetscInt two = 2; 1379 MPI_Comm wcomm = ((PetscObject)pc)->comm; 1380 1381 PetscFunctionBegin; 1382 ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr); 1383 { 1384 /* -pc_gamg_verbose */ 1385 ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG", 1386 "none", pc_gamg->verbose, 1387 &pc_gamg->verbose, NULL);CHKERRQ(ierr); 1388 /* -pc_gamg_repartition */ 1389 ierr = PetscOptionsBool("-pc_gamg_repartition", 1390 "Repartion coarse grids (false)", 1391 "PCGAMGRepartitioning", 1392 pc_gamg->repart, 1393 &pc_gamg->repart, 1394 &flag);CHKERRQ(ierr); 1395 /* -pc_gamg_reuse_interpolation */ 1396 ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation", 1397 "Reuse prolongation operator (true)", 1398 "PCGAMGReuseProl", 1399 pc_gamg->reuse_prol, 1400 &pc_gamg->reuse_prol, 1401 &flag);CHKERRQ(ierr); 1402 /* -pc_gamg_use_agg_gasm */ 1403 ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm", 1404 "Use aggregation agragates for GASM smoother (false)", 1405 "PCGAMGUseASMAggs", 1406 pc_gamg->use_aggs_in_gasm, 1407 &pc_gamg->use_aggs_in_gasm, 1408 &flag);CHKERRQ(ierr); 1409 /* -pc_gamg_process_eq_limit */ 1410 ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", 1411 "Limit (goal) on number of equations per process on coarse grids", 1412 "PCGAMGSetProcEqLim", 1413 pc_gamg->min_eq_proc, 1414 &pc_gamg->min_eq_proc, 1415 &flag);CHKERRQ(ierr); 1416 /* -pc_gamg_coarse_eq_limit */ 1417 ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit", 1418 "Limit on number of equations for the coarse grid", 1419 "PCGAMGSetCoarseEqLim", 1420 pc_gamg->coarse_eq_limit, 1421 &pc_gamg->coarse_eq_limit, 1422 &flag);CHKERRQ(ierr); 1423 /* -pc_gamg_threshold */ 1424 ierr = PetscOptionsReal("-pc_gamg_threshold", 1425 "Relative threshold to use for dropping edges in aggregation graph", 1426 "PCGAMGSetThreshold", 1427 pc_gamg->threshold, 1428 &pc_gamg->threshold, 1429 &flag);CHKERRQ(ierr); 1430 if (flag && pc_gamg->verbose) { 1431 ierr = PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold);CHKERRQ(ierr); 1432 } 1433 1434 ierr = PetscOptionsRealArray("-pc_gamg_eigtarget","Target eigenvalue range as fraction of estimated maximum eigenvalue","PCGAMGSetEigTarget",pc_gamg->eigtarget,&two,NULL);CHKERRQ(ierr); 1435 ierr = PetscOptionsInt("-pc_mg_levels", 1436 "Set number of MG levels", 1437 "PCGAMGSetNlevels", 1438 pc_gamg->Nlevels, 1439 &pc_gamg->Nlevels, 1440 &flag);CHKERRQ(ierr); 1441 } 1442 ierr = PetscOptionsTail();CHKERRQ(ierr); 1443 PetscFunctionReturn(0); 1444 } 1445 1446 /* -------------------------------------------------------------------------- */ 1447 /*MC 1448 PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework. 1449 - This is the entry point to GAMG, registered in pcregis.c 1450 1451 Options Database Keys: 1452 Multigrid options(inherited) 1453 + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1454 . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1455 . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 1456 - -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade 1457 1458 Level: intermediate 1459 1460 Concepts: multigrid 1461 1462 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1463 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1464 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1465 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1466 PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1467 M*/ 1468 EXTERN_C_BEGIN 1469 #undef __FUNCT__ 1470 #define __FUNCT__ "PCCreate_GAMG" 1471 PetscErrorCode PCCreate_GAMG(PC pc) 1472 { 1473 PetscErrorCode ierr; 1474 PC_GAMG *pc_gamg; 1475 PC_MG *mg; 1476 #if defined PETSC_GAMG_USE_LOG 1477 static long count = 0; 1478 #endif 1479 1480 PetscFunctionBegin; 1481 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 1482 ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 1483 ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr); 1484 1485 /* create a supporting struct and attach it to pc */ 1486 ierr = PetscNewLog(pc, PC_GAMG, &pc_gamg);CHKERRQ(ierr); 1487 mg = (PC_MG*)pc->data; 1488 mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 1489 mg->innerctx = pc_gamg; 1490 1491 pc_gamg->setup_count = 0; 1492 /* these should be in subctx but repartitioning needs simple arrays */ 1493 pc_gamg->data_sz = 0; 1494 pc_gamg->data = 0; 1495 1496 /* register AMG type */ 1497 if (!GAMGList) { 1498 ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void (*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr); 1499 ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void (*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr); 1500 } 1501 1502 /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 1503 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 1504 pc->ops->setup = PCSetUp_GAMG; 1505 pc->ops->reset = PCReset_GAMG; 1506 pc->ops->destroy = PCDestroy_GAMG; 1507 1508 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1509 "PCGAMGSetProcEqLim_C", 1510 "PCGAMGSetProcEqLim_GAMG", 1511 PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr); 1512 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1513 "PCGAMGSetCoarseEqLim_C", 1514 "PCGAMGSetCoarseEqLim_GAMG", 1515 PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr); 1516 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1517 "PCGAMGSetRepartitioning_C", 1518 "PCGAMGSetRepartitioning_GAMG", 1519 PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr); 1520 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1521 "PCGAMGSetReuseProl_C", 1522 "PCGAMGSetReuseProl_GAMG", 1523 PCGAMGSetReuseProl_GAMG);CHKERRQ(ierr); 1524 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1525 "PCGAMGSetUseASMAggs_C", 1526 "PCGAMGSetUseASMAggs_GAMG", 1527 PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr); 1528 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1529 "PCGAMGSetThreshold_C", 1530 "PCGAMGSetThreshold_GAMG", 1531 PCGAMGSetThreshold_GAMG);CHKERRQ(ierr); 1532 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc, 1533 "PCGAMGSetType_C", 1534 "PCGAMGSetType_GAMG", 1535 PCGAMGSetType_GAMG);CHKERRQ(ierr); 1536 pc_gamg->repart = PETSC_FALSE; 1537 pc_gamg->reuse_prol = PETSC_TRUE; 1538 pc_gamg->use_aggs_in_gasm = PETSC_FALSE; 1539 pc_gamg->min_eq_proc = 50; 1540 pc_gamg->coarse_eq_limit = 800; 1541 pc_gamg->threshold = 0.001; 1542 pc_gamg->Nlevels = GAMG_MAXLEVELS; 1543 pc_gamg->verbose = 0; 1544 pc_gamg->emax_id = -1; 1545 pc_gamg->eigtarget[0] = 0.05; 1546 pc_gamg->eigtarget[1] = 1.05; 1547 1548 /* private events */ 1549 #if defined PETSC_GAMG_USE_LOG 1550 if (count++ == 0) { 1551 ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr); 1552 ierr = PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr); 1553 /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 1554 /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 1555 /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 1556 ierr = PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr); 1557 ierr = PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr); 1558 ierr = PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr); 1559 ierr = PetscLogEventRegister(" search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr); 1560 ierr = PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr); 1561 ierr = PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr); 1562 ierr = PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr); 1563 ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr); 1564 ierr = PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr); 1565 ierr = PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr); 1566 ierr = PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr); 1567 ierr = PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr); 1568 1569 /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 1570 /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 1571 /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 1572 /* create timer stages */ 1573 #if defined GAMG_STAGES 1574 { 1575 char str[32]; 1576 PetscInt lidx; 1577 sprintf(str,"MG Level %d (finest)",0); 1578 ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr); 1579 for (lidx=1; lidx<9; lidx++) { 1580 sprintf(str,"MG Level %d",lidx); 1581 ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr); 1582 } 1583 } 1584 #endif 1585 } 1586 #endif 1587 /* general events */ 1588 #if defined PETSC_USE_LOG 1589 ierr = PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG);CHKERRQ(ierr); 1590 ierr = PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO);CHKERRQ(ierr); 1591 ierr = PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr); 1592 ierr = PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr); 1593 ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr); 1594 ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr); 1595 ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG);CHKERRQ(ierr); 1596 ierr = PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG);CHKERRQ(ierr); 1597 #endif 1598 1599 /* instantiate derived type */ 1600 ierr = PetscOptionsHead("GAMG options");CHKERRQ(ierr); 1601 { 1602 char tname[256] = GAMGAGG; 1603 ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType",GAMGList, tname, tname, sizeof(tname), NULL);CHKERRQ(ierr); 1604 ierr = PCGAMGSetType(pc, tname);CHKERRQ(ierr); 1605 } 1606 ierr = PetscOptionsTail();CHKERRQ(ierr); 1607 PetscFunctionReturn(0); 1608 } 1609 EXTERN_C_END 1610