Actual source code: util.c
petsc-3.11.4 2019-09-28
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>
7: /*
8: Produces a set of block column indices of the matrix row, one for each block represented in the original row
10: n - the number of block indices in cc[]
11: cc - the block indices (must be large enough to contain the indices)
12: */
13: PETSC_STATIC_INLINE PetscErrorCode MatCollapseRow(Mat Amat,PetscInt row,PetscInt bs,PetscInt *n,PetscInt *cc)
14: {
15: PetscInt cnt = -1,nidx,j;
16: const PetscInt *idx;
20: MatGetRow(Amat,row,&nidx,&idx,NULL);
21: if (nidx) {
22: cnt = 0;
23: cc[cnt] = idx[0]/bs;
24: for (j=1; j<nidx; j++) {
25: if (cc[cnt] < idx[j]/bs) cc[++cnt] = idx[j]/bs;
26: }
27: }
28: MatRestoreRow(Amat,row,&nidx,&idx,NULL);
29: *n = cnt+1;
30: return(0);
31: }
33: /*
34: Produces a set of block column indices of the matrix block row, one for each block represented in the original set of rows
36: ncollapsed - the number of block indices
37: collapsed - the block indices (must be large enough to contain the indices)
38: */
39: PETSC_STATIC_INLINE PetscErrorCode MatCollapseRows(Mat Amat,PetscInt start,PetscInt bs,PetscInt *w0,PetscInt *w1,PetscInt *w2,PetscInt *ncollapsed,PetscInt **collapsed)
40: {
41: PetscInt i,nprev,*cprev = w0,ncur = 0,*ccur = w1,*merged = w2,*cprevtmp;
45: MatCollapseRow(Amat,start,bs,&nprev,cprev);
46: for (i=start+1; i<start+bs; i++) {
47: MatCollapseRow(Amat,i,bs,&ncur,ccur);
48: PetscMergeIntArray(nprev,cprev,ncur,ccur,&nprev,&merged);
49: cprevtmp = cprev; cprev = merged; merged = cprevtmp;
50: }
51: *ncollapsed = nprev;
52: if (collapsed) *collapsed = cprev;
53: return(0);
54: }
57: /* -------------------------------------------------------------------------- */
58: /*
59: PCGAMGCreateGraph - create simple scaled scalar graph from matrix
61: Input Parameter:
62: . Amat - matrix
63: Output Parameter:
64: . a_Gmaat - eoutput scalar graph (symmetric?)
65: */
66: PetscErrorCode PCGAMGCreateGraph(Mat Amat, Mat *a_Gmat)
67: {
69: PetscInt Istart,Iend,Ii,jj,kk,ncols,nloc,NN,MM,bs;
70: MPI_Comm comm;
71: Mat Gmat;
72: MatType mtype;
75: PetscObjectGetComm((PetscObject)Amat,&comm);
76: MatGetOwnershipRange(Amat, &Istart, &Iend);
77: MatGetSize(Amat, &MM, &NN);
78: MatGetBlockSize(Amat, &bs);
79: nloc = (Iend-Istart)/bs;
81: #if defined PETSC_GAMG_USE_LOG
82: PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);
83: #endif
85: if (bs > 1) {
86: const PetscScalar *vals;
87: const PetscInt *idx;
88: PetscInt *d_nnz, *o_nnz,*w0,*w1,*w2;
89: PetscBool ismpiaij,isseqaij;
91: /*
92: Determine the preallocation needed for the scalar matrix derived from the vector matrix.
93: */
95: PetscObjectBaseTypeCompare((PetscObject)Amat,MATSEQAIJ,&isseqaij);
96: PetscObjectBaseTypeCompare((PetscObject)Amat,MATMPIAIJ,&ismpiaij);
97: PetscMalloc2(nloc, &d_nnz,isseqaij ? 0 : nloc, &o_nnz);
99: if (isseqaij) {
100: PetscInt max_d_nnz;
102: /*
103: Determine exact preallocation count for (sequential) scalar matrix
104: */
105: MatSeqAIJGetMaxRowNonzeros(Amat,&max_d_nnz);
106: max_d_nnz = PetscMin(nloc,bs*max_d_nnz);
107: PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2);
108: for (Ii = 0, jj = 0; Ii < Iend; Ii += bs, jj++) {
109: MatCollapseRows(Amat,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL);
110: }
111: PetscFree3(w0,w1,w2);
113: } else if (ismpiaij) {
114: Mat Daij,Oaij;
115: const PetscInt *garray;
116: PetscInt max_d_nnz;
118: MatMPIAIJGetSeqAIJ(Amat,&Daij,&Oaij,&garray);
120: /*
121: Determine exact preallocation count for diagonal block portion of scalar matrix
122: */
123: MatSeqAIJGetMaxRowNonzeros(Daij,&max_d_nnz);
124: max_d_nnz = PetscMin(nloc,bs*max_d_nnz);
125: PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2);
126: for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) {
127: MatCollapseRows(Daij,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL);
128: }
129: PetscFree3(w0,w1,w2);
131: /*
132: Over estimate (usually grossly over), preallocation count for off-diagonal portion of scalar matrix
133: */
134: for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) {
135: o_nnz[jj] = 0;
136: for (kk=0; kk<bs; kk++) { /* rows that get collapsed to a single row */
137: MatGetRow(Oaij,Ii+kk,&ncols,0,0);
138: o_nnz[jj] += ncols;
139: MatRestoreRow(Oaij,Ii+kk,&ncols,0,0);
140: }
141: if (o_nnz[jj] > (NN/bs-nloc)) o_nnz[jj] = NN/bs-nloc;
142: }
144: } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require AIJ matrix type");
146: /* get scalar copy (norms) of matrix */
147: MatGetType(Amat,&mtype);
148: MatCreate(comm, &Gmat);
149: MatSetSizes(Gmat,nloc,nloc,PETSC_DETERMINE,PETSC_DETERMINE);
150: MatSetBlockSizes(Gmat, 1, 1);
151: MatSetType(Gmat, mtype);
152: MatSeqAIJSetPreallocation(Gmat,0,d_nnz);
153: MatMPIAIJSetPreallocation(Gmat,0,d_nnz,0,o_nnz);
154: PetscFree2(d_nnz,o_nnz);
156: for (Ii = Istart; Ii < Iend; Ii++) {
157: PetscInt dest_row = Ii/bs;
158: MatGetRow(Amat,Ii,&ncols,&idx,&vals);
159: for (jj=0; jj<ncols; jj++) {
160: PetscInt dest_col = idx[jj]/bs;
161: PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
162: MatSetValues(Gmat,1,&dest_row,1,&dest_col,&sv,ADD_VALUES);
163: }
164: MatRestoreRow(Amat,Ii,&ncols,&idx,&vals);
165: }
166: MatAssemblyBegin(Gmat,MAT_FINAL_ASSEMBLY);
167: MatAssemblyEnd(Gmat,MAT_FINAL_ASSEMBLY);
168: } else {
169: /* just copy scalar matrix - abs() not taken here but scaled later */
170: MatDuplicate(Amat, MAT_COPY_VALUES, &Gmat);
171: }
173: #if defined PETSC_GAMG_USE_LOG
174: PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);
175: #endif
177: *a_Gmat = Gmat;
178: return(0);
179: }
181: /* -------------------------------------------------------------------------- */
182: /*@C
183: PCGAMGFilterGraph - filter (remove zero and possibly small values from the) graph and make it symmetric if requested
185: Collective on Mat
187: Input Parameter:
188: + a_Gmat - the graph
189: . vfilter - threshold paramter [0,1)
190: - symm - make the result symmetric
192: Level: developer
194: Notes:
195: This is called before graph coarsers are called.
197: .seealso: PCGAMGSetThreshold()
198: @*/
199: PetscErrorCode PCGAMGFilterGraph(Mat *a_Gmat,PetscReal vfilter,PetscBool symm)
200: {
201: PetscErrorCode ierr;
202: PetscInt Istart,Iend,Ii,jj,ncols,nnz0,nnz1, NN, MM, nloc;
203: PetscMPIInt rank;
204: Mat Gmat = *a_Gmat, tGmat, matTrans;
205: MPI_Comm comm;
206: const PetscScalar *vals;
207: const PetscInt *idx;
208: PetscInt *d_nnz, *o_nnz;
209: Vec diag;
210: MatType mtype;
213: #if defined PETSC_GAMG_USE_LOG
214: PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);
215: #endif
216: /* scale Gmat for all values between -1 and 1 */
217: MatCreateVecs(Gmat, &diag, 0);
218: MatGetDiagonal(Gmat, diag);
219: VecReciprocal(diag);
220: VecSqrtAbs(diag);
221: MatDiagonalScale(Gmat, diag, diag);
222: VecDestroy(&diag);
224: if (vfilter < 0.0 && !symm) {
225: /* Just use the provided matrix as the graph but make all values positive */
226: MatInfo info;
227: PetscScalar *avals;
228: PetscBool isaij,ismpiaij;
229: PetscObjectBaseTypeCompare((PetscObject)Gmat,MATSEQAIJ,&isaij);
230: PetscObjectBaseTypeCompare((PetscObject)Gmat,MATMPIAIJ,&ismpiaij);
231: if (!isaij && !ismpiaij) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require (MPI)AIJ matrix type");
232: if (isaij) {
233: MatGetInfo(Gmat,MAT_LOCAL,&info);
234: MatSeqAIJGetArray(Gmat,&avals);
235: for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
236: MatSeqAIJRestoreArray(Gmat,&avals);
237: } else {
238: Mat_MPIAIJ *aij = (Mat_MPIAIJ*)Gmat->data;
239: MatGetInfo(aij->A,MAT_LOCAL,&info);
240: MatSeqAIJGetArray(aij->A,&avals);
241: for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
242: MatSeqAIJRestoreArray(aij->A,&avals);
243: MatGetInfo(aij->B,MAT_LOCAL,&info);
244: MatSeqAIJGetArray(aij->B,&avals);
245: for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
246: MatSeqAIJRestoreArray(aij->B,&avals);
247: }
248: #if defined PETSC_GAMG_USE_LOG
249: PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);
250: #endif
251: return(0);
252: }
254: PetscObjectGetComm((PetscObject)Gmat,&comm);
255: MPI_Comm_rank(comm,&rank);
256: MatGetOwnershipRange(Gmat, &Istart, &Iend);
257: nloc = Iend - Istart;
258: MatGetSize(Gmat, &MM, &NN);
260: if (symm) {
261: MatTranspose(Gmat, MAT_INITIAL_MATRIX, &matTrans);
262: }
264: /* Determine upper bound on nonzeros needed in new filtered matrix */
265: PetscMalloc2(nloc, &d_nnz,nloc, &o_nnz);
266: for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
267: MatGetRow(Gmat,Ii,&ncols,NULL,NULL);
268: d_nnz[jj] = ncols;
269: o_nnz[jj] = ncols;
270: MatRestoreRow(Gmat,Ii,&ncols,NULL,NULL);
271: if (symm) {
272: MatGetRow(matTrans,Ii,&ncols,NULL,NULL);
273: d_nnz[jj] += ncols;
274: o_nnz[jj] += ncols;
275: MatRestoreRow(matTrans,Ii,&ncols,NULL,NULL);
276: }
277: if (d_nnz[jj] > nloc) d_nnz[jj] = nloc;
278: if (o_nnz[jj] > (MM-nloc)) o_nnz[jj] = MM - nloc;
279: }
280: MatGetType(Gmat,&mtype);
281: MatCreate(comm, &tGmat);
282: MatSetSizes(tGmat,nloc,nloc,MM,MM);
283: MatSetBlockSizes(tGmat, 1, 1);
284: MatSetType(tGmat, mtype);
285: MatSeqAIJSetPreallocation(tGmat,0,d_nnz);
286: MatMPIAIJSetPreallocation(tGmat,0,d_nnz,0,o_nnz);
287: PetscFree2(d_nnz,o_nnz);
288: if (symm) {
289: MatDestroy(&matTrans);
290: } else {
291: /* all entries are generated locally so MatAssembly will be slightly faster for large process counts */
292: MatSetOption(tGmat,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
293: }
295: for (Ii = Istart, nnz0 = nnz1 = 0; Ii < Iend; Ii++) {
296: MatGetRow(Gmat,Ii,&ncols,&idx,&vals);
297: for (jj=0; jj<ncols; jj++,nnz0++) {
298: PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
299: if (PetscRealPart(sv) > vfilter) {
300: nnz1++;
301: if (symm) {
302: sv *= 0.5;
303: MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);
304: MatSetValues(tGmat,1,&idx[jj],1,&Ii,&sv,ADD_VALUES);
305: } else {
306: MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);
307: }
308: }
309: }
310: MatRestoreRow(Gmat,Ii,&ncols,&idx,&vals);
311: }
312: MatAssemblyBegin(tGmat,MAT_FINAL_ASSEMBLY);
313: MatAssemblyEnd(tGmat,MAT_FINAL_ASSEMBLY);
315: #if defined PETSC_GAMG_USE_LOG
316: PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);
317: #endif
319: #if defined(PETSC_USE_INFO)
320: {
321: double t1 = (!nnz0) ? 1. : 100.*(double)nnz1/(double)nnz0, t2 = (!nloc) ? 1. : (double)nnz0/(double)nloc;
322: PetscInfo4(*a_Gmat,"\t %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%D)\n",t1,vfilter,t2,MM);
323: }
324: #endif
325: MatDestroy(&Gmat);
326: *a_Gmat = tGmat;
327: return(0);
328: }
330: /* -------------------------------------------------------------------------- */
331: /*
332: PCGAMGGetDataWithGhosts - hacks into Mat MPIAIJ so this must have size > 1
334: Input Parameter:
335: . Gmat - MPIAIJ matrix for scattters
336: . data_sz - number of data terms per node (# cols in output)
337: . data_in[nloc*data_sz] - column oriented data
338: Output Parameter:
339: . a_stride - numbrt of rows of output
340: . a_data_out[stride*data_sz] - output data with ghosts
341: */
342: PetscErrorCode PCGAMGGetDataWithGhosts(Mat Gmat,PetscInt data_sz,PetscReal data_in[],PetscInt *a_stride,PetscReal **a_data_out)
343: {
345: Vec tmp_crds;
346: Mat_MPIAIJ *mpimat = (Mat_MPIAIJ*)Gmat->data;
347: PetscInt nnodes,num_ghosts,dir,kk,jj,my0,Iend,nloc;
348: PetscScalar *data_arr;
349: PetscReal *datas;
350: PetscBool isMPIAIJ;
353: PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPIAIJ);
354: MatGetOwnershipRange(Gmat, &my0, &Iend);
355: nloc = Iend - my0;
356: VecGetLocalSize(mpimat->lvec, &num_ghosts);
357: nnodes = num_ghosts + nloc;
358: *a_stride = nnodes;
359: MatCreateVecs(Gmat, &tmp_crds, 0);
361: PetscMalloc1(data_sz*nnodes, &datas);
362: for (dir=0; dir<data_sz; dir++) {
363: /* set local, and global */
364: for (kk=0; kk<nloc; kk++) {
365: PetscInt gid = my0 + kk;
366: PetscScalar crd = (PetscScalar)data_in[dir*nloc + kk]; /* col oriented */
367: datas[dir*nnodes + kk] = PetscRealPart(crd);
369: VecSetValues(tmp_crds, 1, &gid, &crd, INSERT_VALUES);
370: }
371: VecAssemblyBegin(tmp_crds);
372: VecAssemblyEnd(tmp_crds);
373: /* get ghost datas */
374: VecScatterBegin(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);
375: VecScatterEnd(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);
376: VecGetArray(mpimat->lvec, &data_arr);
377: for (kk=nloc,jj=0;jj<num_ghosts;kk++,jj++) datas[dir*nnodes + kk] = PetscRealPart(data_arr[jj]);
378: VecRestoreArray(mpimat->lvec, &data_arr);
379: }
380: VecDestroy(&tmp_crds);
381: *a_data_out = datas;
382: return(0);
383: }
385: PetscErrorCode PCGAMGHashTableCreate(PetscInt a_size, PCGAMGHashTable *a_tab)
386: {
388: PetscInt kk;
391: a_tab->size = a_size;
392: PetscMalloc2(a_size, &a_tab->table,a_size, &a_tab->data);
393: for (kk=0; kk<a_size; kk++) a_tab->table[kk] = -1;
394: return(0);
395: }
397: PetscErrorCode PCGAMGHashTableDestroy(PCGAMGHashTable *a_tab)
398: {
402: PetscFree2(a_tab->table,a_tab->data);
403: return(0);
404: }
406: PetscErrorCode PCGAMGHashTableAdd(PCGAMGHashTable *a_tab, PetscInt a_key, PetscInt a_data)
407: {
408: PetscInt kk,idx;
411: if (a_key<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Negative key %D.",a_key);
412: for (kk = 0, idx = GAMG_HASH(a_key); kk < a_tab->size; kk++, idx = (idx==(a_tab->size-1)) ? 0 : idx + 1) {
413: if (a_tab->table[idx] == a_key) {
414: /* exists */
415: a_tab->data[idx] = a_data;
416: break;
417: } else if (a_tab->table[idx] == -1) {
418: /* add */
419: a_tab->table[idx] = a_key;
420: a_tab->data[idx] = a_data;
421: break;
422: }
423: }
424: if (kk==a_tab->size) {
425: /* this is not to efficient, waiting until completely full */
426: PetscInt oldsize = a_tab->size, new_size = 2*a_tab->size + 5, *oldtable = a_tab->table, *olddata = a_tab->data;
429: a_tab->size = new_size;
430: PetscMalloc2(a_tab->size, &a_tab->table,a_tab->size, &a_tab->data);
431: for (kk=0;kk<a_tab->size;kk++) a_tab->table[kk] = -1;
432: for (kk=0;kk<oldsize;kk++) {
433: if (oldtable[kk] != -1) {
434: PCGAMGHashTableAdd(a_tab, oldtable[kk], olddata[kk]);
435: }
436: }
437: PetscFree2(oldtable,olddata);
438: PCGAMGHashTableAdd(a_tab, a_key, a_data);
439: }
440: return(0);
441: }