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