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