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