Actual source code: tools.c
petsc-3.4.2 2013-07-02
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>
8: /* -------------------------------------------------------------------------- */
9: /*
10: PCGAMGCreateGraph - create simple scaled scalar graph from matrix
12: Input Parameter:
13: . Amat - matrix
14: Output Parameter:
15: . a_Gmaat - eoutput scalar graph (symmetric?)
16: */
19: PetscErrorCode PCGAMGCreateGraph(const Mat Amat, Mat *a_Gmat)
20: {
22: PetscInt Istart,Iend,Ii,jj,kk,ncols,nloc,NN,MM,bs;
23: PetscMPIInt rank, size;
24: MPI_Comm comm;
25: Mat Gmat;
28: PetscObjectGetComm((PetscObject)Amat,&comm);
29: MPI_Comm_rank(comm,&rank);
30: MPI_Comm_size(comm,&size);
31: MatGetOwnershipRange(Amat, &Istart, &Iend);
32: MatGetSize(Amat, &MM, &NN);
33: MatGetBlockSize(Amat, &bs);
34: nloc = (Iend-Istart)/bs;
36: #if defined PETSC_GAMG_USE_LOG
37: PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);
38: #endif
39: if (bs > 1) {
40: const PetscScalar *vals;
41: const PetscInt *idx;
42: PetscInt *d_nnz, *o_nnz;
43: /* count nnz, there is sparcity in here so this might not be enough */
44: PetscMalloc(nloc*sizeof(PetscInt), &d_nnz);
45: PetscMalloc(nloc*sizeof(PetscInt), &o_nnz);
46: for (Ii = Istart, jj = 0; Ii < Iend; Ii += bs, jj++) {
47: d_nnz[jj] = 0;
48: for (kk=0; kk<bs; kk++) {
49: MatGetRow(Amat,Ii+kk,&ncols,0,0);
50: if (ncols > d_nnz[jj]) {
51: d_nnz[jj] = ncols; /* very pessimistic but could be too low in theory */
52: o_nnz[jj] = ncols;
53: if (d_nnz[jj] > nloc) d_nnz[jj] = nloc;
54: if (o_nnz[jj] > (NN/bs-nloc)) o_nnz[jj] = NN/bs-nloc;
55: }
56: MatRestoreRow(Amat,Ii+kk,&ncols,0,0);
57: }
58: }
60: /* get scalar copy (norms) of matrix -- AIJ specific!!! */
61: MatCreateAIJ(comm, nloc, nloc, PETSC_DETERMINE, PETSC_DETERMINE,0, d_nnz, 0, o_nnz, &Gmat);
63: PetscFree(d_nnz);
64: PetscFree(o_nnz);
66: for (Ii = Istart; Ii < Iend; Ii++) {
67: PetscInt dest_row = Ii/bs;
68: MatGetRow(Amat,Ii,&ncols,&idx,&vals);
69: for (jj=0; jj<ncols; jj++) {
70: PetscInt dest_col = idx[jj]/bs;
71: PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
72: MatSetValues(Gmat,1,&dest_row,1,&dest_col,&sv,ADD_VALUES);
73: }
74: MatRestoreRow(Amat,Ii,&ncols,&idx,&vals);
75: }
76: MatAssemblyBegin(Gmat,MAT_FINAL_ASSEMBLY);
77: MatAssemblyEnd(Gmat,MAT_FINAL_ASSEMBLY);
78: } else {
79: /* just copy scalar matrix - abs() not taken here but scaled later */
80: MatDuplicate(Amat, MAT_COPY_VALUES, &Gmat);
81: }
83: #if defined PETSC_GAMG_USE_LOG
84: PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);
85: #endif
87: *a_Gmat = Gmat;
88: return(0);
89: }
91: /* -------------------------------------------------------------------------- */
92: /*
93: PCGAMGFilterGraph - filter graph and symetrize if needed
95: Input Parameter:
96: . vfilter - threshold paramter [0,1)
97: . symm - symetrize?
98: In/Output Parameter:
99: . a_Gmat - original graph
100: */
103: PetscErrorCode PCGAMGFilterGraph(Mat *a_Gmat,const PetscReal vfilter,const PetscBool symm,const PetscInt verbose)
104: {
105: PetscErrorCode ierr;
106: PetscInt Istart,Iend,Ii,jj,ncols,nnz0,nnz1, NN, MM, nloc;
107: PetscMPIInt rank, size;
108: Mat Gmat = *a_Gmat, tGmat, matTrans;
109: MPI_Comm comm;
110: const PetscScalar *vals;
111: const PetscInt *idx;
112: PetscInt *d_nnz, *o_nnz;
113: Vec diag;
116: PetscObjectGetComm((PetscObject)Gmat,&comm);
117: MPI_Comm_rank(comm,&rank);
118: MPI_Comm_size(comm,&size);
119: MatGetOwnershipRange(Gmat, &Istart, &Iend);
120: nloc = Iend - Istart;
121: MatGetSize(Gmat, &MM, &NN);
122: #if defined PETSC_GAMG_USE_LOG
123: PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);
124: #endif
125: /* scale Gmat so filter works */
126: MatGetVecs(Gmat, &diag, 0);
127: MatGetDiagonal(Gmat, diag);
128: VecReciprocal(diag);
129: VecSqrtAbs(diag);
130: MatDiagonalScale(Gmat, diag, diag);
131: VecDestroy(&diag);
133: if (symm) {
134: MatTranspose(Gmat, MAT_INITIAL_MATRIX, &matTrans);
135: }
137: /* filter - dup zeros out matrix */
138: PetscMalloc(nloc*sizeof(PetscInt), &d_nnz);
139: PetscMalloc(nloc*sizeof(PetscInt), &o_nnz);
140: for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
141: MatGetRow(Gmat,Ii,&ncols,NULL,NULL);
142: d_nnz[jj] = ncols;
143: o_nnz[jj] = ncols;
144: MatRestoreRow(Gmat,Ii,&ncols,NULL,NULL);
145: if (symm) {
146: MatGetRow(matTrans,Ii,&ncols,NULL,NULL);
147: d_nnz[jj] += ncols;
148: o_nnz[jj] += ncols;
149: MatRestoreRow(matTrans,Ii,&ncols,NULL,NULL);
150: }
151: if (d_nnz[jj] > nloc) d_nnz[jj] = nloc;
152: if (o_nnz[jj] > (MM-nloc)) o_nnz[jj] = MM - nloc;
153: }
154: MatCreateAIJ(comm, nloc, nloc, MM, MM, 0, d_nnz, 0, o_nnz, &tGmat);
155: PetscFree(d_nnz);
156: PetscFree(o_nnz);
157: if (symm) {
158: MatDestroy(&matTrans);
159: }
161: for (Ii = Istart, nnz0 = nnz1 = 0; Ii < Iend; Ii++) {
162: MatGetRow(Gmat,Ii,&ncols,&idx,&vals);
163: for (jj=0; jj<ncols; jj++,nnz0++) {
164: PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
165: if (PetscRealPart(sv) > vfilter) {
166: nnz1++;
167: if (symm) {
168: sv *= 0.5;
169: MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);
170: MatSetValues(tGmat,1,&idx[jj],1,&Ii,&sv,ADD_VALUES);
171: } else {
172: MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);
173: }
174: }
175: }
176: MatRestoreRow(Gmat,Ii,&ncols,&idx,&vals);
177: }
178: MatAssemblyBegin(tGmat,MAT_FINAL_ASSEMBLY);
179: MatAssemblyEnd(tGmat,MAT_FINAL_ASSEMBLY);
181: #if defined PETSC_GAMG_USE_LOG
182: PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);
183: #endif
185: if (verbose) {
186: if (verbose == 1) {
187: PetscPrintf(comm,"\t[%d]%s %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%d)\n",rank,__FUNCT__,
188: 100.*(double)nnz1/(double)nnz0,vfilter,(double)nnz0/(double)nloc,MM);
189: } else {
190: PetscInt nnz[2],out[2];
191: nnz[0] = nnz0; nnz[1] = nnz1;
192: MPI_Allreduce(nnz, out, 2, MPIU_INT, MPI_SUM, comm);
193: PetscPrintf(comm,"\t[%d]%s %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%d)\n",rank,__FUNCT__,
194: 100.*(double)out[1]/(double)out[0],vfilter,(double)out[0]/(double)MM,MM);
195: }
196: }
197: MatDestroy(&Gmat);
198: *a_Gmat = tGmat;
199: return(0);
200: }
202: /* -------------------------------------------------------------------------- */
203: /*
204: PCGAMGGetDataWithGhosts - hacks into Mat MPIAIJ so this must have > 1 pe
206: Input Parameter:
207: . Gmat - MPIAIJ matrix for scattters
208: . data_sz - number of data terms per node (# cols in output)
209: . data_in[nloc*data_sz] - column oriented data
210: Output Parameter:
211: . a_stride - numbrt of rows of output
212: . a_data_out[stride*data_sz] - output data with ghosts
213: */
216: PetscErrorCode PCGAMGGetDataWithGhosts(const Mat Gmat,const PetscInt data_sz,const PetscReal data_in[],PetscInt *a_stride,PetscReal **a_data_out)
217: {
219: PetscMPIInt rank,size;
220: MPI_Comm comm;
221: Vec tmp_crds;
222: Mat_MPIAIJ *mpimat = (Mat_MPIAIJ*)Gmat->data;
223: PetscInt nnodes,num_ghosts,dir,kk,jj,my0,Iend,nloc;
224: PetscScalar *data_arr;
225: PetscReal *datas;
226: PetscBool isMPIAIJ;
229: PetscObjectGetComm((PetscObject)Gmat,&comm);
230: PetscObjectTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPIAIJ);
231: MPI_Comm_rank(comm,&rank);
232: MPI_Comm_size(comm,&size);
233: MatGetOwnershipRange(Gmat, &my0, &Iend);
234: nloc = Iend - my0;
235: VecGetLocalSize(mpimat->lvec, &num_ghosts);
236: nnodes = num_ghosts + nloc;
237: *a_stride = nnodes;
238: MatGetVecs(Gmat, &tmp_crds, 0);
240: PetscMalloc(data_sz*nnodes*sizeof(PetscReal), &datas);
241: for (dir=0; dir<data_sz; dir++) {
242: /* set local, and global */
243: for (kk=0; kk<nloc; kk++) {
244: PetscInt gid = my0 + kk;
245: PetscScalar crd = (PetscScalar)data_in[dir*nloc + kk]; /* col oriented */
246: datas[dir*nnodes + kk] = PetscRealPart(crd);
248: VecSetValues(tmp_crds, 1, &gid, &crd, INSERT_VALUES);
249: }
250: VecAssemblyBegin(tmp_crds);
251: VecAssemblyEnd(tmp_crds);
252: /* get ghost datas */
253: VecScatterBegin(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);
254: VecScatterEnd(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);
255: VecGetArray(mpimat->lvec, &data_arr);
256: for (kk=nloc,jj=0;jj<num_ghosts;kk++,jj++) datas[dir*nnodes + kk] = PetscRealPart(data_arr[jj]);
257: VecRestoreArray(mpimat->lvec, &data_arr);
258: }
259: VecDestroy(&tmp_crds);
260: *a_data_out = datas;
261: return(0);
262: }
265: /* hash table stuff - simple, not dymanic, key >= 0, has table
266: *
267: * GAMGTableCreate
268: */
269: /* avoid overflow */
270: #define GAMG_HASH(key) ((7*key)%a_tab->size)
273: PetscErrorCode GAMGTableCreate(PetscInt a_size, GAMGHashTable *a_tab)
274: {
276: PetscInt kk;
279: a_tab->size = a_size;
281: PetscMalloc(a_size*sizeof(PetscInt), &a_tab->table);
282: PetscMalloc(a_size*sizeof(PetscInt), &a_tab->data);
283: for (kk=0; kk<a_size; kk++) a_tab->table[kk] = -1;
284: return(0);
285: }
289: PetscErrorCode GAMGTableDestroy(GAMGHashTable *a_tab)
290: {
294: PetscFree(a_tab->table);
295: PetscFree(a_tab->data);
296: return(0);
297: }
301: PetscErrorCode GAMGTableAdd(GAMGHashTable *a_tab, PetscInt a_key, PetscInt a_data)
302: {
303: PetscInt kk,idx;
306: if (a_key<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Negative key %d.",a_key);
307: for (kk = 0, idx = GAMG_HASH(a_key);
308: kk < a_tab->size;
309: kk++, idx = (idx==(a_tab->size-1)) ? 0 : idx + 1) {
311: if (a_tab->table[idx] == a_key) {
312: /* exists */
313: a_tab->data[idx] = a_data;
314: break;
315: } else if (a_tab->table[idx] == -1) {
316: /* add */
317: a_tab->table[idx] = a_key;
318: a_tab->data[idx] = a_data;
319: break;
320: }
321: }
322: if (kk==a_tab->size) {
323: /* this is not to efficient, waiting until completely full */
324: PetscInt oldsize = a_tab->size, new_size = 2*a_tab->size + 5, *oldtable = a_tab->table, *olddata = a_tab->data;
327: a_tab->size = new_size;
329: PetscMalloc(a_tab->size*sizeof(PetscInt), &a_tab->table);
330: PetscMalloc(a_tab->size*sizeof(PetscInt), &a_tab->data);
332: for (kk=0;kk<a_tab->size;kk++) a_tab->table[kk] = -1;
333: for (kk=0;kk<oldsize;kk++) {
334: if (oldtable[kk] != -1) {
335: GAMGTableAdd(a_tab, oldtable[kk], olddata[kk]);
336: }
337: }
338: PetscFree(oldtable);
339: PetscFree(olddata);
340: GAMGTableAdd(a_tab, a_key, a_data);
341: }
342: return(0);
343: }
347: PetscErrorCode GAMGTableFind(GAMGHashTable *a_tab, PetscInt a_key, PetscInt *a_data)
348: {
349: PetscInt kk,idx;
352: if (a_key<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Negative key %d.",a_key);
353: for (kk = 0, idx = GAMG_HASH(a_key); kk < a_tab->size; kk++, idx = (idx==(a_tab->size-1)) ? 0 : idx + 1) {
354: if (a_tab->table[idx] == a_key) {
355: *a_data = a_tab->data[idx];
356: break;
357: } else if (a_tab->table[idx] == -1) {
358: /* not here */
359: *a_data = -1;
360: break;
361: }
362: }
363: if (kk==a_tab->size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"key %d not found in table",a_key);
364: return(0);
365: }