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