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