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 = PetscInfo2(0,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr); 34 } else { 35 ierr = PetscInfo2(0,"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 && levels > 1) { 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; 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 mg = (PC_MG**)pc->data; 258 } 259 mgtype = mg[0]->am; 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]->residual) { 382 Mat mat; 383 ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 384 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 385 } 386 if (mg[i]->restrct && !mg[i]->interpolate) { 387 ierr = PCMGSetInterpolate(pc,i,mg[i]->restrct);CHKERRQ(ierr); 388 } 389 if (!mg[i]->restrct && mg[i]->interpolate) { 390 ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr); 391 } 392 #if defined(PETSC_USE_DEBUG) 393 if (!mg[i]->restrct || !mg[i]->interpolate) { 394 SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 395 } 396 #endif 397 } 398 for (i=0; i<n-1; i++) { 399 if (!mg[i]->b) { 400 Mat mat; 401 Vec vec; 402 ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 403 ierr = MatGetVecs(mat,&vec,PETSC_NULL);CHKERRQ(ierr); 404 ierr = PCMGSetRhs(pc,i,vec);CHKERRQ(ierr); 405 } 406 if (!mg[i]->r && i) { 407 ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 408 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 409 ierr = VecDestroy(tvec);CHKERRQ(ierr); 410 } 411 if (!mg[i]->x) { 412 ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 413 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 414 ierr = VecDestroy(tvec);CHKERRQ(ierr); 415 } 416 } 417 } 418 419 /* If user did not provide fine grid operators, use those from PC */ 420 /* BUG BUG BUG This will work ONLY the first time called: hence if the user changes 421 the PC matrices between solves PCMG will continue to use first set provided */ 422 ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 423 if (!dA && !dB) { 424 ierr = PetscInfo(pc,"Using outer operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr); 425 ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,uflag);CHKERRQ(ierr); 426 } 427 428 if (mg[0]->galerkin) { 429 Mat B; 430 mg[0]->galerkinused = PETSC_TRUE; 431 /* currently only handle case where mat and pmat are the same on coarser levels */ 432 ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 433 if (!pc->setupcalled) { 434 for (i=n-2; i>-1; i--) { 435 ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 436 ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 437 dB = B; 438 } 439 } else { 440 for (i=n-2; i>-1; i--) { 441 ierr = KSPGetOperators(mg[i]->smoothd,0,&B,0);CHKERRQ(ierr); 442 ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 443 ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 444 dB = B; 445 } 446 } 447 } 448 449 for (i=1; i<n; i++) { 450 if (mg[i]->smoothu == mg[i]->smoothd) { 451 /* if doing only down then initial guess is zero */ 452 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 453 } 454 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 455 ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr); 456 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 457 } 458 for (i=1; i<n; i++) { 459 if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 460 PC uppc,downpc; 461 Mat downmat,downpmat,upmat,uppmat; 462 MatStructure matflag; 463 464 /* check if operators have been set for up, if not use down operators to set them */ 465 ierr = KSPGetPC(mg[i]->smoothu,&uppc);CHKERRQ(ierr); 466 ierr = PCGetOperators(uppc,&upmat,&uppmat,PETSC_NULL);CHKERRQ(ierr); 467 if (!upmat) { 468 ierr = KSPGetPC(mg[i]->smoothd,&downpc);CHKERRQ(ierr); 469 ierr = PCGetOperators(downpc,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 470 ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 471 } 472 473 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 474 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 475 ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr); 476 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 477 } 478 } 479 480 /* 481 If coarse solver is not direct method then DO NOT USE preonly 482 */ 483 ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 484 if (preonly) { 485 ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr); 486 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 487 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 488 ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 489 if (!lu && !redundant && !cholesky) { 490 ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 491 } 492 } 493 494 if (!pc->setupcalled) { 495 if (monitor) { 496 ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr); 497 ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 498 ierr = PetscViewerASCIISetTab(ascii,n);CHKERRQ(ierr); 499 ierr = KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 500 } 501 ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr); 502 } 503 504 if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 505 ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr); 506 if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 507 508 #if defined(PETSC_USE_SOCKET_VIEWER) 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 #endif 524 525 ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr); 526 if (dump) { 527 for (i=1; i<n; i++) { 528 ierr = MatView(mg[i]->restrct,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr); 529 } 530 for (i=0; i<n; i++) { 531 ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 532 ierr = MatView(pc->mat,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr); 533 } 534 } 535 PetscFunctionReturn(0); 536 } 537 538 /* -------------------------------------------------------------------------------------*/ 539 540 #undef __FUNCT__ 541 #define __FUNCT__ "PCMGSetLevels" 542 /*@C 543 PCMGSetLevels - Sets the number of levels to use with MG. 544 Must be called before any other MG routine. 545 546 Collective on PC 547 548 Input Parameters: 549 + pc - the preconditioner context 550 . levels - the number of levels 551 - comms - optional communicators for each level; this is to allow solving the coarser problems 552 on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 553 554 Level: intermediate 555 556 Notes: 557 If the number of levels is one then the multigrid uses the -mg_levels prefix 558 for setting the level options rather than the -mg_coarse prefix. 559 560 .keywords: MG, set, levels, multigrid 561 562 .seealso: PCMGSetType(), PCMGGetLevels() 563 @*/ 564 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 565 { 566 PetscErrorCode ierr; 567 PC_MG **mg=0; 568 569 PetscFunctionBegin; 570 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 571 572 if (pc->data) { 573 SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\ 574 make sure that you call PCMGSetLevels() before KSPSetFromOptions()"); 575 } 576 ierr = PCMGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr); 577 mg[0]->am = PC_MG_MULTIPLICATIVE; 578 pc->data = (void*)mg; 579 pc->ops->applyrichardson = PCApplyRichardson_MG; 580 PetscFunctionReturn(0); 581 } 582 583 #undef __FUNCT__ 584 #define __FUNCT__ "PCMGGetLevels" 585 /*@ 586 PCMGGetLevels - Gets the number of levels to use with MG. 587 588 Not Collective 589 590 Input Parameter: 591 . pc - the preconditioner context 592 593 Output parameter: 594 . levels - the number of levels 595 596 Level: advanced 597 598 .keywords: MG, get, levels, multigrid 599 600 .seealso: PCMGSetLevels() 601 @*/ 602 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels) 603 { 604 PC_MG **mg; 605 606 PetscFunctionBegin; 607 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 608 PetscValidIntPointer(levels,2); 609 610 mg = (PC_MG**)pc->data; 611 *levels = mg[0]->levels; 612 PetscFunctionReturn(0); 613 } 614 615 #undef __FUNCT__ 616 #define __FUNCT__ "PCMGSetType" 617 /*@ 618 PCMGSetType - Determines the form of multigrid to use: 619 multiplicative, additive, full, or the Kaskade algorithm. 620 621 Collective on PC 622 623 Input Parameters: 624 + pc - the preconditioner context 625 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 626 PC_MG_FULL, PC_MG_KASKADE 627 628 Options Database Key: 629 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 630 additive, full, kaskade 631 632 Level: advanced 633 634 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 635 636 .seealso: PCMGSetLevels() 637 @*/ 638 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form) 639 { 640 PC_MG **mg; 641 642 PetscFunctionBegin; 643 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 644 mg = (PC_MG**)pc->data; 645 646 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 647 mg[0]->am = form; 648 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 649 else pc->ops->applyrichardson = 0; 650 PetscFunctionReturn(0); 651 } 652 653 #undef __FUNCT__ 654 #define __FUNCT__ "PCMGSetCycles" 655 /*@ 656 PCMGSetCycles - Sets the type cycles to use. Use PCMGSetCyclesOnLevel() for more 657 complicated cycling. 658 659 Collective on PC 660 661 Input Parameters: 662 + pc - the multigrid context 663 - n - the number of cycles 664 665 Options Database Key: 666 $ -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle. 667 668 Level: advanced 669 670 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 671 672 .seealso: PCMGSetCyclesOnLevel() 673 @*/ 674 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycles(PC pc,PetscInt n) 675 { 676 PC_MG **mg; 677 PetscInt i,levels; 678 679 PetscFunctionBegin; 680 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 681 mg = (PC_MG**)pc->data; 682 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 683 levels = mg[0]->levels; 684 685 for (i=0; i<levels; i++) { 686 mg[i]->cycles = n; 687 } 688 PetscFunctionReturn(0); 689 } 690 691 #undef __FUNCT__ 692 #define __FUNCT__ "PCMGSetGalerkin" 693 /*@ 694 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 695 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 696 697 Collective on PC 698 699 Input Parameters: 700 . pc - the multigrid context 701 702 Options Database Key: 703 $ -pc_mg_galerkin 704 705 Level: intermediate 706 707 .keywords: MG, set, Galerkin 708 709 .seealso: PCMGGetGalerkin() 710 711 @*/ 712 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc) 713 { 714 PC_MG **mg; 715 PetscInt i,levels; 716 717 PetscFunctionBegin; 718 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 719 mg = (PC_MG**)pc->data; 720 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 721 levels = mg[0]->levels; 722 723 for (i=0; i<levels; i++) { 724 mg[i]->galerkin = PETSC_TRUE; 725 } 726 PetscFunctionReturn(0); 727 } 728 729 #undef __FUNCT__ 730 #define __FUNCT__ "PCMGGetGalerkin" 731 /*@ 732 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 733 A_i-1 = r_i * A_i * r_i^t 734 735 Not Collective 736 737 Input Parameter: 738 . pc - the multigrid context 739 740 Output Parameter: 741 . gelerkin - PETSC_TRUE or PETSC_FALSE 742 743 Options Database Key: 744 $ -pc_mg_galerkin 745 746 Level: intermediate 747 748 .keywords: MG, set, Galerkin 749 750 .seealso: PCMGSetGalerkin() 751 752 @*/ 753 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin) 754 { 755 PC_MG **mg; 756 757 PetscFunctionBegin; 758 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 759 mg = (PC_MG**)pc->data; 760 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 761 *galerkin = mg[0]->galerkin; 762 PetscFunctionReturn(0); 763 } 764 765 #undef __FUNCT__ 766 #define __FUNCT__ "PCMGSetNumberSmoothDown" 767 /*@ 768 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 769 use on all levels. Use PCMGGetSmootherDown() to set different 770 pre-smoothing steps on different levels. 771 772 Collective on PC 773 774 Input Parameters: 775 + mg - the multigrid context 776 - n - the number of smoothing steps 777 778 Options Database Key: 779 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 780 781 Level: advanced 782 783 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 784 785 .seealso: PCMGSetNumberSmoothUp() 786 @*/ 787 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(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]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 803 mg[i]->default_smoothd = n; 804 } 805 PetscFunctionReturn(0); 806 } 807 808 #undef __FUNCT__ 809 #define __FUNCT__ "PCMGSetNumberSmoothUp" 810 /*@ 811 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 812 on all levels. Use PCMGGetSmootherUp() to set different numbers of 813 post-smoothing steps on different levels. 814 815 Collective on PC 816 817 Input Parameters: 818 + mg - the multigrid context 819 - n - the number of smoothing steps 820 821 Options Database Key: 822 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 823 824 Level: advanced 825 826 Note: this does not set a value on the coarsest grid, since we assume that 827 there is no separate smooth up on the coarsest grid. 828 829 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 830 831 .seealso: PCMGSetNumberSmoothDown() 832 @*/ 833 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n) 834 { 835 PC_MG **mg; 836 PetscErrorCode ierr; 837 PetscInt i,levels; 838 839 PetscFunctionBegin; 840 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 841 mg = (PC_MG**)pc->data; 842 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 843 levels = mg[0]->levels; 844 845 for (i=1; i<levels; i++) { 846 /* make sure smoother up and down are different */ 847 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 848 ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 849 mg[i]->default_smoothu = n; 850 } 851 PetscFunctionReturn(0); 852 } 853 854 /* ----------------------------------------------------------------------------------------*/ 855 856 /*MC 857 PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional 858 information about the coarser grid matrices and restriction/interpolation operators. 859 860 Options Database Keys: 861 + -pc_mg_levels <nlevels> - number of levels including finest 862 . -pc_mg_cycles 1 or 2 - for V or W-cycle 863 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 864 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 865 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 866 . -pc_mg_log - log information about time spent on each level of the solver 867 . -pc_mg_monitor - print information on the multigrid convergence 868 . -pc_mg_galerkin - use Galerkin process to compute coarser operators 869 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 870 to the Socket viewer for reading from Matlab. 871 872 Notes: 873 874 Level: intermediate 875 876 Concepts: multigrid 877 878 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 879 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycles(), PCMGSetNumberSmoothDown(), 880 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 881 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 882 PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 883 M*/ 884 885 EXTERN_C_BEGIN 886 #undef __FUNCT__ 887 #define __FUNCT__ "PCCreate_MG" 888 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 889 { 890 PetscFunctionBegin; 891 pc->ops->apply = PCApply_MG; 892 pc->ops->setup = PCSetUp_MG; 893 pc->ops->destroy = PCDestroy_MG; 894 pc->ops->setfromoptions = PCSetFromOptions_MG; 895 pc->ops->view = PCView_MG; 896 897 pc->data = (void*)0; 898 PetscFunctionReturn(0); 899 } 900 EXTERN_C_END 901