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