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 73 PetscFunctionBegin; 74 ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 75 ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 76 ierr = MatGetSize(Amat, &MM, &NN);CHKERRQ(ierr); 77 ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 78 nloc = (Iend-Istart)/bs; 79 80 #if defined PETSC_GAMG_USE_LOG 81 ierr = PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr); 82 #endif 83 84 /* TODO GPU: these calls are potentially expensive if matrices are large and we want to use the GPU */ 85 /* A solution consists in providing a new API, MatAIJGetCollapsedAIJ, and each class can provide a fast 86 implementation */ 87 if (bs > 1) { 88 const PetscScalar *vals; 89 const PetscInt *idx; 90 PetscInt *d_nnz, *o_nnz,*w0,*w1,*w2; 91 PetscBool ismpiaij,isseqaij; 92 93 /* 94 Determine the preallocation needed for the scalar matrix derived from the vector matrix. 95 */ 96 97 ierr = PetscObjectBaseTypeCompare((PetscObject)Amat,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 98 ierr = PetscObjectBaseTypeCompare((PetscObject)Amat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 99 ierr = PetscMalloc2(nloc, &d_nnz,isseqaij ? 0 : nloc, &o_nnz);CHKERRQ(ierr); 100 101 if (isseqaij) { 102 PetscInt max_d_nnz; 103 104 /* 105 Determine exact preallocation count for (sequential) scalar matrix 106 */ 107 ierr = MatSeqAIJGetMaxRowNonzeros(Amat,&max_d_nnz);CHKERRQ(ierr); 108 max_d_nnz = PetscMin(nloc,bs*max_d_nnz);CHKERRQ(ierr); 109 ierr = PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2);CHKERRQ(ierr); 110 for (Ii = 0, jj = 0; Ii < Iend; Ii += bs, jj++) { 111 ierr = MatCollapseRows(Amat,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL);CHKERRQ(ierr); 112 } 113 ierr = PetscFree3(w0,w1,w2);CHKERRQ(ierr); 114 115 } else if (ismpiaij) { 116 Mat Daij,Oaij; 117 const PetscInt *garray; 118 PetscInt max_d_nnz; 119 120 ierr = MatMPIAIJGetSeqAIJ(Amat,&Daij,&Oaij,&garray);CHKERRQ(ierr); 121 122 /* 123 Determine exact preallocation count for diagonal block portion of scalar matrix 124 */ 125 ierr = MatSeqAIJGetMaxRowNonzeros(Daij,&max_d_nnz);CHKERRQ(ierr); 126 max_d_nnz = PetscMin(nloc,bs*max_d_nnz);CHKERRQ(ierr); 127 ierr = PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2);CHKERRQ(ierr); 128 for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) { 129 ierr = MatCollapseRows(Daij,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL);CHKERRQ(ierr); 130 } 131 ierr = PetscFree3(w0,w1,w2);CHKERRQ(ierr); 132 133 /* 134 Over estimate (usually grossly over), preallocation count for off-diagonal portion of scalar matrix 135 */ 136 for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) { 137 o_nnz[jj] = 0; 138 for (kk=0; kk<bs; kk++) { /* rows that get collapsed to a single row */ 139 ierr = MatGetRow(Oaij,Ii+kk,&ncols,NULL,NULL);CHKERRQ(ierr); 140 o_nnz[jj] += ncols; 141 ierr = MatRestoreRow(Oaij,Ii+kk,&ncols,NULL,NULL);CHKERRQ(ierr); 142 } 143 if (o_nnz[jj] > (NN/bs-nloc)) o_nnz[jj] = NN/bs-nloc; 144 } 145 146 } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require AIJ matrix type"); 147 148 /* get scalar copy (norms) of matrix */ 149 ierr = MatCreate(comm, &Gmat);CHKERRQ(ierr); 150 ierr = MatSetSizes(Gmat,nloc,nloc,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 151 ierr = MatSetBlockSizes(Gmat, 1, 1);CHKERRQ(ierr); 152 ierr = MatSetType(Gmat, MATAIJ);CHKERRQ(ierr); 153 ierr = MatSeqAIJSetPreallocation(Gmat,0,d_nnz);CHKERRQ(ierr); 154 ierr = MatMPIAIJSetPreallocation(Gmat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 155 ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 156 157 for (Ii = Istart; Ii < Iend; Ii++) { 158 PetscInt dest_row = Ii/bs; 159 ierr = MatGetRow(Amat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr); 160 for (jj=0; jj<ncols; jj++) { 161 PetscInt dest_col = idx[jj]/bs; 162 PetscScalar sv = PetscAbs(PetscRealPart(vals[jj])); 163 ierr = MatSetValues(Gmat,1,&dest_row,1,&dest_col,&sv,ADD_VALUES);CHKERRQ(ierr); 164 } 165 ierr = MatRestoreRow(Amat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr); 166 } 167 ierr = MatAssemblyBegin(Gmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168 ierr = MatAssemblyEnd(Gmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169 } else { 170 /* just copy scalar matrix - abs() not taken here but scaled later */ 171 ierr = MatDuplicate(Amat, MAT_COPY_VALUES, &Gmat);CHKERRQ(ierr); 172 } 173 ierr = MatPropagateSymmetryOptions(Amat, Gmat);CHKERRQ(ierr); 174 175 #if defined PETSC_GAMG_USE_LOG 176 ierr = PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr); 177 #endif 178 179 *a_Gmat = Gmat; 180 PetscFunctionReturn(0); 181 } 182 183 /* -------------------------------------------------------------------------- */ 184 /*@C 185 PCGAMGFilterGraph - filter (remove zero and possibly small values from the) graph and make it symmetric if requested 186 187 Collective on Mat 188 189 Input Parameter: 190 + a_Gmat - the graph 191 . vfilter - threshold parameter [0,1) 192 - symm - make the result symmetric 193 194 Level: developer 195 196 Notes: 197 This is called before graph coarsers are called. 198 199 .seealso: PCGAMGSetThreshold() 200 @*/ 201 PetscErrorCode PCGAMGFilterGraph(Mat *a_Gmat,PetscReal vfilter,PetscBool symm) 202 { 203 PetscErrorCode ierr; 204 PetscInt Istart,Iend,Ii,jj,ncols,nnz0,nnz1, NN, MM, nloc; 205 PetscMPIInt rank; 206 Mat Gmat = *a_Gmat, tGmat, matTrans; 207 MPI_Comm comm; 208 const PetscScalar *vals; 209 const PetscInt *idx; 210 PetscInt *d_nnz, *o_nnz; 211 Vec diag; 212 213 PetscFunctionBegin; 214 #if defined PETSC_GAMG_USE_LOG 215 ierr = PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr); 216 #endif 217 /* scale Gmat for all values between -1 and 1 */ 218 ierr = MatCreateVecs(Gmat, &diag, NULL);CHKERRQ(ierr); 219 ierr = MatGetDiagonal(Gmat, diag);CHKERRQ(ierr); 220 ierr = VecReciprocal(diag);CHKERRQ(ierr); 221 ierr = VecSqrtAbs(diag);CHKERRQ(ierr); 222 ierr = MatDiagonalScale(Gmat, diag, diag);CHKERRQ(ierr); 223 ierr = VecDestroy(&diag);CHKERRQ(ierr); 224 225 /* TODO GPU: optimization proposal, each class provides fast implementation of this 226 procedure via MatAbs API */ 227 if (vfilter < 0.0 && !symm) { 228 /* Just use the provided matrix as the graph but make all values positive */ 229 MatInfo info; 230 PetscScalar *avals; 231 PetscBool isaij,ismpiaij; 232 ierr = PetscObjectBaseTypeCompare((PetscObject)Gmat,MATSEQAIJ,&isaij);CHKERRQ(ierr); 233 ierr = PetscObjectBaseTypeCompare((PetscObject)Gmat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 234 if (!isaij && !ismpiaij) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require (MPI)AIJ matrix type"); 235 if (isaij) { 236 ierr = MatGetInfo(Gmat,MAT_LOCAL,&info);CHKERRQ(ierr); 237 ierr = MatSeqAIJGetArray(Gmat,&avals);CHKERRQ(ierr); 238 for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]); 239 ierr = MatSeqAIJRestoreArray(Gmat,&avals);CHKERRQ(ierr); 240 } else { 241 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)Gmat->data; 242 ierr = MatGetInfo(aij->A,MAT_LOCAL,&info);CHKERRQ(ierr); 243 ierr = MatSeqAIJGetArray(aij->A,&avals);CHKERRQ(ierr); 244 for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]); 245 ierr = MatSeqAIJRestoreArray(aij->A,&avals);CHKERRQ(ierr); 246 ierr = MatGetInfo(aij->B,MAT_LOCAL,&info);CHKERRQ(ierr); 247 ierr = MatSeqAIJGetArray(aij->B,&avals);CHKERRQ(ierr); 248 for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]); 249 ierr = MatSeqAIJRestoreArray(aij->B,&avals);CHKERRQ(ierr); 250 } 251 #if defined PETSC_GAMG_USE_LOG 252 ierr = PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr); 253 #endif 254 PetscFunctionReturn(0); 255 } 256 257 /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold? 258 Also, if the matrix is symmetric, can we skip this 259 operation? It can be very expensive on large matrices. */ 260 ierr = PetscObjectGetComm((PetscObject)Gmat,&comm);CHKERRQ(ierr); 261 ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 262 ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr); 263 nloc = Iend - Istart; 264 ierr = MatGetSize(Gmat, &MM, &NN);CHKERRQ(ierr); 265 266 if (symm) { 267 ierr = MatTranspose(Gmat, MAT_INITIAL_MATRIX, &matTrans);CHKERRQ(ierr); 268 } 269 270 /* Determine upper bound on nonzeros needed in new filtered matrix */ 271 ierr = PetscMalloc2(nloc, &d_nnz,nloc, &o_nnz);CHKERRQ(ierr); 272 for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) { 273 ierr = MatGetRow(Gmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr); 274 d_nnz[jj] = ncols; 275 o_nnz[jj] = ncols; 276 ierr = MatRestoreRow(Gmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr); 277 if (symm) { 278 ierr = MatGetRow(matTrans,Ii,&ncols,NULL,NULL);CHKERRQ(ierr); 279 d_nnz[jj] += ncols; 280 o_nnz[jj] += ncols; 281 ierr = MatRestoreRow(matTrans,Ii,&ncols,NULL,NULL);CHKERRQ(ierr); 282 } 283 if (d_nnz[jj] > nloc) d_nnz[jj] = nloc; 284 if (o_nnz[jj] > (MM-nloc)) o_nnz[jj] = MM - nloc; 285 } 286 ierr = MatCreate(comm, &tGmat);CHKERRQ(ierr); 287 ierr = MatSetSizes(tGmat,nloc,nloc,MM,MM);CHKERRQ(ierr); 288 ierr = MatSetBlockSizes(tGmat, 1, 1);CHKERRQ(ierr); 289 ierr = MatSetType(tGmat, MATAIJ);CHKERRQ(ierr); 290 ierr = MatSeqAIJSetPreallocation(tGmat,0,d_nnz);CHKERRQ(ierr); 291 ierr = MatMPIAIJSetPreallocation(tGmat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 292 ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 293 if (symm) { 294 ierr = MatDestroy(&matTrans);CHKERRQ(ierr); 295 } else { 296 /* all entries are generated locally so MatAssembly will be slightly faster for large process counts */ 297 ierr = MatSetOption(tGmat,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 298 } 299 300 for (Ii = Istart, nnz0 = nnz1 = 0; Ii < Iend; Ii++) { 301 ierr = MatGetRow(Gmat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr); 302 for (jj=0; jj<ncols; jj++,nnz0++) { 303 PetscScalar sv = PetscAbs(PetscRealPart(vals[jj])); 304 if (PetscRealPart(sv) > vfilter) { 305 nnz1++; 306 if (symm) { 307 sv *= 0.5; 308 ierr = MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);CHKERRQ(ierr); 309 ierr = MatSetValues(tGmat,1,&idx[jj],1,&Ii,&sv,ADD_VALUES);CHKERRQ(ierr); 310 } else { 311 ierr = MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);CHKERRQ(ierr); 312 } 313 } 314 } 315 ierr = MatRestoreRow(Gmat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr); 316 } 317 ierr = MatAssemblyBegin(tGmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 318 ierr = MatAssemblyEnd(tGmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 319 if (symm) { 320 ierr = MatSetOption(tGmat,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 321 } else { 322 ierr = MatPropagateSymmetryOptions(Gmat,tGmat);CHKERRQ(ierr); 323 } 324 #if defined PETSC_GAMG_USE_LOG 325 ierr = PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr); 326 #endif 327 328 #if defined(PETSC_USE_INFO) 329 { 330 double t1 = (!nnz0) ? 1. : 100.*(double)nnz1/(double)nnz0, t2 = (!nloc) ? 1. : (double)nnz0/(double)nloc; 331 ierr = PetscInfo4(*a_Gmat,"\t %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%D)\n",t1,vfilter,t2,MM);CHKERRQ(ierr); 332 } 333 #endif 334 ierr = MatDestroy(&Gmat);CHKERRQ(ierr); 335 *a_Gmat = tGmat; 336 PetscFunctionReturn(0); 337 } 338 339 /* -------------------------------------------------------------------------- */ 340 /* 341 PCGAMGGetDataWithGhosts - hacks into Mat MPIAIJ so this must have size > 1 342 343 Input Parameter: 344 . Gmat - MPIAIJ matrix for scattters 345 . data_sz - number of data terms per node (# cols in output) 346 . data_in[nloc*data_sz] - column oriented data 347 Output Parameter: 348 . a_stride - numbrt of rows of output 349 . a_data_out[stride*data_sz] - output data with ghosts 350 */ 351 PetscErrorCode PCGAMGGetDataWithGhosts(Mat Gmat,PetscInt data_sz,PetscReal data_in[],PetscInt *a_stride,PetscReal **a_data_out) 352 { 353 PetscErrorCode ierr; 354 Vec tmp_crds; 355 Mat_MPIAIJ *mpimat = (Mat_MPIAIJ*)Gmat->data; 356 PetscInt nnodes,num_ghosts,dir,kk,jj,my0,Iend,nloc; 357 PetscScalar *data_arr; 358 PetscReal *datas; 359 PetscBool isMPIAIJ; 360 361 PetscFunctionBegin; 362 ierr = PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPIAIJ);CHKERRQ(ierr); 363 ierr = MatGetOwnershipRange(Gmat, &my0, &Iend);CHKERRQ(ierr); 364 nloc = Iend - my0; 365 ierr = VecGetLocalSize(mpimat->lvec, &num_ghosts);CHKERRQ(ierr); 366 nnodes = num_ghosts + nloc; 367 *a_stride = nnodes; 368 ierr = MatCreateVecs(Gmat, &tmp_crds, NULL);CHKERRQ(ierr); 369 370 ierr = PetscMalloc1(data_sz*nnodes, &datas);CHKERRQ(ierr); 371 for (dir=0; dir<data_sz; dir++) { 372 /* set local, and global */ 373 for (kk=0; kk<nloc; kk++) { 374 PetscInt gid = my0 + kk; 375 PetscScalar crd = (PetscScalar)data_in[dir*nloc + kk]; /* col oriented */ 376 datas[dir*nnodes + kk] = PetscRealPart(crd); 377 378 ierr = VecSetValues(tmp_crds, 1, &gid, &crd, INSERT_VALUES);CHKERRQ(ierr); 379 } 380 ierr = VecAssemblyBegin(tmp_crds);CHKERRQ(ierr); 381 ierr = VecAssemblyEnd(tmp_crds);CHKERRQ(ierr); 382 /* get ghost datas */ 383 ierr = VecScatterBegin(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 384 ierr = VecScatterEnd(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 385 ierr = VecGetArray(mpimat->lvec, &data_arr);CHKERRQ(ierr); 386 for (kk=nloc,jj=0;jj<num_ghosts;kk++,jj++) datas[dir*nnodes + kk] = PetscRealPart(data_arr[jj]); 387 ierr = VecRestoreArray(mpimat->lvec, &data_arr);CHKERRQ(ierr); 388 } 389 ierr = VecDestroy(&tmp_crds);CHKERRQ(ierr); 390 *a_data_out = datas; 391 PetscFunctionReturn(0); 392 } 393 394 PetscErrorCode PCGAMGHashTableCreate(PetscInt a_size, PCGAMGHashTable *a_tab) 395 { 396 PetscErrorCode ierr; 397 PetscInt kk; 398 399 PetscFunctionBegin; 400 a_tab->size = a_size; 401 ierr = PetscMalloc2(a_size, &a_tab->table,a_size, &a_tab->data);CHKERRQ(ierr); 402 for (kk=0; kk<a_size; kk++) a_tab->table[kk] = -1; 403 PetscFunctionReturn(0); 404 } 405 406 PetscErrorCode PCGAMGHashTableDestroy(PCGAMGHashTable *a_tab) 407 { 408 PetscErrorCode ierr; 409 410 PetscFunctionBegin; 411 ierr = PetscFree2(a_tab->table,a_tab->data);CHKERRQ(ierr); 412 PetscFunctionReturn(0); 413 } 414 415 PetscErrorCode PCGAMGHashTableAdd(PCGAMGHashTable *a_tab, PetscInt a_key, PetscInt a_data) 416 { 417 PetscInt kk,idx; 418 419 PetscFunctionBegin; 420 if (a_key<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Negative key %D.",a_key); 421 for (kk = 0, idx = GAMG_HASH(a_key); kk < a_tab->size; kk++, idx = (idx==(a_tab->size-1)) ? 0 : idx + 1) { 422 if (a_tab->table[idx] == a_key) { 423 /* exists */ 424 a_tab->data[idx] = a_data; 425 break; 426 } else if (a_tab->table[idx] == -1) { 427 /* add */ 428 a_tab->table[idx] = a_key; 429 a_tab->data[idx] = a_data; 430 break; 431 } 432 } 433 if (kk==a_tab->size) { 434 /* this is not to efficient, waiting until completely full */ 435 PetscInt oldsize = a_tab->size, new_size = 2*a_tab->size + 5, *oldtable = a_tab->table, *olddata = a_tab->data; 436 PetscErrorCode ierr; 437 438 a_tab->size = new_size; 439 ierr = PetscMalloc2(a_tab->size, &a_tab->table,a_tab->size, &a_tab->data);CHKERRQ(ierr); 440 for (kk=0;kk<a_tab->size;kk++) a_tab->table[kk] = -1; 441 for (kk=0;kk<oldsize;kk++) { 442 if (oldtable[kk] != -1) { 443 ierr = PCGAMGHashTableAdd(a_tab, oldtable[kk], olddata[kk]);CHKERRQ(ierr); 444 } 445 } 446 ierr = PetscFree2(oldtable,olddata);CHKERRQ(ierr); 447 ierr = PCGAMGHashTableAdd(a_tab, a_key, a_data);CHKERRQ(ierr); 448 } 449 PetscFunctionReturn(0); 450 } 451