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