Actual source code: gamg.c
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>
6: #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h>
8: #if defined(PETSC_HAVE_CUDA)
9: #include <cuda_runtime.h>
10: #endif
12: #if defined(PETSC_HAVE_HIP)
13: #include <hip/hip_runtime.h>
14: #endif
16: PetscLogEvent petsc_gamg_setup_events[NUM_SET];
17: PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3];
18: PetscLogEvent PC_GAMGGraph_AGG;
19: PetscLogEvent PC_GAMGGraph_GEO;
20: PetscLogEvent PC_GAMGCoarsen_AGG;
21: PetscLogEvent PC_GAMGCoarsen_GEO;
22: PetscLogEvent PC_GAMGProlongator_AGG;
23: PetscLogEvent PC_GAMGProlongator_GEO;
24: PetscLogEvent PC_GAMGOptProlongator_AGG;
26: /* #define GAMG_STAGES */
27: #if defined(GAMG_STAGES)
28: static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
29: #endif
31: static PetscFunctionList GAMGList = NULL;
32: static PetscBool PCGAMGPackageInitialized;
34: /* ----------------------------------------------------------------------------- */
35: PetscErrorCode PCReset_GAMG(PC pc)
36: {
37: PetscErrorCode ierr, level;
38: PC_MG *mg = (PC_MG*)pc->data;
39: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
42: PetscFree(pc_gamg->data);
43: pc_gamg->data_sz = 0;
44: PetscFree(pc_gamg->orig_data);
45: for (level = 0; level < PETSC_MG_MAXLEVELS ; level++) {
46: mg->min_eigen_DinvA[level] = 0;
47: mg->max_eigen_DinvA[level] = 0;
48: }
49: pc_gamg->emin = 0;
50: pc_gamg->emax = 0;
51: return(0);
52: }
54: /* -------------------------------------------------------------------------- */
55: /*
56: PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number
57: of active processors.
59: Input Parameter:
60: . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
61: 'pc_gamg->data_sz' are changed via repartitioning/reduction.
62: . Amat_fine - matrix on this fine (k) level
63: . cr_bs - coarse block size
64: In/Output Parameter:
65: . a_P_inout - prolongation operator to the next level (k-->k-1)
66: . a_nactive_proc - number of active procs
67: Output Parameter:
68: . a_Amat_crs - coarse matrix that is created (k-1)
69: */
71: static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,IS * Pcolumnperm, PetscBool is_last)
72: {
73: PetscErrorCode ierr;
74: PC_MG *mg = (PC_MG*)pc->data;
75: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
76: Mat Cmat,Pold=*a_P_inout;
77: MPI_Comm comm;
78: PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc;
79: PetscInt ncrs_eq,ncrs,f_bs;
82: PetscObjectGetComm((PetscObject)Amat_fine,&comm);
83: MPI_Comm_rank(comm, &rank);
84: MPI_Comm_size(comm, &size);
85: MatGetBlockSize(Amat_fine, &f_bs);
86: PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0);
87: MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);
88: PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0);
90: if (Pcolumnperm) *Pcolumnperm = NULL;
92: /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
93: MatGetLocalSize(Cmat, &ncrs_eq, NULL);
94: if (pc_gamg->data_cell_rows>0) {
95: ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
96: } else {
97: PetscInt bs;
98: MatGetBlockSize(Cmat, &bs);
99: ncrs = ncrs_eq/bs;
100: }
101: /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
102: if (pc_gamg->level_reduction_factors[pc_gamg->current_level] == 0 && PetscDefined(HAVE_CUDA) && pc_gamg->current_level==0) { /* 0 turns reducing to 1 process/device on; do for HIP, etc. */
103: #if defined(PETSC_HAVE_CUDA)
104: PetscShmComm pshmcomm;
105: PetscMPIInt locrank;
106: MPI_Comm loccomm;
107: PetscInt s_nnodes,r_nnodes, new_new_size;
108: cudaError_t cerr;
109: int devCount;
110: PetscShmCommGet(comm,&pshmcomm);
111: PetscShmCommGetMpiShmComm(pshmcomm,&loccomm);
112: MPI_Comm_rank(loccomm, &locrank);
113: s_nnodes = !locrank;
114: MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm);
115: if (size%r_nnodes) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of nodes np=%D nnodes%D",size,r_nnodes);
116: devCount = 0;
117: cerr = cudaGetDeviceCount(&devCount);
118: cudaGetLastError(); /* Reset the last error */
119: if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */
120: new_new_size = r_nnodes * devCount;
121: new_size = new_new_size;
122: PetscInfo5(pc,"Fine grid with Cuda. %D nodes. Change new active set size %d --> %d (devCount=%d #nodes=%D)\n",r_nnodes,nactive,new_size,devCount,r_nnodes);
123: } else {
124: PetscInfo(pc,"With Cuda but no device. Use heuristics.");
125: goto HEURISTIC;
126: }
127: #else
128: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here");
129: #endif
130: } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) {
131: if (nactive%pc_gamg->level_reduction_factors[pc_gamg->current_level]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of active process %D wrt reduction factor %D",nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]);
132: new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level];
133: PetscInfo3(pc,"Manually setting reduction to %d active processes (%d/%D)\n",new_size,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]);
134: } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) {
135: new_size = 1;
136: PetscInfo1(pc,"Force coarsest grid reduction to %d active processes\n",new_size);
137: } else {
138: PetscInt ncrs_eq_glob;
139: #if defined(PETSC_HAVE_CUDA)
140: HEURISTIC:
141: #endif
142: MatGetSize(Cmat, &ncrs_eq_glob, NULL);
143: new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
144: if (!new_size) new_size = 1; /* not likely, posible? */
145: else if (new_size >= nactive) new_size = nactive; /* no change, rare */
146: PetscInfo2(pc,"Coarse grid reduction from %d to %d active processes\n",nactive,new_size);
147: }
148: if (new_size==nactive) {
149: *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
150: if (new_size < size) {
151: /* odd case where multiple coarse grids are on one processor or no coarsening ... */
152: PetscInfo1(pc,"reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",nactive);
153: if (pc_gamg->cpu_pin_coarse_grids) {
154: MatBindToCPU(*a_Amat_crs,PETSC_TRUE);
155: MatBindToCPU(*a_P_inout,PETSC_TRUE);
156: }
157: }
158: /* we know that the grid structure can be reused in MatPtAP */
159: } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
160: PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1;
161: IS is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
162: nloc_old = ncrs_eq/cr_bs;
163: 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);
164: /* get new_size and rfactor */
165: if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
166: /* find factor */
167: if (new_size == 1) rfactor = size; /* don't modify */
168: else {
169: PetscReal best_fact = 0.;
170: jj = -1;
171: for (kk = 1 ; kk <= size ; kk++) {
172: if (!(size%kk)) { /* a candidate */
173: PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
174: if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
175: if (fact > best_fact) {
176: best_fact = fact; jj = kk;
177: }
178: }
179: }
180: if (jj != -1) rfactor = jj;
181: else rfactor = 1; /* a prime */
182: if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
183: else expand_factor = rfactor;
184: }
185: new_size = size/rfactor; /* make new size one that is factor */
186: if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
187: *a_Amat_crs = Cmat;
188: PetscInfo2(pc,"Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",new_size,ncrs_eq);
189: return(0);
190: }
191: }
192: PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);
193: /* make 'is_eq_newproc' */
194: PetscMalloc1(size, &counts);
195: if (pc_gamg->repart) {
196: /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */
197: Mat adj;
198: PetscInfo4(pc,"Repartition: size (active): %d --> %d, %D local equations, using %s process layout\n",*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread");
199: /* get 'adj' */
200: if (cr_bs == 1) {
201: MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);
202: } else {
203: /* make a scalar matrix to partition (no Stokes here) */
204: Mat tMat;
205: PetscInt Istart_crs,Iend_crs,ncols,jj,Ii;
206: const PetscScalar *vals;
207: const PetscInt *idx;
208: PetscInt *d_nnz, *o_nnz, M, N;
209: static PetscInt llev = 0; /* ugly but just used for debugging */
210: MatType mtype;
212: PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);
213: MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);
214: MatGetSize(Cmat, &M, &N);
215: for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
216: MatGetRow(Cmat,Ii,&ncols,NULL,NULL);
217: d_nnz[jj] = ncols/cr_bs;
218: o_nnz[jj] = ncols/cr_bs;
219: MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL);
220: if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
221: if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
222: }
224: MatGetType(Amat_fine,&mtype);
225: MatCreate(comm, &tMat);
226: MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);
227: MatSetType(tMat,mtype);
228: MatSeqAIJSetPreallocation(tMat,0,d_nnz);
229: MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);
230: PetscFree2(d_nnz,o_nnz);
232: for (ii = Istart_crs; ii < Iend_crs; ii++) {
233: PetscInt dest_row = ii/cr_bs;
234: MatGetRow(Cmat,ii,&ncols,&idx,&vals);
235: for (jj = 0; jj < ncols; jj++) {
236: PetscInt dest_col = idx[jj]/cr_bs;
237: PetscScalar v = 1.0;
238: MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);
239: }
240: MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);
241: }
242: MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);
243: MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);
245: if (llev++ == -1) {
246: PetscViewer viewer; char fname[32];
247: PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);
248: PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
249: MatView(tMat, viewer);
250: PetscViewerDestroy(&viewer);
251: }
252: MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);
253: MatDestroy(&tMat);
254: } /* create 'adj' */
256: { /* partition: get newproc_idx */
257: char prefix[256];
258: const char *pcpre;
259: const PetscInt *is_idx;
260: MatPartitioning mpart;
261: IS proc_is;
263: MatPartitioningCreate(comm, &mpart);
264: MatPartitioningSetAdjacency(mpart, adj);
265: PCGetOptionsPrefix(pc, &pcpre);
266: PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
267: PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);
268: MatPartitioningSetFromOptions(mpart);
269: MatPartitioningSetNParts(mpart, new_size);
270: MatPartitioningApply(mpart, &proc_is);
271: MatPartitioningDestroy(&mpart);
273: /* collect IS info */
274: PetscMalloc1(ncrs_eq, &newproc_idx);
275: ISGetIndices(proc_is, &is_idx);
276: for (kk = jj = 0 ; kk < nloc_old ; kk++) {
277: for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
278: newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
279: }
280: }
281: ISRestoreIndices(proc_is, &is_idx);
282: ISDestroy(&proc_is);
283: }
284: MatDestroy(&adj);
286: ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);
287: PetscFree(newproc_idx);
288: } else { /* simple aggregation of parts -- 'is_eq_newproc' */
289: PetscInt targetPE;
290: if (new_size==nactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen");
291: PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);
292: targetPE = (rank/rfactor)*expand_factor;
293: ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);
294: } /* end simple 'is_eq_newproc' */
296: /*
297: Create an index set from the is_eq_newproc index set to indicate the mapping TO
298: */
299: ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);
300: is_eq_num_prim = is_eq_num;
301: /*
302: Determine how many equations/vertices are assigned to each processor
303: */
304: ISPartitioningCount(is_eq_newproc, size, counts);
305: ncrs_eq_new = counts[rank];
306: ISDestroy(&is_eq_newproc);
307: ncrs_new = ncrs_eq_new/cr_bs;
309: PetscFree(counts);
310: /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxilary data into some complex abstracted thing */
311: {
312: Vec src_crd, dest_crd;
313: const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
314: VecScatter vecscat;
315: PetscScalar *array;
316: IS isscat;
317: /* move data (for primal equations only) */
318: /* Create a vector to contain the newly ordered element information */
319: VecCreate(comm, &dest_crd);
320: VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);
321: VecSetType(dest_crd,VECSTANDARD); /* this is needed! */
322: /*
323: There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
324: a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'.
325: */
326: PetscMalloc1(ncrs*node_data_sz, &tidx);
327: ISGetIndices(is_eq_num_prim, &idx);
328: for (ii=0,jj=0; ii<ncrs; ii++) {
329: PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
330: for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
331: }
332: ISRestoreIndices(is_eq_num_prim, &idx);
333: ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);
334: PetscFree(tidx);
335: /*
336: Create a vector to contain the original vertex information for each element
337: */
338: VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);
339: for (jj=0; jj<ndata_cols; jj++) {
340: const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
341: for (ii=0; ii<ncrs; ii++) {
342: for (kk=0; kk<ndata_rows; kk++) {
343: PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
344: PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
345: VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);
346: }
347: }
348: }
349: VecAssemblyBegin(src_crd);
350: VecAssemblyEnd(src_crd);
351: /*
352: Scatter the element vertex information (still in the original vertex ordering)
353: to the correct processor
354: */
355: VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);
356: ISDestroy(&isscat);
357: VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);
358: VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);
359: VecScatterDestroy(&vecscat);
360: VecDestroy(&src_crd);
361: /*
362: Put the element vertex data into a new allocation of the gdata->ele
363: */
364: PetscFree(pc_gamg->data);
365: PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);
367: pc_gamg->data_sz = node_data_sz*ncrs_new;
368: strideNew = ncrs_new*ndata_rows;
370: VecGetArray(dest_crd, &array);
371: for (jj=0; jj<ndata_cols; jj++) {
372: for (ii=0; ii<ncrs_new; ii++) {
373: for (kk=0; kk<ndata_rows; kk++) {
374: PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
375: pc_gamg->data[ix] = PetscRealPart(array[jx]);
376: }
377: }
378: }
379: VecRestoreArray(dest_crd, &array);
380: VecDestroy(&dest_crd);
381: }
382: /* move A and P (columns) with new layout */
383: PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);
384: /*
385: Invert for MatCreateSubMatrix
386: */
387: ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);
388: ISSort(new_eq_indices); /* is this needed? */
389: ISSetBlockSize(new_eq_indices, cr_bs);
390: if (is_eq_num != is_eq_num_prim) {
391: ISDestroy(&is_eq_num_prim); /* could be same as 'is_eq_num' */
392: }
393: if (Pcolumnperm) {
394: PetscObjectReference((PetscObject)new_eq_indices);
395: *Pcolumnperm = new_eq_indices;
396: }
397: ISDestroy(&is_eq_num);
398: PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);
399: PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);
400: /* 'a_Amat_crs' output */
401: {
402: Mat mat;
403: MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);
404: *a_Amat_crs = mat;
405: }
406: MatDestroy(&Cmat);
408: PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);
409: /* prolongator */
410: {
411: IS findices;
412: PetscInt Istart,Iend;
413: Mat Pnew;
415: MatGetOwnershipRange(Pold, &Istart, &Iend);
416: PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);
417: ISCreateStride(comm,Iend-Istart,Istart,1,&findices);
418: ISSetBlockSize(findices,f_bs);
419: MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);
420: ISDestroy(&findices);
421: MatSetOption(Pnew,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);
423: PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);
424: MatDestroy(a_P_inout);
426: /* output - repartitioned */
427: *a_P_inout = Pnew;
428: }
429: ISDestroy(&new_eq_indices);
431: *a_nactive_proc = new_size; /* output */
433: /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
434: if (pc_gamg->cpu_pin_coarse_grids) {
435: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
436: static PetscInt llev = 2;
437: PetscInfo1(pc,"Pinning level %D to the CPU\n",llev++);
438: #endif
439: MatBindToCPU(*a_Amat_crs,PETSC_TRUE);
440: MatBindToCPU(*a_P_inout,PETSC_TRUE);
441: if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
442: Mat A = *a_Amat_crs, P = *a_P_inout;
443: PetscMPIInt size;
444: MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
445: if (size > 1) {
446: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data;
447: VecBindToCPU(a->lvec,PETSC_TRUE);
448: VecBindToCPU(p->lvec,PETSC_TRUE);
449: }
450: }
451: }
452: PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);
453: }
454: return(0);
455: }
457: PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2)
458: {
460: const char *prefix;
461: char addp[32];
462: PC_MG *mg = (PC_MG*)a_pc->data;
463: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
466: PCGetOptionsPrefix(a_pc,&prefix);
467: PetscInfo1(a_pc,"Square Graph on level %D\n",pc_gamg->current_level+1);
468: MatProductCreate(Gmat1,Gmat1,NULL,Gmat2);
469: MatSetOptionsPrefix(*Gmat2,prefix);
470: PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level);
471: MatAppendOptionsPrefix(*Gmat2,addp);
472: if ((*Gmat2)->structurally_symmetric) {
473: MatProductSetType(*Gmat2,MATPRODUCT_AB);
474: } else {
475: MatSetOption(Gmat1,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);
476: MatProductSetType(*Gmat2,MATPRODUCT_AtB);
477: }
478: MatProductSetFromOptions(*Gmat2);
479: PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0);
480: MatProductSymbolic(*Gmat2);
481: PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0);
482: MatProductClear(*Gmat2);
483: /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
484: (*Gmat2)->assembled = PETSC_TRUE;
485: return(0);
486: }
488: /* -------------------------------------------------------------------------- */
489: /*
490: PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
491: by setting data structures and options.
493: Input Parameter:
494: . pc - the preconditioner context
496: */
497: PetscErrorCode PCSetUp_GAMG(PC pc)
498: {
500: PC_MG *mg = (PC_MG*)pc->data;
501: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
502: Mat Pmat = pc->pmat;
503: PetscInt fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS];
504: MPI_Comm comm;
505: PetscMPIInt rank,size,nactivepe;
506: Mat Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS];
507: IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
508: PetscLogDouble nnz0=0.,nnztot=0.;
509: MatInfo info;
510: PetscBool is_last = PETSC_FALSE;
513: PetscObjectGetComm((PetscObject)pc,&comm);
514: MPI_Comm_rank(comm,&rank);
515: MPI_Comm_size(comm,&size);
517: if (pc->setupcalled) {
518: if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) {
519: /* reset everything */
520: PCReset_MG(pc);
521: pc->setupcalled = 0;
522: } else {
523: PC_MG_Levels **mglevels = mg->levels;
524: /* just do Galerkin grids */
525: Mat B,dA,dB;
527: if (pc_gamg->Nlevels > 1) {
528: PetscInt gl;
529: /* currently only handle case where mat and pmat are the same on coarser levels */
530: KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);
531: /* (re)set to get dirty flag */
532: KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);
534: for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) {
535: MatReuse reuse = MAT_INITIAL_MATRIX ;
537: /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
538: KSPGetOperators(mglevels[level]->smoothd,NULL,&B);
539: if (B->product) {
540: if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) {
541: reuse = MAT_REUSE_MATRIX;
542: }
543: }
544: if (reuse == MAT_INITIAL_MATRIX) { MatDestroy(&mglevels[level]->A); }
545: if (reuse == MAT_REUSE_MATRIX) {
546: PetscInfo1(pc,"RAP after first solve, reuse matrix level %D\n",level);
547: } else {
548: PetscInfo1(pc,"RAP after first solve, new matrix level %D\n",level);
549: }
550: PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0);
551: MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B);
552: PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0);
553: if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B;
554: KSPSetOperators(mglevels[level]->smoothd,B,B);
555: dB = B;
556: }
557: }
559: PCSetUp_MG(pc);
560: return(0);
561: }
562: }
564: if (!pc_gamg->data) {
565: if (pc_gamg->orig_data) {
566: MatGetBlockSize(Pmat, &bs);
567: MatGetLocalSize(Pmat, &qq, NULL);
569: pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
570: pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
571: pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
573: PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);
574: for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
575: } else {
576: if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
577: pc_gamg->ops->createdefaultdata(pc,Pmat);
578: }
579: }
581: /* cache original data for reuse */
582: if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
583: PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);
584: for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
585: pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
586: pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
587: }
589: /* get basic dims */
590: MatGetBlockSize(Pmat, &bs);
591: MatGetSize(Pmat, &M, &N);
593: MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info); /* global reduction */
594: nnz0 = info.nz_used;
595: nnztot = info.nz_used;
596: PetscInfo6(pc,"level %D) N=%D, n data rows=%D, n data cols=%D, nnz/row (ave)=%d, np=%D\n",0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(int)(nnz0/(PetscReal)M+0.5),size);
598: /* Get A_i and R_i */
599: for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
600: pc_gamg->current_level = level;
601: if (level >= PETSC_MG_MAXLEVELS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level);
602: level1 = level + 1;
603: PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);
604: #if defined(GAMG_STAGES)
605: PetscLogStagePush(gamg_stages[level]);
606: #endif
607: { /* construct prolongator */
608: Mat Gmat;
609: PetscCoarsenData *agg_lists;
610: Mat Prol11;
612: pc_gamg->ops->graph(pc,Aarr[level], &Gmat);
613: pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);
614: pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);
616: /* could have failed to create new level */
617: if (Prol11) {
618: const char *prefix;
619: char addp[32];
621: /* get new block size of coarse matrices */
622: MatGetBlockSizes(Prol11, NULL, &bs);
624: if (pc_gamg->ops->optprolongator) {
625: /* smooth */
626: pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);
627: }
629: if (pc_gamg->use_aggs_in_asm) {
630: PetscInt bs;
631: MatGetBlockSizes(Prol11, &bs, NULL);
632: PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);
633: }
635: PCGetOptionsPrefix(pc,&prefix);
636: MatSetOptionsPrefix(Prol11,prefix);
637: PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level);
638: MatAppendOptionsPrefix(Prol11,addp);
639: /* Always generate the transpose with CUDA
640: Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
641: MatSetOption(Prol11,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);
642: MatSetFromOptions(Prol11);
643: Parr[level1] = Prol11;
644: } else Parr[level1] = NULL; /* failed to coarsen */
646: MatDestroy(&Gmat);
647: PetscCDDestroy(agg_lists);
648: } /* construct prolongator scope */
649: PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);
650: if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
651: if (!Parr[level1]) { /* failed to coarsen */
652: PetscInfo1(pc,"Stop gridding, level %D\n",level);
653: #if defined(GAMG_STAGES)
654: PetscLogStagePop();
655: #endif
656: break;
657: }
658: PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);
659: MatGetSize(Parr[level1], &M, &N); /* N is next M, a loop test variables */
660: if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?");
661: if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
662: if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
663: pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);
665: PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);
666: MatGetSize(Aarr[level1], &M, &N); /* M is loop test variables */
667: MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);
668: nnztot += info.nz_used;
669: PetscInfo5(pc,"%D) N=%D, n data cols=%D, nnz/row (ave)=%d, %D active pes\n",level1,M,pc_gamg->data_cell_cols,(int)(info.nz_used/(PetscReal)M),nactivepe);
671: #if defined(GAMG_STAGES)
672: PetscLogStagePop();
673: #endif
674: /* stop if one node or one proc -- could pull back for singular problems */
675: if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) {
676: PetscInfo2(pc,"HARD stop of coarsening on level %D. Grid too small: %D block nodes\n",level,M/bs);
677: level++;
678: break;
679: }
680: } /* levels */
681: PetscFree(pc_gamg->data);
683: PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);
684: pc_gamg->Nlevels = level + 1;
685: fine_level = level;
686: PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);
688: if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
689: PetscErrorCode (*savesetfromoptions[PETSC_MG_MAXLEVELS])(PetscOptionItems*,KSP);
691: /* set default smoothers & set operators */
692: for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
693: KSP smoother;
694: PC subpc;
696: PCMGGetSmoother(pc, lidx, &smoother);
697: KSPGetPC(smoother, &subpc);
699: KSPSetNormType(smoother, KSP_NORM_NONE);
700: /* set ops */
701: KSPSetOperators(smoother, Aarr[level], Aarr[level]);
702: PCMGSetInterpolation(pc, lidx, Parr[level+1]);
704: /* set defaults */
705: KSPSetType(smoother, KSPCHEBYSHEV);
707: /* set blocks for ASM smoother that uses the 'aggregates' */
708: if (pc_gamg->use_aggs_in_asm) {
709: PetscInt sz;
710: IS *iss;
712: sz = nASMBlocksArr[level];
713: iss = ASMLocalIDsArr[level];
714: PCSetType(subpc, PCASM);
715: PCASMSetOverlap(subpc, 0);
716: PCASMSetType(subpc,PC_ASM_BASIC);
717: if (!sz) {
718: IS is;
719: ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);
720: PCASMSetLocalSubdomains(subpc, 1, NULL, &is);
721: ISDestroy(&is);
722: } else {
723: PetscInt kk;
724: PCASMSetLocalSubdomains(subpc, sz, NULL, iss);
725: for (kk=0; kk<sz; kk++) {
726: ISDestroy(&iss[kk]);
727: }
728: PetscFree(iss);
729: }
730: ASMLocalIDsArr[level] = NULL;
731: nASMBlocksArr[level] = 0;
732: } else {
733: PCSetType(subpc, PCSOR);
734: }
735: }
736: {
737: /* coarse grid */
738: KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
739: Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
741: PCMGGetSmoother(pc, lidx, &smoother);
742: KSPSetOperators(smoother, Lmat, Lmat);
743: if (!pc_gamg->use_parallel_coarse_grid_solver) {
744: KSPSetNormType(smoother, KSP_NORM_NONE);
745: KSPGetPC(smoother, &subpc);
746: PCSetType(subpc, PCBJACOBI);
747: PCSetUp(subpc);
748: PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);
749: if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
750: KSPGetPC(k2[0],&pc2);
751: PCSetType(pc2, PCLU);
752: PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);
753: KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);
754: KSPSetType(k2[0], KSPPREONLY);
755: }
756: }
758: /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
759: PetscObjectOptionsBegin((PetscObject)pc);
760: PCSetFromOptions_MG(PetscOptionsObject,pc);
761: PetscOptionsEnd();
762: PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);
764: /* setup cheby eigen estimates from SA */
765: if (pc_gamg->use_sa_esteig==1) {
766: for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
767: KSP smoother;
768: PetscBool ischeb;
770: savesetfromoptions[level] = NULL;
771: PCMGGetSmoother(pc, lidx, &smoother);
772: PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb);
773: if (ischeb) {
774: KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data;
776: KSPSetFromOptions(smoother); /* let command line emax override using SA's eigenvalues */
777: if (mg->max_eigen_DinvA[level] > 0 && cheb->emax == 0.) {
778: PC subpc;
779: PetscBool isjac;
780: KSPGetPC(smoother, &subpc);
781: PetscObjectTypeCompare((PetscObject)subpc,PCJACOBI,&isjac);
782: if (isjac && pc_gamg->use_sa_esteig==1) {
783: PetscReal emax,emin;
785: emin = mg->min_eigen_DinvA[level];
786: emax = mg->max_eigen_DinvA[level];
787: PetscInfo4(pc,"PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %D (N=%D) with emax = %g emin = %g\n",level,Aarr[level]->rmap->N,(double)emax,(double)emin);
788: cheb->emin_computed = emin;
789: cheb->emax_computed = emax;
790: KSPChebyshevSetEigenvalues(smoother, cheb->tform[2]*emin + cheb->tform[3]*emax, cheb->tform[0]*emin + cheb->tform[1]*emax);
792: /* We have set the eigenvalues and consumed the transformation values
793: prevent from flagging the recomputation of the eigenvalues again in PCSetUp_MG
794: below when setfromoptions will be called again */
795: savesetfromoptions[level] = smoother->ops->setfromoptions;
796: smoother->ops->setfromoptions = NULL;
797: }
798: }
799: }
800: }
801: }
803: PCSetUp_MG(pc);
805: /* restore Chebyshev smoother for next calls */
806: if (pc_gamg->use_sa_esteig==1) {
807: for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
808: if (savesetfromoptions[level]) {
809: KSP smoother;
810: PCMGGetSmoother(pc, lidx, &smoother);
811: smoother->ops->setfromoptions = savesetfromoptions[level];
812: }
813: }
814: }
816: /* clean up */
817: for (level=1; level<pc_gamg->Nlevels; level++) {
818: MatDestroy(&Parr[level]);
819: MatDestroy(&Aarr[level]);
820: }
821: } else {
822: KSP smoother;
824: PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");
825: PCMGGetSmoother(pc, 0, &smoother);
826: KSPSetOperators(smoother, Aarr[0], Aarr[0]);
827: KSPSetType(smoother, KSPPREONLY);
828: PCSetUp_MG(pc);
829: }
830: return(0);
831: }
833: /* ------------------------------------------------------------------------- */
834: /*
835: PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
836: that was created with PCCreate_GAMG().
838: Input Parameter:
839: . pc - the preconditioner context
841: Application Interface Routine: PCDestroy()
842: */
843: PetscErrorCode PCDestroy_GAMG(PC pc)
844: {
846: PC_MG *mg = (PC_MG*)pc->data;
847: PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx;
850: PCReset_GAMG(pc);
851: if (pc_gamg->ops->destroy) {
852: (*pc_gamg->ops->destroy)(pc);
853: }
854: PetscFree(pc_gamg->ops);
855: PetscFree(pc_gamg->gamg_type_name);
856: PetscFree(pc_gamg);
857: PCDestroy_MG(pc);
858: return(0);
859: }
861: /*@
862: PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
864: Logically Collective on PC
866: Input Parameters:
867: + pc - the preconditioner context
868: - n - the number of equations
870: Options Database Key:
871: . -pc_gamg_process_eq_limit <limit>
873: Notes:
874: GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
875: that has degrees of freedom
877: Level: intermediate
879: .seealso: PCGAMGSetCoarseEqLim(), PCGAMGSetRankReductionFactors()
880: @*/
881: PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n)
882: {
887: PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));
888: return(0);
889: }
891: static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
892: {
893: PC_MG *mg = (PC_MG*)pc->data;
894: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
897: if (n>0) pc_gamg->min_eq_proc = n;
898: return(0);
899: }
901: /*@
902: PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
904: Collective on PC
906: Input Parameters:
907: + pc - the preconditioner context
908: - n - maximum number of equations to aim for
910: Options Database Key:
911: . -pc_gamg_coarse_eq_limit <limit>
913: Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
914: has less than 1000 unknowns.
916: Level: intermediate
918: .seealso: PCGAMGSetProcEqLim(), PCGAMGSetRankReductionFactors()
919: @*/
920: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
921: {
926: PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
927: return(0);
928: }
930: static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
931: {
932: PC_MG *mg = (PC_MG*)pc->data;
933: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
936: if (n>0) pc_gamg->coarse_eq_limit = n;
937: return(0);
938: }
940: /*@
941: PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
943: Collective on PC
945: Input Parameters:
946: + pc - the preconditioner context
947: - n - PETSC_TRUE or PETSC_FALSE
949: Options Database Key:
950: . -pc_gamg_repartition <true,false>
952: Notes:
953: this will generally improve the loading balancing of the work on each level
955: Level: intermediate
957: .seealso: ()
958: @*/
959: PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
960: {
965: PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));
966: return(0);
967: }
969: static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
970: {
971: PC_MG *mg = (PC_MG*)pc->data;
972: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
975: pc_gamg->repart = n;
976: return(0);
977: }
979: /*@
980: PCGAMGSetEstEigKSPMaxIt - Set number of KSP iterations in eigen estimator (for Cheby)
982: Collective on PC
984: Input Parameters:
985: + pc - the preconditioner context
986: - n - number of its
988: Options Database Key:
989: . -pc_gamg_esteig_ksp_max_it <its>
991: Notes:
993: Level: intermediate
995: .seealso: ()
996: @*/
997: PetscErrorCode PCGAMGSetEstEigKSPMaxIt(PC pc, PetscInt n)
998: {
1003: PetscTryMethod(pc,"PCGAMGSetEstEigKSPMaxIt_C",(PC,PetscInt),(pc,n));
1004: return(0);
1005: }
1007: static PetscErrorCode PCGAMGSetEstEigKSPMaxIt_GAMG(PC pc, PetscInt n)
1008: {
1009: PC_MG *mg = (PC_MG*)pc->data;
1010: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1013: pc_gamg->esteig_max_it = n;
1014: return(0);
1015: }
1017: /*@
1018: PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Cheby smoother
1020: Collective on PC
1022: Input Parameters:
1023: + pc - the preconditioner context
1024: - n - number of its
1026: Options Database Key:
1027: . -pc_gamg_use_sa_esteig <true,false>
1029: Notes:
1031: Level: intermediate
1033: .seealso: ()
1034: @*/
1035: PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
1036: {
1041: PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));
1042: return(0);
1043: }
1045: static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
1046: {
1047: PC_MG *mg = (PC_MG*)pc->data;
1048: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1051: pc_gamg->use_sa_esteig = n ? 1 : 0;
1052: return(0);
1053: }
1055: /*@C
1056: PCGAMGSetEstEigKSPType - Set type of KSP in eigen estimator (for Cheby)
1058: Collective on PC
1060: Input Parameters:
1061: + pc - the preconditioner context
1062: - t - ksp type
1064: Options Database Key:
1065: . -pc_gamg_esteig_ksp_type <type>
1067: Notes:
1069: Level: intermediate
1071: .seealso: ()
1072: @*/
1073: PetscErrorCode PCGAMGSetEstEigKSPType(PC pc, char t[])
1074: {
1079: PetscTryMethod(pc,"PCGAMGSetEstEigKSPType_C",(PC,char[]),(pc,t));
1080: return(0);
1081: }
1083: static PetscErrorCode PCGAMGSetEstEigKSPType_GAMG(PC pc, char t[])
1084: {
1086: PC_MG *mg = (PC_MG*)pc->data;
1087: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1090: PetscStrcpy(pc_gamg->esteig_type,t);
1091: return(0);
1092: }
1094: /*@
1095: PCGAMGSetEigenvalues - Set eigenvalues
1097: Collective on PC
1099: Input Parameters:
1100: + pc - the preconditioner context
1101: - emax - max eigenvalue
1102: . emin - min eigenvalue
1104: Options Database Key:
1105: . -pc_gamg_eigenvalues
1107: Level: intermediate
1109: .seealso: PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPType()
1110: @*/
1111: PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin)
1112: {
1117: PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));
1118: return(0);
1119: }
1121: static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin)
1122: {
1123: PC_MG *mg = (PC_MG*)pc->data;
1124: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1127: if (emax <= emin) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin);
1128: if (emax*emin <= 0.0) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin);
1129: pc_gamg->emax = emax;
1130: pc_gamg->emin = emin;
1132: return(0);
1133: }
1135: /*@
1136: PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
1138: Collective on PC
1140: Input Parameters:
1141: + pc - the preconditioner context
1142: - n - PETSC_TRUE or PETSC_FALSE
1144: Options Database Key:
1145: . -pc_gamg_reuse_interpolation <true,false>
1147: Level: intermediate
1149: Notes:
1150: this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1151: rebuilding the preconditioner quicker.
1153: .seealso: ()
1154: @*/
1155: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1156: {
1161: PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
1162: return(0);
1163: }
1165: static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1166: {
1167: PC_MG *mg = (PC_MG*)pc->data;
1168: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1171: pc_gamg->reuse_prol = n;
1172: return(0);
1173: }
1175: /*@
1176: PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner.
1178: Collective on PC
1180: Input Parameters:
1181: + pc - the preconditioner context
1182: - flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
1184: Options Database Key:
1185: . -pc_gamg_asm_use_agg
1187: Level: intermediate
1189: .seealso: ()
1190: @*/
1191: PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1192: {
1197: PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));
1198: return(0);
1199: }
1201: static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1202: {
1203: PC_MG *mg = (PC_MG*)pc->data;
1204: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1207: pc_gamg->use_aggs_in_asm = flg;
1208: return(0);
1209: }
1211: /*@
1212: PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
1214: Collective on PC
1216: Input Parameters:
1217: + pc - the preconditioner context
1218: - flg - PETSC_TRUE to not force coarse grid onto one processor
1220: Options Database Key:
1221: . -pc_gamg_use_parallel_coarse_grid_solver
1223: Level: intermediate
1225: .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids()
1226: @*/
1227: PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1228: {
1233: PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));
1234: return(0);
1235: }
1237: static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1238: {
1239: PC_MG *mg = (PC_MG*)pc->data;
1240: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1243: pc_gamg->use_parallel_coarse_grid_solver = flg;
1244: return(0);
1245: }
1247: /*@
1248: PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU
1250: Collective on PC
1252: Input Parameters:
1253: + pc - the preconditioner context
1254: - flg - PETSC_TRUE to pin coarse grids to CPU
1256: Options Database Key:
1257: . -pc_gamg_cpu_pin_coarse_grids
1259: Level: intermediate
1261: .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve()
1262: @*/
1263: PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1264: {
1269: PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));
1270: return(0);
1271: }
1273: static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1274: {
1275: PC_MG *mg = (PC_MG*)pc->data;
1276: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1279: pc_gamg->cpu_pin_coarse_grids = flg;
1280: return(0);
1281: }
1283: /*@
1284: PCGAMGSetCoarseGridLayoutType - place reduce grids on processors with natural order (compact type)
1286: Collective on PC
1288: Input Parameters:
1289: + pc - the preconditioner context
1290: - flg - Layout type
1292: Options Database Key:
1293: . -pc_gamg_coarse_grid_layout_type
1295: Level: intermediate
1297: .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids()
1298: @*/
1299: PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1300: {
1305: PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));
1306: return(0);
1307: }
1309: static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1310: {
1311: PC_MG *mg = (PC_MG*)pc->data;
1312: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1315: pc_gamg->layout_type = flg;
1316: return(0);
1317: }
1319: /*@
1320: PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use
1322: Not collective on PC
1324: Input Parameters:
1325: + pc - the preconditioner
1326: - n - the maximum number of levels to use
1328: Options Database Key:
1329: . -pc_mg_levels
1331: Level: intermediate
1333: .seealso: ()
1334: @*/
1335: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1336: {
1341: PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
1342: return(0);
1343: }
1345: static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1346: {
1347: PC_MG *mg = (PC_MG*)pc->data;
1348: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1351: pc_gamg->Nlevels = n;
1352: return(0);
1353: }
1355: /*@
1356: PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
1358: Not collective on PC
1360: Input Parameters:
1361: + pc - the preconditioner context
1362: . v - array of threshold values for finest n levels; 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
1363: - n - number of threshold values provided in array
1365: Options Database Key:
1366: . -pc_gamg_threshold <threshold>
1368: Notes:
1369: Increasing the threshold decreases the rate of coarsening. Conversely reducing the threshold increases the rate of coarsening (aggressive coarsening) and thereby reduces the complexity of the coarse grids, and generally results in slower solver converge rates. Reducing coarse grid complexity reduced the complexity of Galerkin coarse grid construction considerably.
1370: Before coarsening or aggregating the graph, GAMG removes small values from the graph with this threshold, and thus reducing the coupling in the graph and a different (perhaps better) coarser set of points.
1372: If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening.
1373: In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold().
1374: If n is greater than the total number of levels, the excess entries in threshold will not be used.
1376: Level: intermediate
1378: .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
1379: @*/
1380: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1381: {
1387: PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));
1388: return(0);
1389: }
1391: static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1392: {
1393: PC_MG *mg = (PC_MG*)pc->data;
1394: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1395: PetscInt i;
1397: for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1398: for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;
1399: return(0);
1400: }
1402: /*@
1403: PCGAMGSetRankReductionFactors - Set manual schedual for process reduction on coarse grids
1405: Collective on PC
1407: Input Parameters:
1408: + pc - the preconditioner context
1409: . v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda
1410: - n - number of values provided in array
1412: Options Database Key:
1413: . -pc_gamg_rank_reduction_factors <factors>
1415: Level: intermediate
1417: .seealso: PCGAMGSetProcEqLim(), PCGAMGSetCoarseEqLim()
1418: @*/
1419: PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1420: {
1426: PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n));
1427: return(0);
1428: }
1430: static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1431: {
1432: PC_MG *mg = (PC_MG*)pc->data;
1433: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1434: PetscInt i;
1436: for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1437: for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
1438: return(0);
1439: }
1441: /*@
1442: PCGAMGSetThresholdScale - Relative threshold reduction at each level
1444: Not collective on PC
1446: Input Parameters:
1447: + pc - the preconditioner context
1448: - scale - the threshold value reduction, ussually < 1.0
1450: Options Database Key:
1451: . -pc_gamg_threshold_scale <v>
1453: Notes:
1454: The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold().
1455: This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold().
1457: Level: advanced
1459: .seealso: PCGAMGSetThreshold()
1460: @*/
1461: PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1462: {
1467: PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));
1468: return(0);
1469: }
1471: static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1472: {
1473: PC_MG *mg = (PC_MG*)pc->data;
1474: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1476: pc_gamg->threshold_scale = v;
1477: return(0);
1478: }
1480: /*@C
1481: PCGAMGSetType - Set solution method
1483: Collective on PC
1485: Input Parameters:
1486: + pc - the preconditioner context
1487: - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1489: Options Database Key:
1490: . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1492: Level: intermediate
1494: .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1495: @*/
1496: PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1497: {
1502: PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1503: return(0);
1504: }
1506: /*@C
1507: PCGAMGGetType - Get solution method
1509: Collective on PC
1511: Input Parameter:
1512: . pc - the preconditioner context
1514: Output Parameter:
1515: . type - the type of algorithm used
1517: Level: intermediate
1519: .seealso: PCGAMGSetType(), PCGAMGType
1520: @*/
1521: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1522: {
1527: PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1528: return(0);
1529: }
1531: static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1532: {
1533: PC_MG *mg = (PC_MG*)pc->data;
1534: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1537: *type = pc_gamg->type;
1538: return(0);
1539: }
1541: static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1542: {
1543: PetscErrorCode ierr,(*r)(PC);
1544: PC_MG *mg = (PC_MG*)pc->data;
1545: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1548: pc_gamg->type = type;
1549: PetscFunctionListFind(GAMGList,type,&r);
1550: if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1551: if (pc_gamg->ops->destroy) {
1552: (*pc_gamg->ops->destroy)(pc);
1553: PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));
1554: pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1555: /* cleaning up common data in pc_gamg - this should disapear someday */
1556: pc_gamg->data_cell_cols = 0;
1557: pc_gamg->data_cell_rows = 0;
1558: pc_gamg->orig_data_cell_cols = 0;
1559: pc_gamg->orig_data_cell_rows = 0;
1560: PetscFree(pc_gamg->data);
1561: pc_gamg->data_sz = 0;
1562: }
1563: PetscFree(pc_gamg->gamg_type_name);
1564: PetscStrallocpy(type,&pc_gamg->gamg_type_name);
1565: (*r)(pc);
1566: return(0);
1567: }
1569: /* -------------------------------------------------------------------------- */
1570: /*
1571: PCMGGetGridComplexity - compute coarse grid complexity of MG hierarchy
1573: Input Parameter:
1574: . pc - the preconditioner context
1576: Output Parameter:
1577: . gc - grid complexity = sum_i(nnz_i) / nnz_0
1579: Level: advanced
1580: */
1581: static PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc)
1582: {
1584: PC_MG *mg = (PC_MG*)pc->data;
1585: PC_MG_Levels **mglevels = mg->levels;
1586: PetscInt lev;
1587: PetscLogDouble nnz0 = 0, sgc = 0;
1588: MatInfo info;
1591: if (!pc->setupcalled) {
1592: *gc = 0;
1593: return(0);
1594: }
1595: if (!mg->nlevels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels");
1596: for (lev=0; lev<mg->nlevels; lev++) {
1597: Mat dB;
1598: KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB);
1599: MatGetInfo(dB,MAT_GLOBAL_SUM,&info); /* global reduction */
1600: sgc += info.nz_used;
1601: if (lev==mg->nlevels-1) nnz0 = info.nz_used;
1602: }
1603: if (nnz0 > 0) *gc = (PetscReal)(sgc/nnz0);
1604: else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number for grid points on finest level is not available");
1605: return(0);
1606: }
1608: static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
1609: {
1610: PetscErrorCode ierr,i;
1611: PC_MG *mg = (PC_MG*)pc->data;
1612: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1613: PetscReal gc=0;
1615: PetscViewerASCIIPrintf(viewer," GAMG specific options\n");
1616: PetscViewerASCIIPrintf(viewer," Threshold for dropping small values in graph on each level =");
1617: for (i=0;i<mg->nlevels; i++) {
1618: PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);
1619: }
1620: PetscViewerASCIIPrintf(viewer,"\n");
1621: PetscViewerASCIIPrintf(viewer," Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);
1622: if (pc_gamg->use_aggs_in_asm) {
1623: PetscViewerASCIIPrintf(viewer," Using aggregates from coarsening process to define subdomains for PCASM\n");
1624: }
1625: if (pc_gamg->use_parallel_coarse_grid_solver) {
1626: PetscViewerASCIIPrintf(viewer," Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");
1627: }
1628: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1629: if (pc_gamg->cpu_pin_coarse_grids) {
1630: /* PetscViewerASCIIPrintf(viewer," Pinning coarse grids to the CPU)\n"); */
1631: }
1632: #endif
1633: /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
1634: /* PetscViewerASCIIPrintf(viewer," Put reduced grids on processes in natural order (ie, 0,1,2...)\n"); */
1635: /* } else { */
1636: /* PetscViewerASCIIPrintf(viewer," Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n"); */
1637: /* } */
1638: if (pc_gamg->ops->view) {
1639: (*pc_gamg->ops->view)(pc,viewer);
1640: }
1641: PCMGGetGridComplexity(pc,&gc);
1642: PetscViewerASCIIPrintf(viewer," Complexity: grid = %g\n",gc);
1643: return(0);
1644: }
1646: PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
1647: {
1649: PC_MG *mg = (PC_MG*)pc->data;
1650: PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1651: PetscBool flag,f2;
1652: MPI_Comm comm;
1653: char prefix[256],tname[32];
1654: PetscInt i,n;
1655: const char *pcpre;
1656: static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL};
1658: PetscObjectGetComm((PetscObject)pc,&comm);
1659: PetscOptionsHead(PetscOptionsObject,"GAMG options");
1660: PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);
1661: if (flag) {
1662: PCGAMGSetType(pc,tname);
1663: }
1664: PetscOptionsFList("-pc_gamg_esteig_ksp_type","Krylov method for eigen estimator","PCGAMGSetEstEigKSPType",KSPList,pc_gamg->esteig_type,tname,sizeof(tname),&flag);
1665: if (flag) {
1666: PCGAMGSetEstEigKSPType(pc,tname);
1667: }
1668: PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);
1669: f2 = PETSC_TRUE;
1670: PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",f2,&f2,&flag);
1671: if (flag) pc_gamg->use_sa_esteig = f2 ? 1 : 0;
1672: PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);
1673: PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation aggregates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL);
1674: PetscOptionsBool("-pc_gamg_use_parallel_coarse_grid_solver","Use parallel coarse grid solver (otherwise put last grid on one process)","PCGAMGSetUseParallelCoarseGridSolve",pc_gamg->use_parallel_coarse_grid_solver,&pc_gamg->use_parallel_coarse_grid_solver,NULL);
1675: PetscOptionsBool("-pc_gamg_cpu_pin_coarse_grids","Pin coarse grids to the CPU","PCGAMGSetCpuPinCoarseGrids",pc_gamg->cpu_pin_coarse_grids,&pc_gamg->cpu_pin_coarse_grids,NULL);
1676: PetscOptionsEnum("-pc_gamg_coarse_grid_layout_type","compact: place reduced grids on processes in natural order; spread: distribute to whole machine for more memory bandwidth","PCGAMGSetCoarseGridLayoutType",LayoutTypes,(PetscEnum)pc_gamg->layout_type,(PetscEnum*)&pc_gamg->layout_type,NULL);
1677: PetscOptionsInt("-pc_gamg_process_eq_limit","Limit (goal) on number of equations per process on coarse grids","PCGAMGSetProcEqLim",pc_gamg->min_eq_proc,&pc_gamg->min_eq_proc,NULL);
1678: PetscOptionsInt("-pc_gamg_esteig_ksp_max_it","Number of iterations of eigen estimator","PCGAMGSetEstEigKSPMaxIt",pc_gamg->esteig_max_it,&pc_gamg->esteig_max_it,NULL);
1679: PetscOptionsInt("-pc_gamg_coarse_eq_limit","Limit on number of equations for the coarse grid","PCGAMGSetCoarseEqLim",pc_gamg->coarse_eq_limit,&pc_gamg->coarse_eq_limit,NULL);
1680: PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL);
1681: n = PETSC_MG_MAXLEVELS;
1682: PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);
1683: if (!flag || n < PETSC_MG_MAXLEVELS) {
1684: if (!flag) n = 1;
1685: i = n;
1686: do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS);
1687: }
1688: n = PETSC_MG_MAXLEVELS;
1689: PetscOptionsIntArray("-pc_gamg_rank_reduction_factors","Manual schedule of coarse grid reduction factors that overrides internal heuristics (0 for first reduction puts one process/device)","PCGAMGSetRankReductionFactors",pc_gamg->level_reduction_factors,&n,&flag);
1690: if (!flag) i = 0;
1691: else i = n;
1692: do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS);
1693: PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);
1694: {
1695: PetscReal eminmax[2] = {0., 0.};
1696: n = 2;
1697: PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag);
1698: if (flag) {
1699: if (n != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
1700: PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);
1701: }
1702: }
1703: /* set options for subtype */
1704: if (pc_gamg->ops->setfromoptions) {(*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);}
1706: PCGetOptionsPrefix(pc, &pcpre);
1707: PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
1708: PetscOptionsTail();
1709: return(0);
1710: }
1712: /* -------------------------------------------------------------------------- */
1713: /*MC
1714: PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
1716: Options Database Keys:
1717: + -pc_gamg_type <type> - one of agg, geo, or classical
1718: . -pc_gamg_repartition <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1719: . -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1720: . -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother
1721: . -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit>
1722: equations on each process that has degrees of freedom
1723: . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
1724: . -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1725: - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1727: Options Database Keys for default Aggregation:
1728: + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1729: . -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1730: - -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1732: Multigrid options:
1733: + -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1734: . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1735: . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1736: - -pc_mg_levels <levels> - Number of levels of multigrid to use.
1738: Notes:
1739: In order to obtain good performance for PCGAMG for vector valued problems you must
1740: Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1741: Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1742: See the Users Manual Chapter 4 for more details
1744: Level: intermediate
1746: .seealso: PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1747: PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetEstEigKSPType()
1748: M*/
1750: PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1751: {
1752: PetscErrorCode ierr,i;
1753: PC_GAMG *pc_gamg;
1754: PC_MG *mg;
1757: /* register AMG type */
1758: PCGAMGInitializePackage();
1760: /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1761: PCSetType(pc, PCMG);
1762: PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);
1764: /* create a supporting struct and attach it to pc */
1765: PetscNewLog(pc,&pc_gamg);
1766: PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);
1767: mg = (PC_MG*)pc->data;
1768: mg->innerctx = pc_gamg;
1770: PetscNewLog(pc,&pc_gamg->ops);
1772: /* these should be in subctx but repartitioning needs simple arrays */
1773: pc_gamg->data_sz = 0;
1774: pc_gamg->data = NULL;
1776: /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1777: pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1778: pc->ops->setup = PCSetUp_GAMG;
1779: pc->ops->reset = PCReset_GAMG;
1780: pc->ops->destroy = PCDestroy_GAMG;
1781: mg->view = PCView_GAMG;
1783: PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);
1784: PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);
1785: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);
1786: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);
1787: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);
1788: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPType_C",PCGAMGSetEstEigKSPType_GAMG);
1789: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPMaxIt_C",PCGAMGSetEstEigKSPMaxIt_GAMG);
1790: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG);
1791: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG);
1792: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);
1793: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);
1794: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);
1795: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);
1796: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);
1797: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);
1798: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG);
1799: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);
1800: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);
1801: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);
1802: PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);
1803: pc_gamg->repart = PETSC_FALSE;
1804: pc_gamg->reuse_prol = PETSC_FALSE;
1805: pc_gamg->use_aggs_in_asm = PETSC_FALSE;
1806: pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1807: pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1808: pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD;
1809: pc_gamg->min_eq_proc = 50;
1810: pc_gamg->coarse_eq_limit = 50;
1811: for (i=0;i<PETSC_MG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1812: pc_gamg->threshold_scale = 1.;
1813: pc_gamg->Nlevels = PETSC_MG_MAXLEVELS;
1814: pc_gamg->current_level = 0; /* don't need to init really */
1815: PetscStrcpy(pc_gamg->esteig_type,NULL);
1816: pc_gamg->esteig_max_it = 10;
1817: pc_gamg->use_sa_esteig = -1;
1818: pc_gamg->emin = 0;
1819: pc_gamg->emax = 0;
1821: pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1823: /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1824: PCGAMGSetType(pc,PCGAMGAGG);
1825: return(0);
1826: }
1828: /*@C
1829: PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
1830: from PCInitializePackage().
1832: Level: developer
1834: .seealso: PetscInitialize()
1835: @*/
1836: PetscErrorCode PCGAMGInitializePackage(void)
1837: {
1839: PetscInt l;
1842: if (PCGAMGPackageInitialized) return(0);
1843: PCGAMGPackageInitialized = PETSC_TRUE;
1844: PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);
1845: PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);
1846: PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);
1847: PetscRegisterFinalize(PCGAMGFinalizePackage);
1849: /* general events */
1850: PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);
1851: PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);
1852: PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1853: PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1854: PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1855: PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1856: PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);
1858: PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1859: PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1860: /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1861: /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1862: /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1863: PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1864: PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1865: PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1866: PetscLogEventRegister(" search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1867: PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1868: PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1869: PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1870: PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1871: PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1872: PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1873: PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1874: PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);
1875: for (l=0;l<PETSC_MG_MAXLEVELS;l++) {
1876: char ename[32];
1878: PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02d",l);
1879: PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]);
1880: PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02d",l);
1881: PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]);
1882: PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02d",l);
1883: PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]);
1884: }
1885: /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1886: /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1887: /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1888: /* create timer stages */
1889: #if defined(GAMG_STAGES)
1890: {
1891: char str[32];
1892: PetscInt lidx;
1893: sprintf(str,"MG Level %d (finest)",0);
1894: PetscLogStageRegister(str, &gamg_stages[0]);
1895: for (lidx=1; lidx<9; lidx++) {
1896: sprintf(str,"MG Level %d",(int)lidx);
1897: PetscLogStageRegister(str, &gamg_stages[lidx]);
1898: }
1899: }
1900: #endif
1901: return(0);
1902: }
1904: /*@C
1905: PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1906: called from PetscFinalize() automatically.
1908: Level: developer
1910: .seealso: PetscFinalize()
1911: @*/
1912: PetscErrorCode PCGAMGFinalizePackage(void)
1913: {
1917: PCGAMGPackageInitialized = PETSC_FALSE;
1918: PetscFunctionListDestroy(&GAMGList);
1919: return(0);
1920: }
1922: /*@C
1923: PCGAMGRegister - Register a PCGAMG implementation.
1925: Input Parameters:
1926: + type - string that will be used as the name of the GAMG type.
1927: - create - function for creating the gamg context.
1929: Level: advanced
1931: .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1932: @*/
1933: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1934: {
1938: PCGAMGInitializePackage();
1939: PetscFunctionListAdd(&GAMGList,type,create);
1940: return(0);
1941: }