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