xref: /petsc/src/ksp/pc/impls/gamg/util.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
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 
7 /*
8    Produces a set of block column indices of the matrix row, one for each block represented in the original row
9 
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;
17   PetscErrorCode ierr;
18 
19   PetscFunctionBegin;
20   ierr = MatGetRow(Amat,row,&nidx,&idx,NULL);CHKERRQ(ierr);
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   ierr = MatRestoreRow(Amat,row,&nidx,&idx,NULL);CHKERRQ(ierr);
29   *n = cnt+1;
30   PetscFunctionReturn(0);
31 }
32 
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
35 
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;
42   PetscErrorCode ierr;
43 
44   PetscFunctionBegin;
45   ierr = MatCollapseRow(Amat,start,bs,&nprev,cprev);CHKERRQ(ierr);
46   for (i=start+1; i<start+bs; i++) {
47     ierr  = MatCollapseRow(Amat,i,bs,&ncur,ccur);CHKERRQ(ierr);
48     ierr  = PetscMergeIntArray(nprev,cprev,ncur,ccur,&nprev,&merged);CHKERRQ(ierr);
49     cprevtmp = cprev; cprev = merged; merged = cprevtmp;
50   }
51   *ncollapsed = nprev;
52   if (collapsed) *collapsed  = cprev;
53   PetscFunctionReturn(0);
54 }
55 
56 
57 /* -------------------------------------------------------------------------- */
58 /*
59    PCGAMGCreateGraph - create simple scaled scalar graph from matrix
60 
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 {
68   PetscErrorCode ierr;
69   PetscInt       Istart,Iend,Ii,jj,kk,ncols,nloc,NN,MM,bs;
70   MPI_Comm       comm;
71   Mat            Gmat;
72   MatType        mtype;
73 
74   PetscFunctionBegin;
75   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
76   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
77   ierr = MatGetSize(Amat, &MM, &NN);CHKERRQ(ierr);
78   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
79   nloc = (Iend-Istart)/bs;
80 
81 #if defined PETSC_GAMG_USE_LOG
82   ierr = PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr);
83 #endif
84 
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;
90 
91     /*
92        Determine the preallocation needed for the scalar matrix derived from the vector matrix.
93     */
94 
95     ierr = PetscObjectBaseTypeCompare((PetscObject)Amat,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
96     ierr = PetscObjectBaseTypeCompare((PetscObject)Amat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
97     ierr = PetscMalloc2(nloc, &d_nnz,isseqaij ? 0 : nloc, &o_nnz);CHKERRQ(ierr);
98 
99     if (isseqaij) {
100       PetscInt       max_d_nnz;
101 
102       /*
103           Determine exact preallocation count for (sequential) scalar matrix
104       */
105       ierr = MatSeqAIJGetMaxRowNonzeros(Amat,&max_d_nnz);CHKERRQ(ierr);
106       max_d_nnz = PetscMin(nloc,bs*max_d_nnz);CHKERRQ(ierr);
107       ierr = PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2);CHKERRQ(ierr);
108       for (Ii = 0, jj = 0; Ii < Iend; Ii += bs, jj++) {
109         ierr = MatCollapseRows(Amat,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL);CHKERRQ(ierr);
110       }
111       ierr = PetscFree3(w0,w1,w2);CHKERRQ(ierr);
112 
113     } else if (ismpiaij) {
114       Mat            Daij,Oaij;
115       const PetscInt *garray;
116       PetscInt       max_d_nnz;
117 
118       ierr = MatMPIAIJGetSeqAIJ(Amat,&Daij,&Oaij,&garray);CHKERRQ(ierr);
119 
120       /*
121           Determine exact preallocation count for diagonal block portion of scalar matrix
122       */
123       ierr = MatSeqAIJGetMaxRowNonzeros(Daij,&max_d_nnz);CHKERRQ(ierr);
124       max_d_nnz = PetscMin(nloc,bs*max_d_nnz);CHKERRQ(ierr);
125       ierr = PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2);CHKERRQ(ierr);
126       for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) {
127         ierr = MatCollapseRows(Daij,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL);CHKERRQ(ierr);
128       }
129       ierr = PetscFree3(w0,w1,w2);CHKERRQ(ierr);
130 
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           ierr = MatGetRow(Oaij,Ii+kk,&ncols,0,0);CHKERRQ(ierr);
138           o_nnz[jj] += ncols;
139           ierr = MatRestoreRow(Oaij,Ii+kk,&ncols,0,0);CHKERRQ(ierr);
140         }
141         if (o_nnz[jj] > (NN/bs-nloc)) o_nnz[jj] = NN/bs-nloc;
142       }
143 
144     } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require AIJ matrix type");
145 
146     /* get scalar copy (norms) of matrix */
147     ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr);
148     ierr = MatCreate(comm, &Gmat);CHKERRQ(ierr);
149     ierr = MatSetSizes(Gmat,nloc,nloc,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
150     ierr = MatSetBlockSizes(Gmat, 1, 1);CHKERRQ(ierr);
151     ierr = MatSetType(Gmat, mtype);CHKERRQ(ierr);
152     ierr = MatSeqAIJSetPreallocation(Gmat,0,d_nnz);CHKERRQ(ierr);
153     ierr = MatMPIAIJSetPreallocation(Gmat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
154     ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
155 
156     for (Ii = Istart; Ii < Iend; Ii++) {
157       PetscInt dest_row = Ii/bs;
158       ierr = MatGetRow(Amat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr);
159       for (jj=0; jj<ncols; jj++) {
160         PetscInt    dest_col = idx[jj]/bs;
161         PetscScalar sv       = PetscAbs(PetscRealPart(vals[jj]));
162         ierr = MatSetValues(Gmat,1,&dest_row,1,&dest_col,&sv,ADD_VALUES);CHKERRQ(ierr);
163       }
164       ierr = MatRestoreRow(Amat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr);
165     }
166     ierr = MatAssemblyBegin(Gmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167     ierr = MatAssemblyEnd(Gmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168   } else {
169     /* just copy scalar matrix - abs() not taken here but scaled later */
170     ierr = MatDuplicate(Amat, MAT_COPY_VALUES, &Gmat);CHKERRQ(ierr);
171   }
172 
173 #if defined PETSC_GAMG_USE_LOG
174   ierr = PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr);
175 #endif
176 
177   *a_Gmat = Gmat;
178   PetscFunctionReturn(0);
179 }
180 
181 /* -------------------------------------------------------------------------- */
182 /*@C
183    PCGAMGFilterGraph - filter (remove zero and possibly small values from the) graph and make it symmetric if requested
184 
185    Collective on Mat
186 
187    Input Parameter:
188 +   a_Gmat - the graph
189 .   vfilter - threshold parameter [0,1)
190 -   symm - make the result symmetric
191 
192    Level: developer
193 
194    Notes:
195     This is called before graph coarsers are called.
196 
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;
211 
212   PetscFunctionBegin;
213 #if defined PETSC_GAMG_USE_LOG
214   ierr = PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr);
215 #endif
216   /* scale Gmat for all values between -1 and 1 */
217   ierr = MatCreateVecs(Gmat, &diag, 0);CHKERRQ(ierr);
218   ierr = MatGetDiagonal(Gmat, diag);CHKERRQ(ierr);
219   ierr = VecReciprocal(diag);CHKERRQ(ierr);
220   ierr = VecSqrtAbs(diag);CHKERRQ(ierr);
221   ierr = MatDiagonalScale(Gmat, diag, diag);CHKERRQ(ierr);
222   ierr = VecDestroy(&diag);CHKERRQ(ierr);
223 
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     ierr = PetscObjectBaseTypeCompare((PetscObject)Gmat,MATSEQAIJ,&isaij);CHKERRQ(ierr);
230     ierr = PetscObjectBaseTypeCompare((PetscObject)Gmat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
231     if (!isaij && !ismpiaij) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require (MPI)AIJ matrix type");
232     if (isaij) {
233       ierr = MatGetInfo(Gmat,MAT_LOCAL,&info);CHKERRQ(ierr);
234       ierr = MatSeqAIJGetArray(Gmat,&avals);CHKERRQ(ierr);
235       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
236       ierr = MatSeqAIJRestoreArray(Gmat,&avals);CHKERRQ(ierr);
237     } else {
238       Mat_MPIAIJ  *aij = (Mat_MPIAIJ*)Gmat->data;
239       ierr = MatGetInfo(aij->A,MAT_LOCAL,&info);CHKERRQ(ierr);
240       ierr = MatSeqAIJGetArray(aij->A,&avals);CHKERRQ(ierr);
241       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
242       ierr = MatSeqAIJRestoreArray(aij->A,&avals);CHKERRQ(ierr);
243       ierr = MatGetInfo(aij->B,MAT_LOCAL,&info);CHKERRQ(ierr);
244       ierr = MatSeqAIJGetArray(aij->B,&avals);CHKERRQ(ierr);
245       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
246       ierr = MatSeqAIJRestoreArray(aij->B,&avals);CHKERRQ(ierr);
247     }
248 #if defined PETSC_GAMG_USE_LOG
249     ierr = PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr);
250 #endif
251     PetscFunctionReturn(0);
252   }
253 
254   ierr = PetscObjectGetComm((PetscObject)Gmat,&comm);CHKERRQ(ierr);
255   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
256   ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr);
257   nloc = Iend - Istart;
258   ierr = MatGetSize(Gmat, &MM, &NN);CHKERRQ(ierr);
259 
260   if (symm) {
261     ierr = MatTranspose(Gmat, MAT_INITIAL_MATRIX, &matTrans);CHKERRQ(ierr);
262   }
263 
264   /* Determine upper bound on nonzeros needed in new filtered matrix */
265   ierr = PetscMalloc2(nloc, &d_nnz,nloc, &o_nnz);CHKERRQ(ierr);
266   for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
267     ierr      = MatGetRow(Gmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
268     d_nnz[jj] = ncols;
269     o_nnz[jj] = ncols;
270     ierr      = MatRestoreRow(Gmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
271     if (symm) {
272       ierr       = MatGetRow(matTrans,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
273       d_nnz[jj] += ncols;
274       o_nnz[jj] += ncols;
275       ierr       = MatRestoreRow(matTrans,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
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   ierr = MatGetType(Gmat,&mtype);CHKERRQ(ierr);
281   ierr = MatCreate(comm, &tGmat);CHKERRQ(ierr);
282   ierr = MatSetSizes(tGmat,nloc,nloc,MM,MM);CHKERRQ(ierr);
283   ierr = MatSetBlockSizes(tGmat, 1, 1);CHKERRQ(ierr);
284   ierr = MatSetType(tGmat, mtype);CHKERRQ(ierr);
285   ierr = MatSeqAIJSetPreallocation(tGmat,0,d_nnz);CHKERRQ(ierr);
286   ierr = MatMPIAIJSetPreallocation(tGmat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
287   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
288   if (symm) {
289     ierr = MatDestroy(&matTrans);CHKERRQ(ierr);
290   } else {
291     /* all entries are generated locally so MatAssembly will be slightly faster for large process counts */
292     ierr = MatSetOption(tGmat,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
293   }
294 
295   for (Ii = Istart, nnz0 = nnz1 = 0; Ii < Iend; Ii++) {
296     ierr = MatGetRow(Gmat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr);
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           ierr = MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);CHKERRQ(ierr);
304           ierr = MatSetValues(tGmat,1,&idx[jj],1,&Ii,&sv,ADD_VALUES);CHKERRQ(ierr);
305         } else {
306           ierr = MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);CHKERRQ(ierr);
307         }
308       }
309     }
310     ierr = MatRestoreRow(Gmat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr);
311   }
312   ierr = MatAssemblyBegin(tGmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
313   ierr = MatAssemblyEnd(tGmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
314 
315 #if defined PETSC_GAMG_USE_LOG
316   ierr = PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr);
317 #endif
318 
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     ierr = PetscInfo4(*a_Gmat,"\t %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%D)\n",t1,vfilter,t2,MM);CHKERRQ(ierr);
323   }
324 #endif
325   ierr    = MatDestroy(&Gmat);CHKERRQ(ierr);
326   *a_Gmat = tGmat;
327   PetscFunctionReturn(0);
328 }
329 
330 /* -------------------------------------------------------------------------- */
331 /*
332    PCGAMGGetDataWithGhosts - hacks into Mat MPIAIJ so this must have size > 1
333 
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 {
344   PetscErrorCode ierr;
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;
351 
352   PetscFunctionBegin;
353   ierr      = PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPIAIJ);CHKERRQ(ierr);
354   ierr      = MatGetOwnershipRange(Gmat, &my0, &Iend);CHKERRQ(ierr);
355   nloc      = Iend - my0;
356   ierr      = VecGetLocalSize(mpimat->lvec, &num_ghosts);CHKERRQ(ierr);
357   nnodes    = num_ghosts + nloc;
358   *a_stride = nnodes;
359   ierr      = MatCreateVecs(Gmat, &tmp_crds, 0);CHKERRQ(ierr);
360 
361   ierr = PetscMalloc1(data_sz*nnodes, &datas);CHKERRQ(ierr);
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);
368 
369       ierr = VecSetValues(tmp_crds, 1, &gid, &crd, INSERT_VALUES);CHKERRQ(ierr);
370     }
371     ierr = VecAssemblyBegin(tmp_crds);CHKERRQ(ierr);
372     ierr = VecAssemblyEnd(tmp_crds);CHKERRQ(ierr);
373     /* get ghost datas */
374     ierr = VecScatterBegin(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
375     ierr = VecScatterEnd(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
376     ierr = VecGetArray(mpimat->lvec, &data_arr);CHKERRQ(ierr);
377     for (kk=nloc,jj=0;jj<num_ghosts;kk++,jj++) datas[dir*nnodes + kk] = PetscRealPart(data_arr[jj]);
378     ierr = VecRestoreArray(mpimat->lvec, &data_arr);CHKERRQ(ierr);
379   }
380   ierr        = VecDestroy(&tmp_crds);CHKERRQ(ierr);
381   *a_data_out = datas;
382   PetscFunctionReturn(0);
383 }
384 
385 PetscErrorCode PCGAMGHashTableCreate(PetscInt a_size, PCGAMGHashTable *a_tab)
386 {
387   PetscErrorCode ierr;
388   PetscInt       kk;
389 
390   PetscFunctionBegin;
391   a_tab->size = a_size;
392   ierr = PetscMalloc2(a_size, &a_tab->table,a_size, &a_tab->data);CHKERRQ(ierr);
393   for (kk=0; kk<a_size; kk++) a_tab->table[kk] = -1;
394   PetscFunctionReturn(0);
395 }
396 
397 PetscErrorCode PCGAMGHashTableDestroy(PCGAMGHashTable *a_tab)
398 {
399   PetscErrorCode ierr;
400 
401   PetscFunctionBegin;
402   ierr = PetscFree2(a_tab->table,a_tab->data);CHKERRQ(ierr);
403   PetscFunctionReturn(0);
404 }
405 
406 PetscErrorCode PCGAMGHashTableAdd(PCGAMGHashTable *a_tab, PetscInt a_key, PetscInt a_data)
407 {
408   PetscInt kk,idx;
409 
410   PetscFunctionBegin;
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;
427     PetscErrorCode ierr;
428 
429     a_tab->size = new_size;
430     ierr = PetscMalloc2(a_tab->size, &a_tab->table,a_tab->size, &a_tab->data);CHKERRQ(ierr);
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         ierr = PCGAMGHashTableAdd(a_tab, oldtable[kk], olddata[kk]);CHKERRQ(ierr);
435        }
436     }
437     ierr = PetscFree2(oldtable,olddata);CHKERRQ(ierr);
438     ierr = PCGAMGHashTableAdd(a_tab, a_key, a_data);CHKERRQ(ierr);
439   }
440   PetscFunctionReturn(0);
441 }
442 
443