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 for (i=0; i<n-1; i++) { 135 if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);} 136 if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);} 137 if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);} 138 if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);} 139 if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);} 140 } 141 142 for (i=0; i<n; i++) { 143 if (mg[i]->smoothd != mg[i]->smoothu) { 144 ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr); 145 } 146 ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr); 147 ierr = PetscFree(mg[i]);CHKERRQ(ierr); 148 } 149 ierr = PetscFree(mg);CHKERRQ(ierr); 150 PetscFunctionReturn(0); 151 } 152 153 154 155 EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**); 156 EXTERN PetscErrorCode PCMGFCycle_Private(PC_MG**); 157 EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**); 158 159 /* 160 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 161 or full cycle of multigrid. 162 163 Note: 164 A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 165 */ 166 #undef __FUNCT__ 167 #define __FUNCT__ "PCApply_MG" 168 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 169 { 170 PC_MG **mg = (PC_MG**)pc->data; 171 PetscErrorCode ierr; 172 PetscInt levels = mg[0]->levels; 173 174 PetscFunctionBegin; 175 mg[levels-1]->b = b; 176 mg[levels-1]->x = x; 177 if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) { 178 Vec tvec; 179 ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr); 180 ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr); 181 ierr = VecDestroy(tvec);CHKERRQ(ierr); 182 } 183 if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 184 ierr = VecSet(x,0.0);CHKERRQ(ierr); 185 ierr = PCMGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr); 186 } 187 else if (mg[0]->am == PC_MG_ADDITIVE) { 188 ierr = PCMGACycle_Private(mg);CHKERRQ(ierr); 189 } 190 else if (mg[0]->am == PC_MG_KASKADE) { 191 ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr); 192 } 193 else { 194 ierr = PCMGFCycle_Private(mg);CHKERRQ(ierr); 195 } 196 PetscFunctionReturn(0); 197 } 198 199 #undef __FUNCT__ 200 #define __FUNCT__ "PCApplyRichardson_MG" 201 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its) 202 { 203 PC_MG **mg = (PC_MG**)pc->data; 204 PetscErrorCode ierr; 205 PetscInt levels = mg[0]->levels; 206 PetscTruth converged = PETSC_FALSE; 207 208 PetscFunctionBegin; 209 mg[levels-1]->b = b; 210 mg[levels-1]->x = x; 211 212 mg[levels-1]->rtol = rtol; 213 mg[levels-1]->abstol = abstol; 214 mg[levels-1]->dtol = dtol; 215 if (rtol) { 216 /* compute initial residual norm for relative convergence test */ 217 PetscReal rnorm; 218 ierr = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr); 219 ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 220 mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol); 221 } else if (abstol) { 222 mg[levels-1]->ttol = abstol; 223 } else { 224 mg[levels-1]->ttol = 0.0; 225 } 226 227 while (its-- && !converged) { 228 ierr = PCMGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr); 229 } 230 PetscFunctionReturn(0); 231 } 232 233 #undef __FUNCT__ 234 #define __FUNCT__ "PCSetFromOptions_MG" 235 PetscErrorCode PCSetFromOptions_MG(PC pc) 236 { 237 PetscErrorCode ierr; 238 PetscInt m,levels = 1; 239 PetscTruth flg; 240 PC_MG **mg = (PC_MG**)pc->data; 241 PCMGType mgtype; 242 243 PetscFunctionBegin; 244 245 ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 246 if (!pc->data) { 247 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 248 ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 249 mg = (PC_MG**)pc->data; 250 } 251 mgtype = mg[0]->am; 252 ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","PCMGSetCycles",1,&m,&flg);CHKERRQ(ierr); 253 if (flg) { 254 ierr = PCMGSetCycles(pc,m);CHKERRQ(ierr); 255 } 256 ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr); 257 if (flg) { 258 ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 259 } 260 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 261 if (flg) { 262 ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 263 } 264 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 265 if (flg) { 266 ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 267 } 268 ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 269 if (flg) {ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);} 270 ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr); 271 if (flg) { 272 PetscInt i; 273 char eventname[128]; 274 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 275 levels = mg[0]->levels; 276 for (i=0; i<levels; i++) { 277 sprintf(eventname,"MSetup Level %d",(int)i); 278 ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr); 279 sprintf(eventname,"MGSolve Level %d to 0",(int)i); 280 ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr); 281 } 282 } 283 ierr = PetscOptionsTail();CHKERRQ(ierr); 284 PetscFunctionReturn(0); 285 } 286 287 const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 288 289 #undef __FUNCT__ 290 #define __FUNCT__ "PCView_MG" 291 static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 292 { 293 PC_MG **mg = (PC_MG**)pc->data; 294 PetscErrorCode ierr; 295 PetscInt levels = mg[0]->levels,i; 296 PetscTruth iascii; 297 298 PetscFunctionBegin; 299 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 300 if (iascii) { 301 ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%D, pre-smooths=%D, post-smooths=%D\n", 302 PCMGTypes[mg[0]->am],levels,mg[0]->cycles,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr); 303 if (mg[0]->galerkin) { 304 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 305 } 306 for (i=0; i<levels; i++) { 307 if (!i) { 308 ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 309 } else { 310 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 311 } 312 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 313 ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr); 314 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 315 if (i && mg[i]->smoothd == mg[i]->smoothu) { 316 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 317 } else if (i){ 318 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 319 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 320 ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr); 321 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 322 } 323 } 324 } else { 325 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 326 } 327 PetscFunctionReturn(0); 328 } 329 330 /* 331 Calls setup for the KSP on each level 332 */ 333 #undef __FUNCT__ 334 #define __FUNCT__ "PCSetUp_MG" 335 static PetscErrorCode PCSetUp_MG(PC pc) 336 { 337 PC_MG **mg = (PC_MG**)pc->data; 338 PetscErrorCode ierr; 339 PetscInt i,n = mg[0]->levels; 340 PC cpc; 341 PetscTruth preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset; 342 PetscViewerASCIIMonitor ascii; 343 PetscViewer viewer = PETSC_NULL; 344 MPI_Comm comm; 345 Mat dA,dB; 346 MatStructure uflag; 347 Vec tvec; 348 349 PetscFunctionBegin; 350 351 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 352 /* so use those from global PC */ 353 /* Is this what we always want? What if user wants to keep old one? */ 354 ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 355 ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr); 356 if (!opsset || cpc->setupcalled == 2) { 357 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); 358 ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 359 } 360 361 if (mg[0]->galerkin) { 362 Mat B; 363 mg[0]->galerkinused = PETSC_TRUE; 364 /* currently only handle case where mat and pmat are the same on coarser levels */ 365 ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 366 if (!pc->setupcalled) { 367 for (i=n-2; i>-1; i--) { 368 ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 369 ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 370 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 371 dB = B; 372 } 373 ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr); 374 } else { 375 for (i=n-2; i>-1; i--) { 376 ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 377 ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 378 ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 379 dB = B; 380 } 381 } 382 } 383 384 if (!pc->setupcalled) { 385 ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr); 386 387 for (i=0; i<n; i++) { 388 if (monitor) { 389 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr); 390 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 391 ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 392 } 393 ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr); 394 } 395 for (i=1; i<n; i++) { 396 if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) { 397 if (monitor) { 398 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr); 399 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 400 ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 401 } 402 ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr); 403 } 404 } 405 for (i=1; i<n; i++) { 406 if (!mg[i]->residual) { 407 Mat mat; 408 ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 409 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 410 } 411 if (mg[i]->restrct && !mg[i]->interpolate) { 412 ierr = PCMGSetInterpolate(pc,i,mg[i]->restrct);CHKERRQ(ierr); 413 } 414 if (!mg[i]->restrct && mg[i]->interpolate) { 415 ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr); 416 } 417 #if defined(PETSC_USE_DEBUG) 418 if (!mg[i]->restrct || !mg[i]->interpolate) { 419 SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 420 } 421 #endif 422 } 423 for (i=0; i<n-1; i++) { 424 if (!mg[i]->b) { 425 Vec *vec; 426 ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 427 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 428 ierr = PetscFree(vec);CHKERRQ(ierr); 429 } 430 if (!mg[i]->r && i) { 431 ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 432 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 433 ierr = VecDestroy(tvec);CHKERRQ(ierr); 434 } 435 if (!mg[i]->x) { 436 ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 437 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 438 ierr = VecDestroy(tvec);CHKERRQ(ierr); 439 } 440 } 441 } 442 443 444 for (i=1; i<n; i++) { 445 if (mg[i]->smoothu == mg[i]->smoothd) { 446 /* if doing only down then initial guess is zero */ 447 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 448 } 449 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 450 ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr); 451 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 452 } 453 for (i=1; i<n; i++) { 454 if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 455 Mat downmat,downpmat; 456 MatStructure matflag; 457 PetscTruth opsset; 458 459 /* check if operators have been set for up, if not use down operators to set them */ 460 ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 461 if (!opsset) { 462 ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 463 ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 464 } 465 466 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 467 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 468 ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr); 469 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 470 } 471 } 472 473 /* 474 If coarse solver is not direct method then DO NOT USE preonly 475 */ 476 ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 477 if (preonly) { 478 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 479 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 480 ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 481 if (!lu && !redundant && !cholesky) { 482 ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 483 } 484 } 485 486 if (!pc->setupcalled) { 487 if (monitor) { 488 ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr); 489 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 490 ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);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 /* 500 Dump the interpolation/restriction matrices plus the 501 Jacobian/stiffness on each level. This allows Matlab users to 502 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 503 504 Only support one or the other at the same time. 505 */ 506 #if defined(PETSC_USE_SOCKET_VIEWER) 507 ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr); 508 if (dump) { 509 viewer = PETSC_VIEWER_SOCKET_(pc->comm); 510 } 511 #endif 512 ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr); 513 if (dump) { 514 viewer = PETSC_VIEWER_BINARY_(pc->comm); 515 } 516 517 if (viewer) { 518 for (i=1; i<n; i++) { 519 ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr); 520 } 521 for (i=0; i<n; i++) { 522 ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 523 ierr = MatView(pc->mat,viewer);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=0; 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 693 Options Database Key: 694 $ -pc_mg_galerkin 695 696 Level: intermediate 697 698 .keywords: MG, set, Galerkin 699 700 .seealso: PCMGGetGalerkin() 701 702 @*/ 703 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc) 704 { 705 PC_MG **mg; 706 PetscInt i,levels; 707 708 PetscFunctionBegin; 709 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 710 mg = (PC_MG**)pc->data; 711 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 712 levels = mg[0]->levels; 713 714 for (i=0; i<levels; i++) { 715 mg[i]->galerkin = PETSC_TRUE; 716 } 717 PetscFunctionReturn(0); 718 } 719 720 #undef __FUNCT__ 721 #define __FUNCT__ "PCMGGetGalerkin" 722 /*@ 723 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 724 A_i-1 = r_i * A_i * r_i^t 725 726 Not Collective 727 728 Input Parameter: 729 . pc - the multigrid context 730 731 Output Parameter: 732 . gelerkin - PETSC_TRUE or PETSC_FALSE 733 734 Options Database Key: 735 $ -pc_mg_galerkin 736 737 Level: intermediate 738 739 .keywords: MG, set, Galerkin 740 741 .seealso: PCMGSetGalerkin() 742 743 @*/ 744 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin) 745 { 746 PC_MG **mg; 747 748 PetscFunctionBegin; 749 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 750 mg = (PC_MG**)pc->data; 751 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 752 *galerkin = mg[0]->galerkin; 753 PetscFunctionReturn(0); 754 } 755 756 #undef __FUNCT__ 757 #define __FUNCT__ "PCMGSetNumberSmoothDown" 758 /*@ 759 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 760 use on all levels. Use PCMGGetSmootherDown() to set different 761 pre-smoothing steps on different levels. 762 763 Collective on PC 764 765 Input Parameters: 766 + mg - the multigrid context 767 - n - the number of smoothing steps 768 769 Options Database Key: 770 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 771 772 Level: advanced 773 774 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 775 776 .seealso: PCMGSetNumberSmoothUp() 777 @*/ 778 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n) 779 { 780 PC_MG **mg; 781 PetscErrorCode ierr; 782 PetscInt i,levels; 783 784 PetscFunctionBegin; 785 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 786 mg = (PC_MG**)pc->data; 787 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 788 levels = mg[0]->levels; 789 790 for (i=1; i<levels; i++) { 791 /* make sure smoother up and down are different */ 792 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 793 ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 794 mg[i]->default_smoothd = n; 795 } 796 PetscFunctionReturn(0); 797 } 798 799 #undef __FUNCT__ 800 #define __FUNCT__ "PCMGSetNumberSmoothUp" 801 /*@ 802 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 803 on all levels. Use PCMGGetSmootherUp() to set different numbers of 804 post-smoothing steps on different levels. 805 806 Collective on PC 807 808 Input Parameters: 809 + mg - the multigrid context 810 - n - the number of smoothing steps 811 812 Options Database Key: 813 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 814 815 Level: advanced 816 817 Note: this does not set a value on the coarsest grid, since we assume that 818 there is no separate smooth up on the coarsest grid. 819 820 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 821 822 .seealso: PCMGSetNumberSmoothDown() 823 @*/ 824 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n) 825 { 826 PC_MG **mg; 827 PetscErrorCode ierr; 828 PetscInt i,levels; 829 830 PetscFunctionBegin; 831 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 832 mg = (PC_MG**)pc->data; 833 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 834 levels = mg[0]->levels; 835 836 for (i=1; i<levels; i++) { 837 /* make sure smoother up and down are different */ 838 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 839 ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 840 mg[i]->default_smoothu = n; 841 } 842 PetscFunctionReturn(0); 843 } 844 845 /* ----------------------------------------------------------------------------------------*/ 846 847 /*MC 848 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 849 information about the coarser grid matrices and restriction/interpolation operators. 850 851 Options Database Keys: 852 + -pc_mg_levels <nlevels> - number of levels including finest 853 . -pc_mg_cycles 1 or 2 - for V or W-cycle 854 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 855 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 856 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 857 . -pc_mg_log - log information about time spent on each level of the solver 858 . -pc_mg_monitor - print information on the multigrid convergence 859 . -pc_mg_galerkin - use Galerkin process to compute coarser operators 860 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 861 to the Socket viewer for reading from Matlab. 862 863 Notes: 864 865 Level: intermediate 866 867 Concepts: multigrid 868 869 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 870 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycles(), PCMGSetNumberSmoothDown(), 871 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 872 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 873 PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 874 M*/ 875 876 EXTERN_C_BEGIN 877 #undef __FUNCT__ 878 #define __FUNCT__ "PCCreate_MG" 879 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 880 { 881 PetscFunctionBegin; 882 pc->ops->apply = PCApply_MG; 883 pc->ops->setup = PCSetUp_MG; 884 pc->ops->destroy = PCDestroy_MG; 885 pc->ops->setfromoptions = PCSetFromOptions_MG; 886 pc->ops->view = PCView_MG; 887 888 pc->data = (void*)0; 889 PetscFunctionReturn(0); 890 } 891 EXTERN_C_END 892