1 #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 2 #include <petsc-private/kspimpl.h> 3 #include <petscsf.h> 4 5 PetscFunctionList PCGAMGClassicalProlongatorList = NULL; 6 PetscBool PCGAMGClassicalPackageInitialized = PETSC_FALSE; 7 8 typedef struct { 9 PetscReal interp_threshold; /* interpolation threshold */ 10 char prolongtype[256]; 11 PetscInt nsmooths; /* number of jacobi smoothings on the prolongator */ 12 } PC_GAMG_Classical; 13 14 #undef __FUNCT__ 15 #define __FUNCT__ "PCGAMGClassicalSetType" 16 /*@C 17 PCGAMGClassicalSetType - Sets the type of classical interpolation to use 18 19 Collective on PC 20 21 Input Parameters: 22 . pc - the preconditioner context 23 24 Options Database Key: 25 . -pc_gamg_classical_type 26 27 Level: intermediate 28 29 .seealso: () 30 @*/ 31 PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type) 32 { 33 PetscErrorCode ierr; 34 35 PetscFunctionBegin; 36 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 37 ierr = PetscTryMethod(pc,"PCGAMGClassicalSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 38 PetscFunctionReturn(0); 39 } 40 41 #undef __FUNCT__ 42 #define __FUNCT__ "PCGAMGClassicalGetType" 43 /*@C 44 PCGAMGClassicalGetType - Gets the type of classical interpolation to use 45 46 Collective on PC 47 48 Input Parameter: 49 . pc - the preconditioner context 50 51 Output Parameter: 52 . type - the type used 53 54 Level: intermediate 55 56 .seealso: () 57 @*/ 58 PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type) 59 { 60 PetscErrorCode ierr; 61 62 PetscFunctionBegin; 63 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 64 ierr = PetscUseMethod(pc,"PCGAMGClassicalGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 65 PetscFunctionReturn(0); 66 } 67 68 #undef __FUNCT__ 69 #define __FUNCT__ "PCGAMGClassicalSetType_GAMG" 70 static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type) 71 { 72 PetscErrorCode ierr; 73 PC_MG *mg = (PC_MG*)pc->data; 74 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 75 PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 76 77 PetscFunctionBegin; 78 ierr = PetscStrcpy(cls->prolongtype,type);CHKERRQ(ierr); 79 PetscFunctionReturn(0); 80 } 81 82 #undef __FUNCT__ 83 #define __FUNCT__ "PCGAMGClassicalGetType_GAMG" 84 static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type) 85 { 86 PC_MG *mg = (PC_MG*)pc->data; 87 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 88 PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 89 90 PetscFunctionBegin; 91 *type = cls->prolongtype; 92 PetscFunctionReturn(0); 93 } 94 95 #undef __FUNCT__ 96 #define __FUNCT__ "PCGAMGGraph_Classical" 97 PetscErrorCode PCGAMGGraph_Classical(PC pc,const Mat A,Mat *G) 98 { 99 PetscInt s,f,n,idx,lidx,gidx; 100 PetscInt r,c,ncols; 101 const PetscInt *rcol; 102 const PetscScalar *rval; 103 PetscInt *gcol; 104 PetscScalar *gval; 105 PetscReal rmax; 106 PetscInt cmax = 0; 107 PC_MG *mg; 108 PC_GAMG *gamg; 109 PetscErrorCode ierr; 110 PetscInt *gsparse,*lsparse; 111 PetscScalar *Amax; 112 MatType mtype; 113 114 PetscFunctionBegin; 115 mg = (PC_MG *)pc->data; 116 gamg = (PC_GAMG *)mg->innerctx; 117 118 ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr); 119 n=f-s; 120 ierr = PetscMalloc1(n,&lsparse);CHKERRQ(ierr); 121 ierr = PetscMalloc1(n,&gsparse);CHKERRQ(ierr); 122 ierr = PetscMalloc1(n,&Amax);CHKERRQ(ierr); 123 124 for (r = 0;r < n;r++) { 125 lsparse[r] = 0; 126 gsparse[r] = 0; 127 } 128 129 for (r = s;r < f;r++) { 130 /* determine the maximum off-diagonal in each row */ 131 rmax = 0.; 132 ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 133 for (c = 0; c < ncols; c++) { 134 if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) { 135 rmax = PetscRealPart(-rval[c]); 136 } 137 } 138 Amax[r-s] = rmax; 139 if (ncols > cmax) cmax = ncols; 140 lidx = 0; 141 gidx = 0; 142 /* create the local and global sparsity patterns */ 143 for (c = 0; c < ncols; c++) { 144 if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 145 if (rcol[c] < f && rcol[c] >= s) { 146 lidx++; 147 } else { 148 gidx++; 149 } 150 } 151 } 152 ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 153 lsparse[r-s] = lidx; 154 gsparse[r-s] = gidx; 155 } 156 ierr = PetscMalloc1(cmax,&gval);CHKERRQ(ierr); 157 ierr = PetscMalloc1(cmax,&gcol);CHKERRQ(ierr); 158 159 ierr = MatCreate(PetscObjectComm((PetscObject)A),G); CHKERRQ(ierr); 160 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 161 ierr = MatSetType(*G,mtype);CHKERRQ(ierr); 162 ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 163 ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr); 164 ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr); 165 for (r = s;r < f;r++) { 166 ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 167 idx = 0; 168 for (c = 0; c < ncols; c++) { 169 /* classical strength of connection */ 170 if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 171 gcol[idx] = rcol[c]; 172 gval[idx] = rval[c]; 173 idx++; 174 } 175 } 176 ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr); 177 ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 178 } 179 ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 180 ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 181 182 ierr = PetscFree(gval);CHKERRQ(ierr); 183 ierr = PetscFree(gcol);CHKERRQ(ierr); 184 ierr = PetscFree(lsparse);CHKERRQ(ierr); 185 ierr = PetscFree(gsparse);CHKERRQ(ierr); 186 ierr = PetscFree(Amax);CHKERRQ(ierr); 187 PetscFunctionReturn(0); 188 } 189 190 191 #undef __FUNCT__ 192 #define __FUNCT__ "PCGAMGCoarsen_Classical" 193 PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists) 194 { 195 PetscErrorCode ierr; 196 MatCoarsen crs; 197 MPI_Comm fcomm = ((PetscObject)pc)->comm; 198 199 PetscFunctionBegin; 200 201 202 /* construct the graph if necessary */ 203 if (!G) { 204 SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening"); 205 } 206 207 ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr); 208 ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 209 ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr); 210 ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr); 211 ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 212 ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr); 213 ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 214 215 PetscFunctionReturn(0); 216 } 217 218 #undef __FUNCT__ 219 #define __FUNCT__ "PCGAMGProlongator_Classical_Direct" 220 PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P) 221 { 222 PetscErrorCode ierr; 223 PC_MG *mg = (PC_MG*)pc->data; 224 PC_GAMG *gamg = (PC_GAMG*)mg->innerctx; 225 PetscBool iscoarse,isMPIAIJ,isSEQAIJ; 226 PetscInt fn,cn,fs,fe,cs,ce,i,j,ncols,col,row_f,row_c,cmax=0,idx,noff; 227 PetscInt *lcid,*gcid,*lsparse,*gsparse,*colmap,*pcols; 228 const PetscInt *rcol; 229 PetscReal *Amax_pos,*Amax_neg; 230 PetscScalar g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta,pij; 231 PetscScalar *pvals; 232 const PetscScalar *rval; 233 Mat lA,gA=NULL; 234 MatType mtype; 235 Vec C,lvec; 236 PetscLayout clayout; 237 PetscSF sf; 238 Mat_MPIAIJ *mpiaij; 239 240 PetscFunctionBegin; 241 ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 242 fn = fe-fs; 243 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 244 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ);CHKERRQ(ierr); 245 if (!isMPIAIJ && !isSEQAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix"); 246 if (isMPIAIJ) { 247 mpiaij = (Mat_MPIAIJ*)A->data; 248 lA = mpiaij->A; 249 gA = mpiaij->B; 250 lvec = mpiaij->lvec; 251 ierr = VecGetSize(lvec,&noff);CHKERRQ(ierr); 252 colmap = mpiaij->garray; 253 ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 254 ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 255 ierr = PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap);CHKERRQ(ierr); 256 ierr = PetscMalloc1(noff,&gcid);CHKERRQ(ierr); 257 } else { 258 lA = A; 259 } 260 ierr = PetscMalloc1(fn,&lsparse);CHKERRQ(ierr); 261 ierr = PetscMalloc1(fn,&gsparse);CHKERRQ(ierr); 262 ierr = PetscMalloc1(fn,&lcid);CHKERRQ(ierr); 263 ierr = PetscMalloc1(fn,&Amax_pos);CHKERRQ(ierr); 264 ierr = PetscMalloc1(fn,&Amax_neg);CHKERRQ(ierr); 265 266 /* count the number of coarse unknowns */ 267 cn = 0; 268 for (i=0;i<fn;i++) { 269 /* filter out singletons */ 270 ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 271 lcid[i] = -1; 272 if (!iscoarse) { 273 cn++; 274 } 275 } 276 277 /* create the coarse vector */ 278 ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C);CHKERRQ(ierr); 279 ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr); 280 281 cn = 0; 282 for (i=0;i<fn;i++) { 283 ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 284 if (!iscoarse) { 285 lcid[i] = cs+cn; 286 cn++; 287 } else { 288 lcid[i] = -1; 289 } 290 } 291 292 if (gA) { 293 ierr = PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr); 294 ierr = PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr); 295 } 296 297 /* determine the biggest off-diagonal entries in each row */ 298 for (i=fs;i<fe;i++) { 299 Amax_pos[i-fs] = 0.; 300 Amax_neg[i-fs] = 0.; 301 ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 302 for(j=0;j<ncols;j++){ 303 if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]); 304 if ((PetscRealPart(rval[j]) > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]); 305 } 306 if (ncols > cmax) cmax = ncols; 307 ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 308 } 309 ierr = PetscMalloc1(cmax,&pcols);CHKERRQ(ierr); 310 ierr = PetscMalloc1(cmax,&pvals);CHKERRQ(ierr); 311 ierr = VecDestroy(&C);CHKERRQ(ierr); 312 313 /* count the on and off processor sparsity patterns for the prolongator */ 314 for (i=0;i<fn;i++) { 315 /* on */ 316 lsparse[i] = 0; 317 gsparse[i] = 0; 318 if (lcid[i] >= 0) { 319 lsparse[i] = 1; 320 gsparse[i] = 0; 321 } else { 322 ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 323 for (j = 0;j < ncols;j++) { 324 col = rcol[j]; 325 if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 326 lsparse[i] += 1; 327 } 328 } 329 ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 330 /* off */ 331 if (gA) { 332 ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 333 for (j = 0; j < ncols; j++) { 334 col = rcol[j]; 335 if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 336 gsparse[i] += 1; 337 } 338 } 339 ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 340 } 341 } 342 } 343 344 /* preallocate and create the prolongator */ 345 ierr = MatCreate(PetscObjectComm((PetscObject)A),P); CHKERRQ(ierr); 346 ierr = MatGetType(G,&mtype);CHKERRQ(ierr); 347 ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 348 349 ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 350 ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 351 ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 352 353 /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */ 354 for (i = 0;i < fn;i++) { 355 /* determine on or off */ 356 row_f = i + fs; 357 row_c = lcid[i]; 358 if (row_c >= 0) { 359 pij = 1.; 360 ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr); 361 } else { 362 g_pos = 0.; 363 g_neg = 0.; 364 a_pos = 0.; 365 a_neg = 0.; 366 diag = 0.; 367 368 /* local connections */ 369 ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 370 for (j = 0; j < ncols; j++) { 371 col = rcol[j]; 372 if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 373 if (PetscRealPart(rval[j]) > 0.) { 374 g_pos += rval[j]; 375 } else { 376 g_neg += rval[j]; 377 } 378 } 379 if (col != i) { 380 if (PetscRealPart(rval[j]) > 0.) { 381 a_pos += rval[j]; 382 } else { 383 a_neg += rval[j]; 384 } 385 } else { 386 diag = rval[j]; 387 } 388 } 389 ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 390 391 /* ghosted connections */ 392 if (gA) { 393 ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 394 for (j = 0; j < ncols; j++) { 395 col = rcol[j]; 396 if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 397 if (PetscRealPart(rval[j]) > 0.) { 398 g_pos += rval[j]; 399 } else { 400 g_neg += rval[j]; 401 } 402 } 403 if (PetscRealPart(rval[j]) > 0.) { 404 a_pos += rval[j]; 405 } else { 406 a_neg += rval[j]; 407 } 408 } 409 ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 410 } 411 412 if (g_neg == 0.) { 413 alpha = 0.; 414 } else { 415 alpha = -a_neg/g_neg; 416 } 417 418 if (g_pos == 0.) { 419 diag += a_pos; 420 beta = 0.; 421 } else { 422 beta = -a_pos/g_pos; 423 } 424 if (diag == 0.) { 425 invdiag = 0.; 426 } else invdiag = 1. / diag; 427 /* on */ 428 ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 429 idx = 0; 430 for (j = 0;j < ncols;j++) { 431 col = rcol[j]; 432 if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 433 row_f = i + fs; 434 row_c = lcid[col]; 435 /* set the values for on-processor ones */ 436 if (PetscRealPart(rval[j]) < 0.) { 437 pij = rval[j]*alpha*invdiag; 438 } else { 439 pij = rval[j]*beta*invdiag; 440 } 441 if (PetscAbsScalar(pij) != 0.) { 442 pvals[idx] = pij; 443 pcols[idx] = row_c; 444 idx++; 445 } 446 } 447 } 448 ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 449 /* off */ 450 if (gA) { 451 ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 452 for (j = 0; j < ncols; j++) { 453 col = rcol[j]; 454 if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 455 row_f = i + fs; 456 row_c = gcid[col]; 457 /* set the values for on-processor ones */ 458 if (PetscRealPart(rval[j]) < 0.) { 459 pij = rval[j]*alpha*invdiag; 460 } else { 461 pij = rval[j]*beta*invdiag; 462 } 463 if (PetscAbsScalar(pij) != 0.) { 464 pvals[idx] = pij; 465 pcols[idx] = row_c; 466 idx++; 467 } 468 } 469 } 470 ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 471 } 472 ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr); 473 } 474 } 475 476 ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 477 ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 478 479 ierr = PetscFree(lsparse);CHKERRQ(ierr); 480 ierr = PetscFree(gsparse);CHKERRQ(ierr); 481 ierr = PetscFree(pcols);CHKERRQ(ierr); 482 ierr = PetscFree(pvals);CHKERRQ(ierr); 483 ierr = PetscFree(Amax_pos);CHKERRQ(ierr); 484 ierr = PetscFree(Amax_neg);CHKERRQ(ierr); 485 ierr = PetscFree(lcid);CHKERRQ(ierr); 486 if (gA) { 487 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 488 ierr = PetscFree(gcid);CHKERRQ(ierr); 489 } 490 491 PetscFunctionReturn(0); 492 } 493 494 #undef __FUNCT__ 495 #define __FUNCT__ "PCGAMGTruncateProlongator_Private" 496 PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P) 497 { 498 PetscInt j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax; 499 PetscErrorCode ierr; 500 const PetscScalar *pval; 501 const PetscInt *pcol; 502 PetscScalar *pnval; 503 PetscInt *pncol; 504 PetscInt ncols; 505 Mat Pnew; 506 PetscInt *lsparse,*gsparse; 507 PetscReal pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg; 508 PC_MG *mg = (PC_MG*)pc->data; 509 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 510 PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 511 512 PetscFunctionBegin; 513 /* trim and rescale with reallocation */ 514 ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr); 515 ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr); 516 pn = pf-ps; 517 pcn = pcf-pcs; 518 ierr = PetscMalloc1(pn,&lsparse);CHKERRQ(ierr); 519 ierr = PetscMalloc1(pn,&gsparse);CHKERRQ(ierr); 520 /* allocate */ 521 cmax = 0; 522 for (i=ps;i<pf;i++) { 523 lsparse[i-ps] = 0; 524 gsparse[i-ps] = 0; 525 ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 526 if (ncols > cmax) { 527 cmax = ncols; 528 } 529 pmax_pos = 0.; 530 pmax_neg = 0.; 531 for (j=0;j<ncols;j++) { 532 if (PetscRealPart(pval[j]) > pmax_pos) { 533 pmax_pos = PetscRealPart(pval[j]); 534 } else if (PetscRealPart(pval[j]) < pmax_neg) { 535 pmax_neg = PetscRealPart(pval[j]); 536 } 537 } 538 for (j=0;j<ncols;j++) { 539 if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 540 if (pcol[j] >= pcs && pcol[j] < pcf) { 541 lsparse[i-ps]++; 542 } else { 543 gsparse[i-ps]++; 544 } 545 } 546 } 547 ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 548 } 549 550 ierr = PetscMalloc1(cmax,&pnval);CHKERRQ(ierr); 551 ierr = PetscMalloc1(cmax,&pncol);CHKERRQ(ierr); 552 553 ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr); 554 ierr = MatSetType(Pnew, MATAIJ);CHKERRQ(ierr); 555 ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 556 ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr); 557 ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr); 558 559 for (i=ps;i<pf;i++) { 560 ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 561 pmax_pos = 0.; 562 pmax_neg = 0.; 563 for (j=0;j<ncols;j++) { 564 if (PetscRealPart(pval[j]) > pmax_pos) { 565 pmax_pos = PetscRealPart(pval[j]); 566 } else if (PetscRealPart(pval[j]) < pmax_neg) { 567 pmax_neg = PetscRealPart(pval[j]); 568 } 569 } 570 pthresh_pos = 0.; 571 pthresh_neg = 0.; 572 ptot_pos = 0.; 573 ptot_neg = 0.; 574 for (j=0;j<ncols;j++) { 575 if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) { 576 pthresh_pos += PetscRealPart(pval[j]); 577 } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) { 578 pthresh_neg += PetscRealPart(pval[j]); 579 } 580 if (PetscRealPart(pval[j]) > 0.) { 581 ptot_pos += PetscRealPart(pval[j]); 582 } else { 583 ptot_neg += PetscRealPart(pval[j]); 584 } 585 } 586 if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos; 587 if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg; 588 idx=0; 589 for (j=0;j<ncols;j++) { 590 if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) { 591 pnval[idx] = ptot_pos*pval[j]; 592 pncol[idx] = pcol[j]; 593 idx++; 594 } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 595 pnval[idx] = ptot_neg*pval[j]; 596 pncol[idx] = pcol[j]; 597 idx++; 598 } 599 } 600 ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 601 ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr); 602 } 603 604 ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 605 ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 606 ierr = MatDestroy(P);CHKERRQ(ierr); 607 608 *P = Pnew; 609 ierr = PetscFree(lsparse);CHKERRQ(ierr); 610 ierr = PetscFree(gsparse);CHKERRQ(ierr); 611 ierr = PetscFree(pncol);CHKERRQ(ierr); 612 ierr = PetscFree(pnval);CHKERRQ(ierr); 613 PetscFunctionReturn(0); 614 } 615 616 #undef __FUNCT__ 617 #define __FUNCT__ "PCGAMGProlongator_Classical_Standard" 618 PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P) 619 { 620 PetscErrorCode ierr; 621 Mat lA,*lAs; 622 MatType mtype; 623 Vec cv; 624 PetscInt *gcid,*lcid,*lsparse,*gsparse,*picol; 625 PetscInt fs,fe,cs,ce,nl,i,j,k,li,lni,ci,ncols,maxcols,fn,cn,cid; 626 PetscMPIInt size; 627 const PetscInt *lidx,*icol,*gidx; 628 PetscBool iscoarse; 629 PetscScalar vi,pentry,pjentry; 630 PetscScalar *pcontrib,*pvcol; 631 const PetscScalar *vcol; 632 PetscReal diag,jdiag,jwttotal; 633 PetscInt pncols; 634 PetscSF sf; 635 PetscLayout clayout; 636 IS lis; 637 638 PetscFunctionBegin; 639 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 640 ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 641 fn = fe-fs; 642 ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr); 643 if (size > 1) { 644 ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 645 /* increase the overlap by two to get neighbors of neighbors */ 646 ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr); 647 ierr = ISSort(lis);CHKERRQ(ierr); 648 /* get the local part of A */ 649 ierr = MatGetSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs);CHKERRQ(ierr); 650 lA = lAs[0]; 651 /* build an SF out of it */ 652 ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr); 653 ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 654 ierr = ISGetIndices(lis,&lidx);CHKERRQ(ierr); 655 ierr = PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx);CHKERRQ(ierr); 656 ierr = ISRestoreIndices(lis,&lidx);CHKERRQ(ierr); 657 } else { 658 lA = A; 659 nl = fn; 660 } 661 /* create a communication structure for the overlapped portion and transmit coarse indices */ 662 ierr = PetscMalloc1(fn,&lsparse);CHKERRQ(ierr); 663 ierr = PetscMalloc1(fn,&gsparse);CHKERRQ(ierr); 664 ierr = PetscMalloc1(nl,&pcontrib);CHKERRQ(ierr); 665 /* create coarse vector */ 666 cn = 0; 667 for (i=0;i<fn;i++) { 668 ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr); 669 if (!iscoarse) { 670 cn++; 671 } 672 } 673 ierr = PetscMalloc1(fn,&gcid);CHKERRQ(ierr); 674 ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr); 675 ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr); 676 cn = 0; 677 for (i=0;i<fn;i++) { 678 ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 679 if (!iscoarse) { 680 gcid[i] = cs+cn; 681 cn++; 682 } else { 683 gcid[i] = -1; 684 } 685 } 686 if (size > 1) { 687 ierr = PetscMalloc1(nl,&lcid);CHKERRQ(ierr); 688 ierr = PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr); 689 ierr = PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr); 690 } else { 691 lcid = gcid; 692 } 693 /* count to preallocate the prolongator */ 694 ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr); 695 maxcols = 0; 696 /* count the number of unique contributing coarse cells for each fine */ 697 for (i=0;i<nl;i++) { 698 pcontrib[i] = 0.; 699 ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 700 if (gidx[i] >= fs && gidx[i] < fe) { 701 li = gidx[i] - fs; 702 lsparse[li] = 0; 703 gsparse[li] = 0; 704 cid = lcid[i]; 705 if (cid >= 0) { 706 lsparse[li] = 1; 707 } else { 708 for (j=0;j<ncols;j++) { 709 if (lcid[icol[j]] >= 0) { 710 pcontrib[icol[j]] = 1.; 711 } else { 712 ci = icol[j]; 713 ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 714 ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 715 for (k=0;k<ncols;k++) { 716 if (lcid[icol[k]] >= 0) { 717 pcontrib[icol[k]] = 1.; 718 } 719 } 720 ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 721 ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 722 } 723 } 724 for (j=0;j<ncols;j++) { 725 if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 726 lni = lcid[icol[j]]; 727 if (lni >= cs && lni < ce) { 728 lsparse[li]++; 729 } else { 730 gsparse[li]++; 731 } 732 pcontrib[icol[j]] = 0.; 733 } else { 734 ci = icol[j]; 735 ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 736 ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 737 for (k=0;k<ncols;k++) { 738 if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 739 lni = lcid[icol[k]]; 740 if (lni >= cs && lni < ce) { 741 lsparse[li]++; 742 } else { 743 gsparse[li]++; 744 } 745 pcontrib[icol[k]] = 0.; 746 } 747 } 748 ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 749 ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 750 } 751 } 752 } 753 if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li]; 754 } 755 ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 756 } 757 ierr = PetscMalloc1(maxcols,&picol);CHKERRQ(ierr); 758 ierr = PetscMalloc1(maxcols,&pvcol);CHKERRQ(ierr); 759 ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr); 760 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 761 ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 762 ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 763 ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 764 ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 765 for (i=0;i<nl;i++) { 766 diag = 0.; 767 if (gidx[i] >= fs && gidx[i] < fe) { 768 li = gidx[i] - fs; 769 pncols=0; 770 cid = lcid[i]; 771 if (cid >= 0) { 772 pncols = 1; 773 picol[0] = cid; 774 pvcol[0] = 1.; 775 } else { 776 ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 777 for (j=0;j<ncols;j++) { 778 pentry = vcol[j]; 779 if (lcid[icol[j]] >= 0) { 780 /* coarse neighbor */ 781 pcontrib[icol[j]] += pentry; 782 } else if (icol[j] != i) { 783 /* the neighbor is a strongly connected fine node */ 784 ci = icol[j]; 785 vi = vcol[j]; 786 ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 787 ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 788 jwttotal=0.; 789 jdiag = 0.; 790 for (k=0;k<ncols;k++) { 791 if (ci == icol[k]) { 792 jdiag = PetscRealPart(vcol[k]); 793 } 794 } 795 for (k=0;k<ncols;k++) { 796 if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 797 pjentry = vcol[k]; 798 jwttotal += PetscRealPart(pjentry); 799 } 800 } 801 if (jwttotal != 0.) { 802 jwttotal = PetscRealPart(vi)/jwttotal; 803 for (k=0;k<ncols;k++) { 804 if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 805 pjentry = vcol[k]*jwttotal; 806 pcontrib[icol[k]] += pjentry; 807 } 808 } 809 } else { 810 diag += PetscRealPart(vi); 811 } 812 ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 813 ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 814 } else { 815 diag += PetscRealPart(vcol[j]); 816 } 817 } 818 if (diag != 0.) { 819 diag = 1./diag; 820 for (j=0;j<ncols;j++) { 821 if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 822 /* the neighbor is a coarse node */ 823 if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) { 824 lni = lcid[icol[j]]; 825 pvcol[pncols] = -pcontrib[icol[j]]*diag; 826 picol[pncols] = lni; 827 pncols++; 828 } 829 pcontrib[icol[j]] = 0.; 830 } else { 831 /* the neighbor is a strongly connected fine node */ 832 ci = icol[j]; 833 ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 834 ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 835 for (k=0;k<ncols;k++) { 836 if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 837 if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) { 838 lni = lcid[icol[k]]; 839 pvcol[pncols] = -pcontrib[icol[k]]*diag; 840 picol[pncols] = lni; 841 pncols++; 842 } 843 pcontrib[icol[k]] = 0.; 844 } 845 } 846 ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 847 ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 848 } 849 pcontrib[icol[j]] = 0.; 850 } 851 ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 852 } 853 } 854 ci = gidx[i]; 855 li = gidx[i] - fs; 856 if (pncols > 0) { 857 ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr); 858 } 859 } 860 } 861 ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr); 862 ierr = PetscFree(pcontrib);CHKERRQ(ierr); 863 ierr = PetscFree(picol);CHKERRQ(ierr); 864 ierr = PetscFree(pvcol);CHKERRQ(ierr); 865 ierr = PetscFree(lsparse);CHKERRQ(ierr); 866 ierr = PetscFree(gsparse);CHKERRQ(ierr); 867 ierr = ISDestroy(&lis);CHKERRQ(ierr); 868 ierr = PetscFree(gcid);CHKERRQ(ierr); 869 if (size > 1) { 870 ierr = PetscFree(lcid);CHKERRQ(ierr); 871 ierr = MatDestroyMatrices(1,&lAs);CHKERRQ(ierr); 872 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 873 } 874 ierr = VecDestroy(&cv);CHKERRQ(ierr); 875 ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 876 ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 877 PetscFunctionReturn(0); 878 } 879 880 #undef __FUNCT__ 881 #define __FUNCT__ "PCGAMGOptProl_Classical_Jacobi" 882 PetscErrorCode PCGAMGOptProl_Classical_Jacobi(PC pc,const Mat A,Mat *P) 883 { 884 885 PetscErrorCode ierr; 886 PetscInt f,s,n,cf,cs,i,idx; 887 PetscInt *coarserows; 888 PetscInt ncols; 889 const PetscInt *pcols; 890 const PetscScalar *pvals; 891 Mat Pnew; 892 Vec diag; 893 PC_MG *mg = (PC_MG*)pc->data; 894 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 895 PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 896 897 PetscFunctionBegin; 898 if (cls->nsmooths == 0) { 899 ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 900 PetscFunctionReturn(0); 901 } 902 ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr); 903 n = f-s; 904 ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr); 905 ierr = PetscMalloc1(n,&coarserows);CHKERRQ(ierr); 906 /* identify the rows corresponding to coarse unknowns */ 907 idx = 0; 908 for (i=s;i<f;i++) { 909 ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 910 /* assume, for now, that it's a coarse unknown if it has a single unit entry */ 911 if (ncols == 1) { 912 if (pvals[0] == 1.) { 913 coarserows[idx] = i; 914 idx++; 915 } 916 } 917 ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 918 } 919 ierr = MatCreateVecs(A,&diag,0);CHKERRQ(ierr); 920 ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr); 921 ierr = VecReciprocal(diag);CHKERRQ(ierr); 922 for (i=0;i<cls->nsmooths;i++) { 923 ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr); 924 ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr); 925 ierr = MatDiagonalScale(Pnew,diag,0);CHKERRQ(ierr); 926 ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 927 ierr = MatDestroy(P);CHKERRQ(ierr); 928 *P = Pnew; 929 Pnew = NULL; 930 } 931 ierr = VecDestroy(&diag);CHKERRQ(ierr); 932 ierr = PetscFree(coarserows);CHKERRQ(ierr); 933 ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 934 PetscFunctionReturn(0); 935 } 936 937 #undef __FUNCT__ 938 #define __FUNCT__ "PCGAMGProlongator_Classical" 939 PetscErrorCode PCGAMGProlongator_Classical(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P) 940 { 941 PetscErrorCode ierr; 942 PetscErrorCode (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*); 943 PC_MG *mg = (PC_MG*)pc->data; 944 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 945 PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 946 947 PetscFunctionBegin; 948 ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr); 949 if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type"); 950 ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr); 951 PetscFunctionReturn(0); 952 } 953 954 #undef __FUNCT__ 955 #define __FUNCT__ "PCGAMGDestroy_Classical" 956 PetscErrorCode PCGAMGDestroy_Classical(PC pc) 957 { 958 PetscErrorCode ierr; 959 PC_MG *mg = (PC_MG*)pc->data; 960 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 961 962 PetscFunctionBegin; 963 ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 964 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr); 965 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",NULL);CHKERRQ(ierr); 966 PetscFunctionReturn(0); 967 } 968 969 #undef __FUNCT__ 970 #define __FUNCT__ "PCGAMGSetFromOptions_Classical" 971 PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc) 972 { 973 PC_MG *mg = (PC_MG*)pc->data; 974 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 975 PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 976 char tname[256]; 977 PetscErrorCode ierr; 978 PetscBool flg; 979 980 PetscFunctionBegin; 981 ierr = PetscOptionsHead("GAMG-Classical options");CHKERRQ(ierr); 982 ierr = PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation","PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr); 983 if (flg) { 984 ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr); 985 } 986 ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr); 987 ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr); 988 ierr = PetscOptionsTail();CHKERRQ(ierr); 989 PetscFunctionReturn(0); 990 } 991 992 #undef __FUNCT__ 993 #define __FUNCT__ "PCGAMGSetData_Classical" 994 PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A) 995 { 996 PC_MG *mg = (PC_MG*)pc->data; 997 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 998 999 PetscFunctionBegin; 1000 /* no data for classical AMG */ 1001 pc_gamg->data = NULL; 1002 pc_gamg->data_cell_cols = 0; 1003 pc_gamg->data_cell_rows = 0; 1004 pc_gamg->data_sz = 0; 1005 PetscFunctionReturn(0); 1006 } 1007 1008 1009 #undef __FUNCT__ 1010 #define __FUNCT__ "PCGAMGClassicalFinalizePackage" 1011 PetscErrorCode PCGAMGClassicalFinalizePackage(void) 1012 { 1013 PetscErrorCode ierr; 1014 1015 PetscFunctionBegin; 1016 PCGAMGClassicalPackageInitialized = PETSC_FALSE; 1017 ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr); 1018 PetscFunctionReturn(0); 1019 } 1020 1021 #undef __FUNCT__ 1022 #define __FUNCT__ "PCGAMGClassicalInitializePackage" 1023 PetscErrorCode PCGAMGClassicalInitializePackage(void) 1024 { 1025 PetscErrorCode ierr; 1026 1027 PetscFunctionBegin; 1028 if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0); 1029 ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr); 1030 ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr); 1031 ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr); 1032 PetscFunctionReturn(0); 1033 } 1034 1035 /* -------------------------------------------------------------------------- */ 1036 /* 1037 PCCreateGAMG_Classical 1038 1039 */ 1040 #undef __FUNCT__ 1041 #define __FUNCT__ "PCCreateGAMG_Classical" 1042 PetscErrorCode PCCreateGAMG_Classical(PC pc) 1043 { 1044 PetscErrorCode ierr; 1045 PC_MG *mg = (PC_MG*)pc->data; 1046 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1047 PC_GAMG_Classical *pc_gamg_classical; 1048 1049 PetscFunctionBegin; 1050 ierr = PCGAMGClassicalInitializePackage(); 1051 if (pc_gamg->subctx) { 1052 /* call base class */ 1053 ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr); 1054 } 1055 1056 /* create sub context for SA */ 1057 ierr = PetscNewLog(pc,&pc_gamg_classical);CHKERRQ(ierr); 1058 pc_gamg->subctx = pc_gamg_classical; 1059 pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 1060 /* reset does not do anything; setup not virtual */ 1061 1062 /* set internal function pointers */ 1063 pc_gamg->ops->destroy = PCGAMGDestroy_Classical; 1064 pc_gamg->ops->graph = PCGAMGGraph_Classical; 1065 pc_gamg->ops->coarsen = PCGAMGCoarsen_Classical; 1066 pc_gamg->ops->prolongator = PCGAMGProlongator_Classical; 1067 pc_gamg->ops->optprol = PCGAMGOptProl_Classical_Jacobi; 1068 pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 1069 1070 pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical; 1071 pc_gamg_classical->interp_threshold = 0.2; 1072 pc_gamg_classical->nsmooths = 0; 1073 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr); 1074 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",PCGAMGClassicalGetType_GAMG);CHKERRQ(ierr); 1075 ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr); 1076 PetscFunctionReturn(0); 1077 } 1078