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 = (PetscInt) 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); /* pre-smooth */ 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); /* post smooth */ 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 = PC_MG_CYCLE_V; 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 mg[i]->cyclesperpcapply = 1; 121 } 122 *result = mg; 123 PetscFunctionReturn(0); 124 } 125 126 #undef __FUNCT__ 127 #define __FUNCT__ "PCDestroy_MG" 128 static PetscErrorCode PCDestroy_MG(PC pc) 129 { 130 PC_MG **mg = (PC_MG**)pc->data; 131 PetscErrorCode ierr; 132 PetscInt i,n = mg[0]->levels; 133 134 PetscFunctionBegin; 135 for (i=0; i<n-1; i++) { 136 if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);} 137 if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);} 138 if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);} 139 if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);} 140 if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);} 141 } 142 143 for (i=0; i<n; i++) { 144 if (mg[i]->smoothd != mg[i]->smoothu) { 145 ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr); 146 } 147 ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr); 148 ierr = PetscFree(mg[i]);CHKERRQ(ierr); 149 } 150 ierr = PetscFree(mg);CHKERRQ(ierr); 151 PetscFunctionReturn(0); 152 } 153 154 155 156 EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**); 157 EXTERN PetscErrorCode PCMGFCycle_Private(PC_MG**); 158 EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**); 159 160 /* 161 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 162 or full cycle of multigrid. 163 164 Note: 165 A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 166 */ 167 #undef __FUNCT__ 168 #define __FUNCT__ "PCApply_MG" 169 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 170 { 171 PC_MG **mg = (PC_MG**)pc->data; 172 PetscErrorCode ierr; 173 PetscInt levels = mg[0]->levels,i; 174 175 PetscFunctionBegin; 176 mg[levels-1]->b = b; 177 mg[levels-1]->x = x; 178 if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) { 179 Vec tvec; 180 ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr); 181 ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr); 182 ierr = VecDestroy(tvec);CHKERRQ(ierr); 183 } 184 if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 185 ierr = VecSet(x,0.0);CHKERRQ(ierr); 186 for (i=0; i<mg[0]->cyclesperpcapply; i++) { 187 ierr = PCMGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr); 188 } 189 } 190 else if (mg[0]->am == PC_MG_ADDITIVE) { 191 ierr = PCMGACycle_Private(mg);CHKERRQ(ierr); 192 } 193 else if (mg[0]->am == PC_MG_KASKADE) { 194 ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr); 195 } 196 else { 197 ierr = PCMGFCycle_Private(mg);CHKERRQ(ierr); 198 } 199 PetscFunctionReturn(0); 200 } 201 202 #undef __FUNCT__ 203 #define __FUNCT__ "PCApplyRichardson_MG" 204 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its) 205 { 206 PC_MG **mg = (PC_MG**)pc->data; 207 PetscErrorCode ierr; 208 PetscInt levels = mg[0]->levels; 209 PetscTruth converged = PETSC_FALSE; 210 211 PetscFunctionBegin; 212 mg[levels-1]->b = b; 213 mg[levels-1]->x = x; 214 215 mg[levels-1]->rtol = rtol; 216 mg[levels-1]->abstol = abstol; 217 mg[levels-1]->dtol = dtol; 218 if (rtol) { 219 /* compute initial residual norm for relative convergence test */ 220 PetscReal rnorm; 221 ierr = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr); 222 ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 223 mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol); 224 } else if (abstol) { 225 mg[levels-1]->ttol = abstol; 226 } else { 227 mg[levels-1]->ttol = 0.0; 228 } 229 230 while (its-- && !converged) { 231 ierr = PCMGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr); 232 } 233 PetscFunctionReturn(0); 234 } 235 236 #undef __FUNCT__ 237 #define __FUNCT__ "PCSetFromOptions_MG" 238 PetscErrorCode PCSetFromOptions_MG(PC pc) 239 { 240 PetscErrorCode ierr; 241 PetscInt m,levels = 1,cycles; 242 PetscTruth flg; 243 PC_MG **mg = (PC_MG**)pc->data; 244 PCMGType mgtype=PC_MG_MULTIPLICATIVE; 245 PCMGCycleType mgctype; 246 247 PetscFunctionBegin; 248 ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 249 if (!pc->data) { 250 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 251 ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 252 mg = (PC_MG**)pc->data; 253 } 254 mgctype = (PCMGCycleType)mg[0]->cycles; 255 ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 256 if (flg) { 257 ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 258 }; 259 ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr); 260 if (flg) { 261 ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 262 } 263 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 264 if (flg) { 265 ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 266 } 267 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 268 if (flg) { 269 ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 270 } 271 ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 272 if (flg) { 273 ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 274 } 275 if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 276 ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 277 if (flg) { 278 ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 279 } 280 } 281 ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr); 282 if (flg) { 283 PetscInt i; 284 char eventname[128]; 285 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 286 levels = mg[0]->levels; 287 for (i=0; i<levels; i++) { 288 sprintf(eventname,"MSetup Level %d",(int)i); 289 ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr); 290 sprintf(eventname,"MGSolve Level %d to 0",(int)i); 291 ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr); 292 } 293 } 294 ierr = PetscOptionsTail();CHKERRQ(ierr); 295 PetscFunctionReturn(0); 296 } 297 298 const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 299 const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 300 301 #undef __FUNCT__ 302 #define __FUNCT__ "PCView_MG" 303 static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 304 { 305 PC_MG **mg = (PC_MG**)pc->data; 306 PetscErrorCode ierr; 307 PetscInt levels = mg[0]->levels,i; 308 PetscTruth iascii; 309 310 PetscFunctionBegin; 311 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 312 if (iascii) { 313 ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s, pre-smooths=%D, post-smooths=%D\n", 314 PCMGTypes[mg[0]->am],levels,(mg[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w", 315 mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr); 316 if (mg[0]->galerkin) { 317 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 318 } 319 for (i=0; i<levels; i++) { 320 if (!i) { 321 ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 322 } else { 323 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 324 } 325 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 326 ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr); 327 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 328 if (i && mg[i]->smoothd == mg[i]->smoothu) { 329 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 330 } else if (i){ 331 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 332 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 333 ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr); 334 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 335 } 336 } 337 } else { 338 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 339 } 340 PetscFunctionReturn(0); 341 } 342 343 /* 344 Calls setup for the KSP on each level 345 */ 346 #undef __FUNCT__ 347 #define __FUNCT__ "PCSetUp_MG" 348 static PetscErrorCode PCSetUp_MG(PC pc) 349 { 350 PC_MG **mg = (PC_MG**)pc->data; 351 PetscErrorCode ierr; 352 PetscInt i,n = mg[0]->levels; 353 PC cpc; 354 PetscTruth preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset; 355 PetscViewerASCIIMonitor ascii; 356 PetscViewer viewer = PETSC_NULL; 357 MPI_Comm comm; 358 Mat dA,dB; 359 MatStructure uflag; 360 Vec tvec; 361 362 PetscFunctionBegin; 363 364 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 365 /* so use those from global PC */ 366 /* Is this what we always want? What if user wants to keep old one? */ 367 ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 368 ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr); 369 if (!opsset || cpc->setupcalled == 2) { 370 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); 371 ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 372 } 373 374 if (mg[0]->galerkin) { 375 Mat B; 376 mg[0]->galerkinused = PETSC_TRUE; 377 /* currently only handle case where mat and pmat are the same on coarser levels */ 378 ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 379 if (!pc->setupcalled) { 380 for (i=n-2; i>-1; i--) { 381 ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 382 ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 383 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 384 dB = B; 385 } 386 ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr); 387 } else { 388 for (i=n-2; i>-1; i--) { 389 ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 390 ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 391 ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 392 dB = B; 393 } 394 } 395 } 396 397 if (!pc->setupcalled) { 398 ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr); 399 400 for (i=0; i<n; i++) { 401 if (monitor) { 402 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr); 403 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 404 ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 405 } 406 ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr); 407 } 408 for (i=1; i<n; i++) { 409 if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) { 410 if (monitor) { 411 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr); 412 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 413 ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 414 } 415 ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr); 416 } 417 } 418 for (i=1; i<n; i++) { 419 if (!mg[i]->residual) { 420 Mat mat; 421 ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 422 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 423 } 424 if (mg[i]->restrct && !mg[i]->interpolate) { 425 ierr = PCMGSetInterpolate(pc,i,mg[i]->restrct);CHKERRQ(ierr); 426 } 427 if (!mg[i]->restrct && mg[i]->interpolate) { 428 ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr); 429 } 430 #if defined(PETSC_USE_DEBUG) 431 if (!mg[i]->restrct || !mg[i]->interpolate) { 432 SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 433 } 434 #endif 435 } 436 for (i=0; i<n-1; i++) { 437 if (!mg[i]->b) { 438 Vec *vec; 439 ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 440 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 441 ierr = PetscFree(vec);CHKERRQ(ierr); 442 } 443 if (!mg[i]->r && i) { 444 ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 445 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 446 ierr = VecDestroy(tvec);CHKERRQ(ierr); 447 } 448 if (!mg[i]->x) { 449 ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 450 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 451 ierr = VecDestroy(tvec);CHKERRQ(ierr); 452 } 453 } 454 } 455 456 457 for (i=1; i<n; i++) { 458 if (mg[i]->smoothu == mg[i]->smoothd) { 459 /* if doing only down then initial guess is zero */ 460 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 461 } 462 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 463 ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr); 464 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 465 } 466 for (i=1; i<n; i++) { 467 if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 468 Mat downmat,downpmat; 469 MatStructure matflag; 470 PetscTruth opsset; 471 472 /* check if operators have been set for up, if not use down operators to set them */ 473 ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 474 if (!opsset) { 475 ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 476 ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 477 } 478 479 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 480 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 481 ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr); 482 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 483 } 484 } 485 486 /* 487 If coarse solver is not direct method then DO NOT USE preonly 488 */ 489 ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 490 if (preonly) { 491 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 492 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 493 ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 494 if (!lu && !redundant && !cholesky) { 495 ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 496 } 497 } 498 499 if (!pc->setupcalled) { 500 if (monitor) { 501 ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr); 502 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 503 ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 504 } 505 ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr); 506 } 507 508 if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 509 ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr); 510 if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 511 512 /* 513 Dump the interpolation/restriction matrices plus the 514 Jacobian/stiffness on each level. This allows Matlab users to 515 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 516 517 Only support one or the other at the same time. 518 */ 519 #if defined(PETSC_USE_SOCKET_VIEWER) 520 ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr); 521 if (dump) { 522 viewer = PETSC_VIEWER_SOCKET_(pc->comm); 523 } 524 #endif 525 ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr); 526 if (dump) { 527 viewer = PETSC_VIEWER_BINARY_(pc->comm); 528 } 529 530 if (viewer) { 531 for (i=1; i<n; i++) { 532 ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr); 533 } 534 for (i=0; i<n; i++) { 535 ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 536 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 537 } 538 } 539 PetscFunctionReturn(0); 540 } 541 542 /* -------------------------------------------------------------------------------------*/ 543 544 #undef __FUNCT__ 545 #define __FUNCT__ "PCMGSetLevels" 546 /*@C 547 PCMGSetLevels - Sets the number of levels to use with MG. 548 Must be called before any other MG routine. 549 550 Collective on PC 551 552 Input Parameters: 553 + pc - the preconditioner context 554 . levels - the number of levels 555 - comms - optional communicators for each level; this is to allow solving the coarser problems 556 on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 557 558 Level: intermediate 559 560 Notes: 561 If the number of levels is one then the multigrid uses the -mg_levels prefix 562 for setting the level options rather than the -mg_coarse prefix. 563 564 .keywords: MG, set, levels, multigrid 565 566 .seealso: PCMGSetType(), PCMGGetLevels() 567 @*/ 568 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 569 { 570 PetscErrorCode ierr; 571 PC_MG **mg=0; 572 573 PetscFunctionBegin; 574 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 575 576 if (pc->data) { 577 SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\ 578 make sure that you call PCMGSetLevels() before KSPSetFromOptions()"); 579 } 580 ierr = PCMGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr); 581 mg[0]->am = PC_MG_MULTIPLICATIVE; 582 pc->data = (void*)mg; 583 pc->ops->applyrichardson = PCApplyRichardson_MG; 584 PetscFunctionReturn(0); 585 } 586 587 #undef __FUNCT__ 588 #define __FUNCT__ "PCMGGetLevels" 589 /*@ 590 PCMGGetLevels - Gets the number of levels to use with MG. 591 592 Not Collective 593 594 Input Parameter: 595 . pc - the preconditioner context 596 597 Output parameter: 598 . levels - the number of levels 599 600 Level: advanced 601 602 .keywords: MG, get, levels, multigrid 603 604 .seealso: PCMGSetLevels() 605 @*/ 606 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels) 607 { 608 PC_MG **mg; 609 610 PetscFunctionBegin; 611 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 612 PetscValidIntPointer(levels,2); 613 614 mg = (PC_MG**)pc->data; 615 *levels = mg[0]->levels; 616 PetscFunctionReturn(0); 617 } 618 619 #undef __FUNCT__ 620 #define __FUNCT__ "PCMGSetType" 621 /*@ 622 PCMGSetType - Determines the form of multigrid to use: 623 multiplicative, additive, full, or the Kaskade algorithm. 624 625 Collective on PC 626 627 Input Parameters: 628 + pc - the preconditioner context 629 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 630 PC_MG_FULL, PC_MG_KASKADE 631 632 Options Database Key: 633 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 634 additive, full, kaskade 635 636 Level: advanced 637 638 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 639 640 .seealso: PCMGSetLevels() 641 @*/ 642 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form) 643 { 644 PC_MG **mg; 645 646 PetscFunctionBegin; 647 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 648 mg = (PC_MG**)pc->data; 649 650 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 651 mg[0]->am = form; 652 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 653 else pc->ops->applyrichardson = 0; 654 PetscFunctionReturn(0); 655 } 656 657 #undef __FUNCT__ 658 #define __FUNCT__ "PCMGSetCycleType" 659 /*@ 660 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 661 complicated cycling. 662 663 Collective on PC 664 665 Input Parameters: 666 + pc - the multigrid context 667 - PC_MG_CYCLE_V or PC_MG_CYCLE_W 668 669 Options Database Key: 670 $ -pc_mg_cycle_type v or w 671 672 Level: advanced 673 674 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 675 676 .seealso: PCMGSetCycleTypeOnLevel() 677 @*/ 678 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n) 679 { 680 PC_MG **mg; 681 PetscInt i,levels; 682 683 PetscFunctionBegin; 684 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 685 mg = (PC_MG**)pc->data; 686 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 687 levels = mg[0]->levels; 688 689 for (i=0; i<levels; i++) { 690 mg[i]->cycles = n; 691 } 692 PetscFunctionReturn(0); 693 } 694 695 #undef __FUNCT__ 696 #define __FUNCT__ "PCMGMultiplicativeSetCycles" 697 /*@ 698 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 699 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 700 701 Collective on PC 702 703 Input Parameters: 704 + pc - the multigrid context 705 - n - number of cycles (default is 1) 706 707 Options Database Key: 708 $ -pc_mg_multiplicative_cycles n 709 710 Level: advanced 711 712 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 713 714 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 715 716 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 717 @*/ 718 PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 719 { 720 PC_MG **mg; 721 PetscInt i,levels; 722 723 PetscFunctionBegin; 724 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 725 mg = (PC_MG**)pc->data; 726 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 727 levels = mg[0]->levels; 728 729 for (i=0; i<levels; i++) { 730 mg[i]->cyclesperpcapply = n; 731 } 732 PetscFunctionReturn(0); 733 } 734 735 #undef __FUNCT__ 736 #define __FUNCT__ "PCMGSetGalerkin" 737 /*@ 738 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 739 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 740 741 Collective on PC 742 743 Input Parameters: 744 . pc - the multigrid context 745 746 Options Database Key: 747 $ -pc_mg_galerkin 748 749 Level: intermediate 750 751 .keywords: MG, set, Galerkin 752 753 .seealso: PCMGGetGalerkin() 754 755 @*/ 756 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc) 757 { 758 PC_MG **mg; 759 PetscInt i,levels; 760 761 PetscFunctionBegin; 762 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 763 mg = (PC_MG**)pc->data; 764 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 765 levels = mg[0]->levels; 766 767 for (i=0; i<levels; i++) { 768 mg[i]->galerkin = PETSC_TRUE; 769 } 770 PetscFunctionReturn(0); 771 } 772 773 #undef __FUNCT__ 774 #define __FUNCT__ "PCMGGetGalerkin" 775 /*@ 776 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 777 A_i-1 = r_i * A_i * r_i^t 778 779 Not Collective 780 781 Input Parameter: 782 . pc - the multigrid context 783 784 Output Parameter: 785 . gelerkin - PETSC_TRUE or PETSC_FALSE 786 787 Options Database Key: 788 $ -pc_mg_galerkin 789 790 Level: intermediate 791 792 .keywords: MG, set, Galerkin 793 794 .seealso: PCMGSetGalerkin() 795 796 @*/ 797 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin) 798 { 799 PC_MG **mg; 800 801 PetscFunctionBegin; 802 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 803 mg = (PC_MG**)pc->data; 804 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 805 *galerkin = mg[0]->galerkin; 806 PetscFunctionReturn(0); 807 } 808 809 #undef __FUNCT__ 810 #define __FUNCT__ "PCMGSetNumberSmoothDown" 811 /*@ 812 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 813 use on all levels. Use PCMGGetSmootherDown() to set different 814 pre-smoothing steps on different levels. 815 816 Collective on PC 817 818 Input Parameters: 819 + mg - the multigrid context 820 - n - the number of smoothing steps 821 822 Options Database Key: 823 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 824 825 Level: advanced 826 827 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 828 829 .seealso: PCMGSetNumberSmoothUp() 830 @*/ 831 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n) 832 { 833 PC_MG **mg; 834 PetscErrorCode ierr; 835 PetscInt i,levels; 836 837 PetscFunctionBegin; 838 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 839 mg = (PC_MG**)pc->data; 840 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 841 levels = mg[0]->levels; 842 843 for (i=1; i<levels; i++) { 844 /* make sure smoother up and down are different */ 845 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 846 ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 847 mg[i]->default_smoothd = n; 848 } 849 PetscFunctionReturn(0); 850 } 851 852 #undef __FUNCT__ 853 #define __FUNCT__ "PCMGSetNumberSmoothUp" 854 /*@ 855 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 856 on all levels. Use PCMGGetSmootherUp() to set different numbers of 857 post-smoothing steps on different levels. 858 859 Collective on PC 860 861 Input Parameters: 862 + mg - the multigrid context 863 - n - the number of smoothing steps 864 865 Options Database Key: 866 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 867 868 Level: advanced 869 870 Note: this does not set a value on the coarsest grid, since we assume that 871 there is no separate smooth up on the coarsest grid. 872 873 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 874 875 .seealso: PCMGSetNumberSmoothDown() 876 @*/ 877 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n) 878 { 879 PC_MG **mg; 880 PetscErrorCode ierr; 881 PetscInt i,levels; 882 883 PetscFunctionBegin; 884 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 885 mg = (PC_MG**)pc->data; 886 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 887 levels = mg[0]->levels; 888 889 for (i=1; i<levels; i++) { 890 /* make sure smoother up and down are different */ 891 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 892 ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 893 mg[i]->default_smoothu = n; 894 } 895 PetscFunctionReturn(0); 896 } 897 898 /* ----------------------------------------------------------------------------------------*/ 899 900 /*MC 901 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 902 information about the coarser grid matrices and restriction/interpolation operators. 903 904 Options Database Keys: 905 + -pc_mg_levels <nlevels> - number of levels including finest 906 . -pc_mg_cycles v or w 907 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 908 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 909 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 910 . -pc_mg_log - log information about time spent on each level of the solver 911 . -pc_mg_monitor - print information on the multigrid convergence 912 . -pc_mg_galerkin - use Galerkin process to compute coarser operators 913 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 914 to the Socket viewer for reading from Matlab. 915 916 Notes: 917 918 Level: intermediate 919 920 Concepts: multigrid 921 922 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 923 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 924 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 925 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 926 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 927 M*/ 928 929 EXTERN_C_BEGIN 930 #undef __FUNCT__ 931 #define __FUNCT__ "PCCreate_MG" 932 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 933 { 934 PetscFunctionBegin; 935 pc->ops->apply = PCApply_MG; 936 pc->ops->setup = PCSetUp_MG; 937 pc->ops->destroy = PCDestroy_MG; 938 pc->ops->setfromoptions = PCSetFromOptions_MG; 939 pc->ops->view = PCView_MG; 940 941 pc->data = (void*)0; 942 PetscFunctionReturn(0); 943 } 944 EXTERN_C_END 945