1 #define PETSCKSP_DLL 2 3 /* 4 Defines the multigrid preconditioner interface. 5 */ 6 #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 7 8 9 /* 10 MGMCycle_Private - Given an MG structure created with MGCreate() runs 11 one multiplicative cycle down through the levels and 12 back up. 13 14 Input Parameter: 15 . mg - structure created with MGCreate(). 16 */ 17 #undef __FUNCT__ 18 #define __FUNCT__ "MGMCycle_Private" 19 PetscErrorCode MGMCycle_Private(MG *mglevels,PetscTruth *converged) 20 { 21 MG mg = *mglevels,mgc; 22 PetscErrorCode ierr; 23 PetscInt cycles = mg->cycles; 24 PetscScalar zero = 0.0; 25 26 PetscFunctionBegin; 27 if (converged) *converged = PETSC_FALSE; 28 29 if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);} 30 ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr); 31 if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);} 32 if (mg->level) { /* not the coarsest grid */ 33 ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr); 34 35 /* if on finest level and have convergence criteria set */ 36 if (mg->level == mg->levels-1 && mg->ttol) { 37 PetscReal rnorm; 38 ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr); 39 if (rnorm <= mg->ttol) { 40 *converged = PETSC_TRUE; 41 if (rnorm < mg->abstol) { 42 ierr = PetscLogInfo((0,"MGMCycle_Private:Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",rnorm,mg->abstol));CHKERRQ(ierr); 43 } else { 44 ierr = PetscLogInfo((0,"MGMCycle_Private:Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",rnorm,mg->ttol));CHKERRQ(ierr); 45 } 46 PetscFunctionReturn(0); 47 } 48 } 49 50 mgc = *(mglevels - 1); 51 ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr); 52 ierr = VecSet(&zero,mgc->x);CHKERRQ(ierr); 53 while (cycles--) { 54 ierr = MGMCycle_Private(mglevels-1,converged);CHKERRQ(ierr); 55 } 56 ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr); 57 if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);} 58 ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr); 59 if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);} 60 } 61 PetscFunctionReturn(0); 62 } 63 64 /* 65 MGCreate_Private - Creates a MG structure for use with the 66 multigrid code. Level 0 is the coarsest. (But the 67 finest level is stored first in the array). 68 69 */ 70 #undef __FUNCT__ 71 #define __FUNCT__ "MGCreate_Private" 72 static PetscErrorCode MGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,MG **result) 73 { 74 MG *mg; 75 PetscErrorCode ierr; 76 PetscInt i; 77 PetscMPIInt size; 78 char *prefix; 79 PC ipc; 80 81 PetscFunctionBegin; 82 ierr = PetscMalloc(levels*sizeof(MG),&mg);CHKERRQ(ierr); 83 ierr = PetscLogObjectMemory(pc,levels*(sizeof(MG)+sizeof(struct _MG)));CHKERRQ(ierr); 84 85 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 86 87 for (i=0; i<levels; i++) { 88 ierr = PetscNew(struct _MG,&mg[i]);CHKERRQ(ierr); 89 mg[i]->level = i; 90 mg[i]->levels = levels; 91 mg[i]->cycles = 1; 92 mg[i]->galerkin = PETSC_FALSE; 93 mg[i]->galerkinused = PETSC_FALSE; 94 mg[i]->default_smoothu = 1; 95 mg[i]->default_smoothd = 1; 96 97 if (comms) comm = comms[i]; 98 ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr); 99 ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr); 100 ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr); 101 102 /* do special stuff for coarse grid */ 103 if (!i && levels > 1) { 104 ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 105 106 /* coarse solve is (redundant) LU by default */ 107 ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 108 ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr); 109 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 110 if (size > 1) { 111 ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 112 ierr = PCRedundantGetPC(ipc,&ipc);CHKERRQ(ierr); 113 } 114 ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 115 116 } else { 117 char tprefix[128]; 118 sprintf(tprefix,"mg_levels_%d_",(int)i); 119 ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr); 120 } 121 ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr); 122 mg[i]->smoothu = mg[i]->smoothd; 123 mg[i]->rtol = 0.0; 124 mg[i]->abstol = 0.0; 125 mg[i]->dtol = 0.0; 126 mg[i]->ttol = 0.0; 127 mg[i]->eventsetup = 0; 128 mg[i]->eventsolve = 0; 129 } 130 *result = mg; 131 PetscFunctionReturn(0); 132 } 133 134 #undef __FUNCT__ 135 #define __FUNCT__ "PCDestroy_MG" 136 static PetscErrorCode PCDestroy_MG(PC pc) 137 { 138 MG *mg = (MG*)pc->data; 139 PetscErrorCode ierr; 140 PetscInt i,n = mg[0]->levels; 141 142 PetscFunctionBegin; 143 if (mg[0]->galerkinused) { 144 Mat B; 145 for (i=0; i<n-1; i++) { 146 ierr = KSPGetOperators(mg[i]->smoothd,0,&B,0);CHKERRQ(ierr); 147 ierr = MatDestroy(B);CHKERRQ(ierr); 148 } 149 } 150 151 for (i=0; i<n-1; i++) { 152 if (mg[i]->r) {ierr = VecDestroy(mg[i]->r);CHKERRQ(ierr);} 153 if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);} 154 if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);} 155 if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);} 156 if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);} 157 } 158 if (mg[n-1]->r) {ierr = VecDestroy(mg[n-1]->r);CHKERRQ(ierr);} 159 160 for (i=0; i<n; i++) { 161 if (mg[i]->smoothd != mg[i]->smoothu) { 162 ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr); 163 } 164 ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr); 165 ierr = PetscFree(mg[i]);CHKERRQ(ierr); 166 } 167 ierr = PetscFree(mg);CHKERRQ(ierr); 168 PetscFunctionReturn(0); 169 } 170 171 172 173 EXTERN PetscErrorCode MGACycle_Private(MG*); 174 EXTERN PetscErrorCode MGFCycle_Private(MG*); 175 EXTERN PetscErrorCode MGKCycle_Private(MG*); 176 177 /* 178 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 179 or full cycle of multigrid. 180 181 Note: 182 A simple wrapper which calls MGMCycle(),MGACycle(), or MGFCycle(). 183 */ 184 #undef __FUNCT__ 185 #define __FUNCT__ "PCApply_MG" 186 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 187 { 188 MG *mg = (MG*)pc->data; 189 PetscScalar zero = 0.0; 190 PetscErrorCode ierr; 191 PetscInt levels = mg[0]->levels; 192 193 PetscFunctionBegin; 194 mg[levels-1]->b = b; 195 mg[levels-1]->x = x; 196 if (mg[0]->am == MGMULTIPLICATIVE) { 197 ierr = VecSet(&zero,x);CHKERRQ(ierr); 198 ierr = MGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr); 199 } 200 else if (mg[0]->am == MGADDITIVE) { 201 ierr = MGACycle_Private(mg);CHKERRQ(ierr); 202 } 203 else if (mg[0]->am == MGKASKADE) { 204 ierr = MGKCycle_Private(mg);CHKERRQ(ierr); 205 } 206 else { 207 ierr = MGFCycle_Private(mg);CHKERRQ(ierr); 208 } 209 PetscFunctionReturn(0); 210 } 211 212 #undef __FUNCT__ 213 #define __FUNCT__ "PCApplyRichardson_MG" 214 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its) 215 { 216 MG *mg = (MG*)pc->data; 217 PetscErrorCode ierr; 218 PetscInt levels = mg[0]->levels; 219 PetscTruth converged = PETSC_FALSE; 220 221 PetscFunctionBegin; 222 mg[levels-1]->b = b; 223 mg[levels-1]->x = x; 224 225 mg[levels-1]->rtol = rtol; 226 mg[levels-1]->abstol = abstol; 227 mg[levels-1]->dtol = dtol; 228 if (rtol) { 229 /* compute initial residual norm for relative convergence test */ 230 PetscReal rnorm; 231 ierr = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr); 232 ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 233 mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol); 234 } else if (abstol) { 235 mg[levels-1]->ttol = abstol; 236 } else { 237 mg[levels-1]->ttol = 0.0; 238 } 239 240 while (its-- && !converged) { 241 ierr = MGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr); 242 } 243 PetscFunctionReturn(0); 244 } 245 246 #undef __FUNCT__ 247 #define __FUNCT__ "PCSetFromOptions_MG" 248 static PetscErrorCode PCSetFromOptions_MG(PC pc) 249 { 250 PetscErrorCode ierr; 251 PetscInt indx,m,levels = 1; 252 PetscTruth flg; 253 const char *type[] = {"additive","multiplicative","full","cascade","kascade"}; 254 255 PetscFunctionBegin; 256 257 ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 258 if (!pc->data) { 259 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","MGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 260 ierr = MGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 261 } 262 ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 263 if (flg) { 264 ierr = MGSetCycles(pc,m);CHKERRQ(ierr); 265 } 266 ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","MGSetGalerkin",&flg);CHKERRQ(ierr); 267 if (flg) { 268 ierr = MGSetGalerkin(pc);CHKERRQ(ierr); 269 } 270 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 271 if (flg) { 272 ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 273 } 274 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 275 if (flg) { 276 ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 277 } 278 ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr); 279 if (flg) { 280 MGType mg = (MGType) 0; 281 switch (indx) { 282 case 0: 283 mg = MGADDITIVE; 284 break; 285 case 1: 286 mg = MGMULTIPLICATIVE; 287 break; 288 case 2: 289 mg = MGFULL; 290 break; 291 case 3: 292 mg = MGKASKADE; 293 break; 294 case 4: 295 mg = MGKASKADE; 296 break; 297 } 298 ierr = MGSetType(pc,mg);CHKERRQ(ierr); 299 } 300 ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr); 301 if (flg) { 302 MG *mg = (MG*)pc->data; 303 PetscInt i; 304 char eventname[128]; 305 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 306 levels = mg[0]->levels; 307 for (i=0; i<levels; i++) { 308 sprintf(eventname,"MSetup Level %d",(int)i); 309 ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr); 310 sprintf(eventname,"MGSolve Level %d to 0",(int)i); 311 ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr); 312 } 313 } 314 ierr = PetscOptionsTail();CHKERRQ(ierr); 315 PetscFunctionReturn(0); 316 } 317 318 #undef __FUNCT__ 319 #define __FUNCT__ "PCView_MG" 320 static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 321 { 322 MG *mg = (MG*)pc->data; 323 PetscErrorCode ierr; 324 PetscInt levels = mg[0]->levels,i; 325 const char *cstring; 326 PetscTruth iascii; 327 328 PetscFunctionBegin; 329 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 330 if (iascii) { 331 if (mg[0]->am == MGMULTIPLICATIVE) cstring = "multiplicative"; 332 else if (mg[0]->am == MGADDITIVE) cstring = "additive"; 333 else if (mg[0]->am == MGFULL) cstring = "full"; 334 else if (mg[0]->am == MGKASKADE) cstring = "Kaskade"; 335 else cstring = "unknown"; 336 ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%D, pre-smooths=%D, post-smooths=%D\n", 337 cstring,levels,mg[0]->cycles,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr); 338 if (mg[0]->galerkin) { 339 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 340 } 341 for (i=0; i<levels; i++) { 342 if (!i) { 343 ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 344 } else { 345 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 346 } 347 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 348 ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr); 349 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 350 if (i && mg[i]->smoothd == mg[i]->smoothu) { 351 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 352 } else if (i){ 353 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 354 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 355 ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr); 356 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 357 } 358 } 359 } else { 360 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 361 } 362 PetscFunctionReturn(0); 363 } 364 365 /* 366 Calls setup for the KSP on each level 367 */ 368 #undef __FUNCT__ 369 #define __FUNCT__ "PCSetUp_MG" 370 static PetscErrorCode PCSetUp_MG(PC pc) 371 { 372 MG *mg = (MG*)pc->data; 373 PetscErrorCode ierr; 374 PetscInt i,n = mg[0]->levels; 375 PC cpc; 376 PetscTruth preonly,lu,redundant,monitor = PETSC_FALSE,dump; 377 PetscViewer ascii; 378 MPI_Comm comm; 379 Mat dA,dB; 380 MatStructure uflag; 381 382 PetscFunctionBegin; 383 if (!pc->setupcalled) { 384 ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr); 385 386 for (i=0; i<n; i++) { 387 if (monitor) { 388 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr); 389 ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 390 ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr); 391 ierr = KSPSetMonitor(mg[i]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 392 } 393 ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr); 394 } 395 for (i=1; i<n; i++) { 396 if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 397 if (monitor) { 398 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr); 399 ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 400 ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr); 401 ierr = KSPSetMonitor(mg[i]->smoothu,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 402 } 403 ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr); 404 } 405 } 406 for (i=1; i<n; i++) { 407 if (mg[i]->restrct && !mg[i]->interpolate) { 408 ierr = MGSetInterpolate(pc,i,mg[i]->restrct);CHKERRQ(ierr); 409 } 410 if (!mg[i]->restrct && mg[i]->interpolate) { 411 ierr = MGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr); 412 } 413 #if defined(PETSC_USE_DEBUG) 414 if (!mg[i]->restrct || !mg[i]->interpolate) { 415 SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 416 } 417 #endif 418 } 419 } 420 421 /* If user did not provide fine grid operators, use those from PC */ 422 /* BUG BUG BUG This will work ONLY the first time called: hence if the user changes 423 the PC matrices between solves PCMG will continue to use first set provided */ 424 ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 425 if (!dA && !dB) { 426 ierr = PetscLogInfo((pc,"PCSetUp_MG: Using outer operators to define finest grid operator \n because MGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n")); 427 ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,uflag);CHKERRQ(ierr); 428 } 429 430 if (mg[0]->galerkin) { 431 mg[0]->galerkinused = PETSC_TRUE; 432 /* currently only handle case where mat and pmat are the same on coarser levels */ 433 Mat B; 434 ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 435 if (!pc->setupcalled) { 436 for (i=n-2; i>-1; i--) { 437 ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 438 ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 439 dB = B; 440 } 441 } else { 442 for (i=n-2; i>-1; i--) { 443 ierr = KSPGetOperators(mg[i]->smoothd,0,&B,0);CHKERRQ(ierr); 444 ierr = MatPtAP(dB,mg[i]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 445 ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 446 dB = B; 447 } 448 } 449 } 450 451 for (i=1; i<n; i++) { 452 if (mg[i]->smoothu == mg[i]->smoothd) { 453 /* if doing only down then initial guess is zero */ 454 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 455 } 456 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 457 ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr); 458 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 459 } 460 for (i=1; i<n; i++) { 461 if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 462 PC uppc,downpc; 463 Mat downmat,downpmat,upmat,uppmat; 464 MatStructure matflag; 465 466 /* check if operators have been set for up, if not use down operators to set them */ 467 ierr = KSPGetPC(mg[i]->smoothu,&uppc);CHKERRQ(ierr); 468 ierr = PCGetOperators(uppc,&upmat,&uppmat,PETSC_NULL);CHKERRQ(ierr); 469 if (!upmat) { 470 ierr = KSPGetPC(mg[i]->smoothd,&downpc);CHKERRQ(ierr); 471 ierr = PCGetOperators(downpc,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 472 ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 473 } 474 475 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 476 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 477 ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr); 478 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 479 } 480 } 481 482 /* 483 If coarse solver is not direct method then DO NOT USE preonly 484 */ 485 ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 486 if (preonly) { 487 ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr); 488 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 489 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 490 if (!lu && !redundant) { 491 ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 492 } 493 } 494 495 if (!pc->setupcalled) { 496 if (monitor) { 497 ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr); 498 ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 499 ierr = PetscViewerASCIISetTab(ascii,n);CHKERRQ(ierr); 500 ierr = KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 501 } 502 ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr); 503 } 504 505 if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 506 ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr); 507 if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 508 509 /* 510 Dump the interpolation/restriction matrices to matlab plus the 511 Jacobian/stiffness on each level. This allows Matlab users to 512 easily check if the Galerkin condition A_c = R A_f R^T is satisfied */ 513 ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr); 514 if (dump) { 515 for (i=1; i<n; i++) { 516 ierr = MatView(mg[i]->restrct,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr); 517 } 518 for (i=0; i<n; i++) { 519 ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 520 ierr = MatView(pc->mat,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr); 521 } 522 } 523 ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr); 524 if (dump) { 525 for (i=1; i<n; i++) { 526 ierr = MatView(mg[i]->restrct,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr); 527 } 528 for (i=0; i<n; i++) { 529 ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 530 ierr = MatView(pc->mat,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr); 531 } 532 } 533 PetscFunctionReturn(0); 534 } 535 536 /* -------------------------------------------------------------------------------------*/ 537 538 #undef __FUNCT__ 539 #define __FUNCT__ "MGSetLevels" 540 /*@C 541 MGSetLevels - Sets the number of levels to use with MG. 542 Must be called before any other MG routine. 543 544 Collective on PC 545 546 Input Parameters: 547 + pc - the preconditioner context 548 . levels - the number of levels 549 - comms - optional communicators for each level; this is to allow solving the coarser problems 550 on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 551 552 Level: intermediate 553 554 Notes: 555 If the number of levels is one then the multigrid uses the -mg_levels prefix 556 for setting the level options rather than the -mg_coarse prefix. 557 558 .keywords: MG, set, levels, multigrid 559 560 .seealso: MGSetType(), MGGetLevels() 561 @*/ 562 PetscErrorCode PETSCKSP_DLLEXPORT MGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 563 { 564 PetscErrorCode ierr; 565 MG *mg; 566 567 PetscFunctionBegin; 568 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 569 570 if (pc->data) { 571 SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\ 572 make sure that you call MGSetLevels() before KSPSetFromOptions()"); 573 } 574 ierr = MGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr); 575 mg[0]->am = MGMULTIPLICATIVE; 576 pc->data = (void*)mg; 577 pc->ops->applyrichardson = PCApplyRichardson_MG; 578 PetscFunctionReturn(0); 579 } 580 581 #undef __FUNCT__ 582 #define __FUNCT__ "MGGetLevels" 583 /*@ 584 MGGetLevels - Gets the number of levels to use with MG. 585 586 Not Collective 587 588 Input Parameter: 589 . pc - the preconditioner context 590 591 Output parameter: 592 . levels - the number of levels 593 594 Level: advanced 595 596 .keywords: MG, get, levels, multigrid 597 598 .seealso: MGSetLevels() 599 @*/ 600 PetscErrorCode PETSCKSP_DLLEXPORT MGGetLevels(PC pc,PetscInt *levels) 601 { 602 MG *mg; 603 604 PetscFunctionBegin; 605 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 606 PetscValidIntPointer(levels,2); 607 608 mg = (MG*)pc->data; 609 *levels = mg[0]->levels; 610 PetscFunctionReturn(0); 611 } 612 613 #undef __FUNCT__ 614 #define __FUNCT__ "MGSetType" 615 /*@ 616 MGSetType - Determines the form of multigrid to use: 617 multiplicative, additive, full, or the Kaskade algorithm. 618 619 Collective on PC 620 621 Input Parameters: 622 + pc - the preconditioner context 623 - form - multigrid form, one of MGMULTIPLICATIVE, MGADDITIVE, 624 MGFULL, MGKASKADE 625 626 Options Database Key: 627 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 628 additive, full, kaskade 629 630 Level: advanced 631 632 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 633 634 .seealso: MGSetLevels() 635 @*/ 636 PetscErrorCode PETSCKSP_DLLEXPORT MGSetType(PC pc,MGType form) 637 { 638 MG *mg; 639 640 PetscFunctionBegin; 641 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 642 mg = (MG*)pc->data; 643 644 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 645 mg[0]->am = form; 646 if (form == MGMULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 647 else pc->ops->applyrichardson = 0; 648 PetscFunctionReturn(0); 649 } 650 651 #undef __FUNCT__ 652 #define __FUNCT__ "MGSetCycles" 653 /*@ 654 MGSetCycles - Sets the type cycles to use. Use MGSetCyclesOnLevel() for more 655 complicated cycling. 656 657 Collective on PC 658 659 Input Parameters: 660 + pc - the multigrid context 661 - n - the number of cycles 662 663 Options Database Key: 664 $ -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle. 665 666 Level: advanced 667 668 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 669 670 .seealso: MGSetCyclesOnLevel() 671 @*/ 672 PetscErrorCode PETSCKSP_DLLEXPORT MGSetCycles(PC pc,PetscInt n) 673 { 674 MG *mg; 675 PetscInt i,levels; 676 677 PetscFunctionBegin; 678 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 679 mg = (MG*)pc->data; 680 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 681 levels = mg[0]->levels; 682 683 for (i=0; i<levels; i++) { 684 mg[i]->cycles = n; 685 } 686 PetscFunctionReturn(0); 687 } 688 689 #undef __FUNCT__ 690 #define __FUNCT__ "MGSetGalerkin" 691 /*@ 692 MGSetGalerkin - Causes the coarser grid matrices to be computed from the 693 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 694 695 Collective on PC 696 697 Input Parameters: 698 + pc - the multigrid context 699 - n - the number of cycles 700 701 Options Database Key: 702 $ -pc_mg_galerkin 703 704 Level: intermediate 705 706 .keywords: MG, set, Galerkin 707 708 @*/ 709 PetscErrorCode PETSCKSP_DLLEXPORT MGSetGalerkin(PC pc) 710 { 711 MG *mg; 712 PetscInt i,levels; 713 714 PetscFunctionBegin; 715 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 716 mg = (MG*)pc->data; 717 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 718 levels = mg[0]->levels; 719 720 for (i=0; i<levels; i++) { 721 mg[i]->galerkin = PETSC_TRUE; 722 } 723 PetscFunctionReturn(0); 724 } 725 726 #undef __FUNCT__ 727 #define __FUNCT__ "MGCheck" 728 /*@ 729 MGCheck - Checks that all components of the MG structure have 730 been set. 731 732 Collective on PC 733 734 Input Parameters: 735 . mg - the MG structure 736 737 Level: advanced 738 739 .keywords: MG, check, set, multigrid 740 @*/ 741 PetscErrorCode PETSCKSP_DLLEXPORT MGCheck(PC pc) 742 { 743 MG *mg; 744 PetscInt i,n,count = 0; 745 746 PetscFunctionBegin; 747 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 748 mg = (MG*)pc->data; 749 750 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 751 752 n = mg[0]->levels; 753 754 for (i=1; i<n; i++) { 755 if (!mg[i]->restrct) { 756 (*PetscErrorPrintf)("No restrict set level %D \n",n-i); count++; 757 } 758 if (!mg[i]->interpolate) { 759 (*PetscErrorPrintf)("No interpolate set level %D \n",n-i); count++; 760 } 761 if (!mg[i]->r) { 762 (*PetscErrorPrintf)("No r set level %D \n",n-i); count++; 763 } 764 if (!mg[i-1]->x) { 765 (*PetscErrorPrintf)("No x set level %D \n",n-i); count++; 766 } 767 if (!mg[i-1]->b) { 768 (*PetscErrorPrintf)("No b set level %D \n",n-i); count++; 769 } 770 } 771 PetscFunctionReturn(count); 772 } 773 774 775 #undef __FUNCT__ 776 #define __FUNCT__ "MGSetNumberSmoothDown" 777 /*@ 778 MGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 779 use on all levels. Use MGGetSmootherDown() to set different 780 pre-smoothing steps on different levels. 781 782 Collective on PC 783 784 Input Parameters: 785 + mg - the multigrid context 786 - n - the number of smoothing steps 787 788 Options Database Key: 789 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 790 791 Level: advanced 792 793 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 794 795 .seealso: MGSetNumberSmoothUp() 796 @*/ 797 PetscErrorCode PETSCKSP_DLLEXPORT MGSetNumberSmoothDown(PC pc,PetscInt n) 798 { 799 MG *mg; 800 PetscErrorCode ierr; 801 PetscInt i,levels; 802 803 PetscFunctionBegin; 804 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 805 mg = (MG*)pc->data; 806 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 807 levels = mg[0]->levels; 808 809 for (i=0; i<levels; i++) { 810 /* make sure smoother up and down are different */ 811 ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 812 ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 813 mg[i]->default_smoothd = n; 814 } 815 PetscFunctionReturn(0); 816 } 817 818 #undef __FUNCT__ 819 #define __FUNCT__ "MGSetNumberSmoothUp" 820 /*@ 821 MGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 822 on all levels. Use MGGetSmootherUp() to set different numbers of 823 post-smoothing steps on different levels. 824 825 Collective on PC 826 827 Input Parameters: 828 + mg - the multigrid context 829 - n - the number of smoothing steps 830 831 Options Database Key: 832 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 833 834 Level: advanced 835 836 Note: this does not set a value on the coarsest grid, since we assume that 837 there is no seperate smooth up on the coarsest grid. 838 839 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 840 841 .seealso: MGSetNumberSmoothDown() 842 @*/ 843 PetscErrorCode PETSCKSP_DLLEXPORT MGSetNumberSmoothUp(PC pc,PetscInt n) 844 { 845 MG *mg; 846 PetscErrorCode ierr; 847 PetscInt i,levels; 848 849 PetscFunctionBegin; 850 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 851 mg = (MG*)pc->data; 852 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 853 levels = mg[0]->levels; 854 855 for (i=1; i<levels; i++) { 856 /* make sure smoother up and down are different */ 857 ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 858 ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 859 mg[i]->default_smoothu = n; 860 } 861 PetscFunctionReturn(0); 862 } 863 864 /* ----------------------------------------------------------------------------------------*/ 865 866 /*MC 867 PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional 868 information about the coarser grid matrices and restriction/interpolation operators. 869 870 Options Database Keys: 871 + -pc_mg_levels <nlevels> - number of levels including finest 872 . -pc_mg_cycles 1 or 2 - for V or W-cycle 873 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 874 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 875 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 876 . -pc_mg_log - log information about time spent on each level of the solver 877 . -pc_mg_monitor - print information on the multigrid convergence 878 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 879 to the Socket viewer for reading from Matlab. 880 881 Notes: 882 883 Level: intermediate 884 885 Concepts: multigrid 886 887 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 888 MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(), 889 MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(), 890 MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(), 891 MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR() 892 M*/ 893 894 EXTERN_C_BEGIN 895 #undef __FUNCT__ 896 #define __FUNCT__ "PCCreate_MG" 897 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 898 { 899 PetscFunctionBegin; 900 pc->ops->apply = PCApply_MG; 901 pc->ops->setup = PCSetUp_MG; 902 pc->ops->destroy = PCDestroy_MG; 903 pc->ops->setfromoptions = PCSetFromOptions_MG; 904 pc->ops->view = PCView_MG; 905 906 pc->data = (void*)0; 907 PetscFunctionReturn(0); 908 } 909 EXTERN_C_END 910