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