1 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2 3 static const SNESMSType SNESMSDefault = SNESMSM62; 4 static PetscBool SNESMSRegisterAllCalled; 5 static PetscBool SNESMSPackageInitialized; 6 7 typedef struct _SNESMSTableau *SNESMSTableau; 8 struct _SNESMSTableau { 9 char *name; 10 PetscInt nstages; /* Number of stages */ 11 PetscInt nregisters; /* Number of registers */ 12 PetscReal stability; /* Scaled stability region */ 13 PetscReal *gamma; /* Coefficients of 3S* method */ 14 PetscReal *delta; /* Coefficients of 3S* method */ 15 PetscReal *betasub; /* Subdiagonal of beta in Shu-Osher form */ 16 }; 17 18 typedef struct _SNESMSTableauLink *SNESMSTableauLink; 19 struct _SNESMSTableauLink { 20 struct _SNESMSTableau tab; 21 SNESMSTableauLink next; 22 }; 23 static SNESMSTableauLink SNESMSTableauList; 24 25 typedef struct { 26 SNESMSTableau tableau; /* Tableau in low-storage form */ 27 PetscReal damping; /* Damping parameter, like length of (pseudo) time step */ 28 PetscBool norms; /* Compute norms, usually only for monitoring purposes */ 29 } SNES_MS; 30 31 #undef __FUNCT__ 32 #define __FUNCT__ "SNESMSRegisterAll" 33 /*@C 34 SNESMSRegisterAll - Registers all of the multi-stage methods in SNESMS 35 36 Not Collective, but should be called by all processes which will need the schemes to be registered 37 38 Level: advanced 39 40 .keywords: SNES, SNESMS, register, all 41 42 .seealso: SNESMSRegisterDestroy() 43 @*/ 44 PetscErrorCode SNESMSRegisterAll(void) 45 { 46 PetscErrorCode ierr; 47 48 PetscFunctionBegin; 49 if (SNESMSRegisterAllCalled) PetscFunctionReturn(0); 50 SNESMSRegisterAllCalled = PETSC_TRUE; 51 52 { 53 PetscReal gamma[3][1] = {{1.0},{0.0},{0.0}}; 54 PetscReal delta[1] = {0.0}; 55 PetscReal betasub[1] = {1.0}; 56 ierr = SNESMSRegister(SNESMSEULER,1,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 57 } 58 59 { 60 PetscReal gamma[3][6] = { 61 {0.0000000000000000E+00, -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02, -1.4204296130641869E-01}, 62 {1.0000000000000000E+00, 1.1111025767083920E+00, 5.6150921583923230E-01, 7.4151723494934041E-01, 3.1714538168600587E-01, 4.6479276238548706E-01}, 63 {0.0000000000000000E+00, 0.0000000000000000E+00, 0.0000000000000000E+00, 6.7968174970583317E-01, -4.1755042846051737E-03, -1.9115668129923846E-01} 64 }; 65 PetscReal delta[6] = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00}; 66 PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01}; 67 ierr = SNESMSRegister(SNESMSM62,6,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 68 } 69 { 70 PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}}; 71 PetscReal delta[4] = {0,0,0,0}; 72 PetscReal betasub[4] = {0.25, 0.5, 0.55, 1.0}; 73 ierr = SNESMSRegister(SNESMSJAMESON83,4,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 74 } 75 { /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */ 76 PetscReal gamma[3][2] = {{0,0},{0,0},{1,1}}; 77 PetscReal delta[2] = {0,0}; 78 PetscReal betasub[2] = {0.3333,1.0}; 79 ierr = SNESMSRegister(SNESMSVLTP21,2,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 80 } 81 { /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */ 82 PetscReal gamma[3][3] = {{0,0,0},{0,0,0},{1,1,1}}; 83 PetscReal delta[3] = {0,0,0}; 84 PetscReal betasub[3] = {0.1481,0.4000,1.0}; 85 ierr = SNESMSRegister(SNESMSVLTP31,3,3,1.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 86 } 87 { /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */ 88 PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}}; 89 PetscReal delta[4] = {0,0,0,0}; 90 PetscReal betasub[4] = {0.0833,0.2069,0.4265,1.0}; 91 ierr = SNESMSRegister(SNESMSVLTP41,4,3,2.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 92 } 93 { /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */ 94 PetscReal gamma[3][5] = {{0,0,0,0,0},{0,0,0,0,0},{1,1,1,1,1}}; 95 PetscReal delta[5] = {0,0,0,0,0}; 96 PetscReal betasub[5] = {0.0533,0.1263,0.2375,0.4414,1.0}; 97 ierr = SNESMSRegister(SNESMSVLTP51,5,3,2.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 98 } 99 { /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */ 100 PetscReal gamma[3][6] = {{0,0,0,0,0,0},{0,0,0,0,0,0},{1,1,1,1,1,1}}; 101 PetscReal delta[6] = {0,0,0,0,0,0}; 102 PetscReal betasub[6] = {0.0370,0.0851,0.1521,0.2562,0.4512,1.0}; 103 ierr = SNESMSRegister(SNESMSVLTP61,6,3,3.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 104 } 105 PetscFunctionReturn(0); 106 } 107 108 #undef __FUNCT__ 109 #define __FUNCT__ "SNESMSRegisterDestroy" 110 /*@C 111 SNESMSRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister(). 112 113 Not Collective 114 115 Level: advanced 116 117 .keywords: TSRosW, register, destroy 118 .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic() 119 @*/ 120 PetscErrorCode SNESMSRegisterDestroy(void) 121 { 122 PetscErrorCode ierr; 123 SNESMSTableauLink link; 124 125 PetscFunctionBegin; 126 while ((link = SNESMSTableauList)) { 127 SNESMSTableau t = &link->tab; 128 SNESMSTableauList = link->next; 129 ierr = PetscFree3(t->gamma,t->delta,t->betasub);CHKERRQ(ierr); 130 ierr = PetscFree(t->name);CHKERRQ(ierr); 131 ierr = PetscFree(link);CHKERRQ(ierr); 132 } 133 SNESMSRegisterAllCalled = PETSC_FALSE; 134 PetscFunctionReturn(0); 135 } 136 137 #undef __FUNCT__ 138 #define __FUNCT__ "SNESMSInitializePackage" 139 /*@C 140 SNESMSInitializePackage - This function initializes everything in the SNESMS package. It is called 141 from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to SNESCreate_MS() 142 when using static libraries. 143 144 Input Parameter: 145 path - The dynamic library path, or PETSC_NULL 146 147 Level: developer 148 149 .keywords: SNES, SNESMS, initialize, package 150 .seealso: PetscInitialize() 151 @*/ 152 PetscErrorCode SNESMSInitializePackage(const char path[]) 153 { 154 PetscErrorCode ierr; 155 156 PetscFunctionBegin; 157 if (SNESMSPackageInitialized) PetscFunctionReturn(0); 158 SNESMSPackageInitialized = PETSC_TRUE; 159 ierr = SNESMSRegisterAll();CHKERRQ(ierr); 160 ierr = PetscRegisterFinalize(SNESMSFinalizePackage);CHKERRQ(ierr); 161 PetscFunctionReturn(0); 162 } 163 164 #undef __FUNCT__ 165 #define __FUNCT__ "SNESMSFinalizePackage" 166 /*@C 167 SNESMSFinalizePackage - This function destroys everything in the SNESMS package. It is 168 called from PetscFinalize(). 169 170 Level: developer 171 172 .keywords: Petsc, destroy, package 173 .seealso: PetscFinalize() 174 @*/ 175 PetscErrorCode SNESMSFinalizePackage(void) 176 { 177 PetscErrorCode ierr; 178 179 PetscFunctionBegin; 180 SNESMSPackageInitialized = PETSC_FALSE; 181 ierr = SNESMSRegisterDestroy();CHKERRQ(ierr); 182 PetscFunctionReturn(0); 183 } 184 185 #undef __FUNCT__ 186 #define __FUNCT__ "SNESMSRegister" 187 /*@C 188 SNESMSRegister - register a multistage scheme 189 190 Not Collective, but the same schemes should be registered on all processes on which they will be used 191 192 Input Parameters: 193 + name - identifier for method 194 . nstages - number of stages 195 . nregisters - number of registers used by low-storage implementation 196 . gamma - coefficients, see Ketcheson's paper 197 . delta - coefficients, see Ketcheson's paper 198 - betasub - subdiagonal of Shu-Osher form 199 200 Notes: 201 The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 202 203 Level: advanced 204 205 .keywords: SNES, register 206 207 .seealso: SNESMS 208 @*/ 209 PetscErrorCode SNESMSRegister(const SNESMSType name,PetscInt nstages,PetscInt nregisters,PetscReal stability,const PetscReal gamma[],const PetscReal delta[],const PetscReal betasub[]) 210 { 211 PetscErrorCode ierr; 212 SNESMSTableauLink link; 213 SNESMSTableau t; 214 215 PetscFunctionBegin; 216 PetscValidCharPointer(name,1); 217 if (nstages < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have at least one stage"); 218 if (nregisters != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 3-register form"); 219 PetscValidPointer(gamma,4); 220 PetscValidPointer(delta,5); 221 PetscValidPointer(betasub,6); 222 223 ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 224 ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 225 t = &link->tab; 226 ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 227 t->nstages = nstages; 228 t->nregisters = nregisters; 229 t->stability = stability; 230 ierr = PetscMalloc3(nstages*nregisters,PetscReal,&t->gamma,nstages,PetscReal,&t->delta,nstages,PetscReal,&t->betasub);CHKERRQ(ierr); 231 ierr = PetscMemcpy(t->gamma,gamma,nstages*nregisters*sizeof(PetscReal));CHKERRQ(ierr); 232 ierr = PetscMemcpy(t->delta,delta,nstages*sizeof(PetscReal));CHKERRQ(ierr); 233 ierr = PetscMemcpy(t->betasub,betasub,nstages*sizeof(PetscReal));CHKERRQ(ierr); 234 235 link->next = SNESMSTableauList; 236 SNESMSTableauList = link; 237 PetscFunctionReturn(0); 238 } 239 240 #undef __FUNCT__ 241 #define __FUNCT__ "SNESMSStep_3Sstar" 242 /* 243 X - initial state, updated in-place. 244 F - residual, computed at the initial X on input 245 */ 246 static PetscErrorCode SNESMSStep_3Sstar(SNES snes,Vec X,Vec F) 247 { 248 PetscErrorCode ierr; 249 SNES_MS *ms = (SNES_MS*)snes->data; 250 SNESMSTableau t = ms->tableau; 251 const PetscReal *gamma = t->gamma,*delta = t->delta,*betasub = t->betasub; 252 Vec S1,S2,S3,Y; 253 PetscInt i,nstages = t->nstages;; 254 255 256 PetscFunctionBegin; 257 Y = snes->work[0]; 258 S1 = X; 259 S2 = snes->work[1]; 260 S3 = snes->work[2]; 261 ierr = VecZeroEntries(S2);CHKERRQ(ierr); 262 ierr = VecCopy(X,S3);CHKERRQ(ierr); 263 for (i=0; i<nstages; i++) { 264 Vec Ss[4] = {S1,S2,S3,Y}; 265 PetscScalar scoeff[4] = {gamma[0*nstages+i]-1.0,gamma[1*nstages+i],gamma[2*nstages+i],-betasub[i]*ms->damping}; 266 ierr = VecAXPY(S2,delta[i],S1);CHKERRQ(ierr); 267 if (i>0) { 268 ierr = SNESComputeFunction(snes,S1,F);CHKERRQ(ierr); 269 if (snes->domainerror) { 270 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 271 PetscFunctionReturn(0); 272 } 273 } 274 ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 275 ierr = VecMAXPY(S1,4,scoeff,Ss);CHKERRQ(ierr); 276 } 277 PetscFunctionReturn(0); 278 } 279 280 #undef __FUNCT__ 281 #define __FUNCT__ "SNESSolve_MS" 282 static PetscErrorCode SNESSolve_MS(SNES snes) 283 { 284 SNES_MS *ms = (SNES_MS*)snes->data; 285 Vec X = snes->vec_sol,F = snes->vec_func; 286 PetscReal fnorm; 287 MatStructure mstruct; 288 PetscInt i; 289 PetscErrorCode ierr; 290 291 PetscFunctionBegin; 292 snes->reason = SNES_CONVERGED_ITERATING; 293 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 294 snes->iter = 0; 295 snes->norm = 0.; 296 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 297 if (!snes->vec_func_init_set) { 298 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 299 if (snes->domainerror) { 300 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 301 PetscFunctionReturn(0); 302 } 303 } else { 304 snes->vec_func_init_set = PETSC_FALSE; 305 } 306 if (snes->jacobian) { /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */ 307 ierr = SNESComputeJacobian(snes,snes->vec_sol,&snes->jacobian,&snes->jacobian_pre,&mstruct);CHKERRQ(ierr); 308 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,mstruct);CHKERRQ(ierr); 309 } 310 if (ms->norms) { 311 if (!snes->norm_init_set) { 312 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 313 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 314 } else { 315 fnorm = snes->norm_init; 316 snes->norm_init_set = PETSC_FALSE; 317 } 318 /* Monitor convergence */ 319 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 320 snes->iter = 0; 321 snes->norm = fnorm; 322 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 323 SNESLogConvHistory(snes,snes->norm,0); 324 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 325 326 /* set parameter for default relative tolerance convergence test */ 327 snes->ttol = fnorm*snes->rtol; 328 /* Test for convergence */ 329 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 330 if (snes->reason) PetscFunctionReturn(0); 331 } 332 333 /* Call general purpose update function */ 334 if (snes->ops->update) { 335 ierr = (*snes->ops->update)(snes,snes->iter);CHKERRQ(ierr); 336 } 337 for (i = 0; i < snes->max_its; i++) { 338 ierr = SNESMSStep_3Sstar(snes,X,F);CHKERRQ(ierr); 339 340 if (i+1 < snes->max_its || ms->norms) { 341 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 342 if (snes->domainerror) { 343 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 344 PetscFunctionReturn(0); 345 } 346 } 347 348 if (ms->norms) { 349 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 350 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 351 /* Monitor convergence */ 352 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 353 snes->iter = i+1; 354 snes->norm = fnorm; 355 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 356 SNESLogConvHistory(snes,snes->norm,0); 357 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 358 359 /* Test for convergence */ 360 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 361 if (snes->reason) PetscFunctionReturn(0); 362 } 363 364 /* Call general purpose update function */ 365 if (snes->ops->update) { 366 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 367 } 368 } 369 if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 370 PetscFunctionReturn(0); 371 } 372 373 #undef __FUNCT__ 374 #define __FUNCT__ "SNESSetUp_MS" 375 static PetscErrorCode SNESSetUp_MS(SNES snes) 376 { 377 SNES_MS * ms = (SNES_MS *)snes->data; 378 PetscErrorCode ierr; 379 380 PetscFunctionBegin; 381 if (!ms->tableau) {ierr = SNESMSSetType(snes,SNESMSDefault);CHKERRQ(ierr);} 382 ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr); 383 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 384 PetscFunctionReturn(0); 385 } 386 387 #undef __FUNCT__ 388 #define __FUNCT__ "SNESReset_MS" 389 static PetscErrorCode SNESReset_MS(SNES snes) 390 { 391 392 PetscFunctionBegin; 393 PetscFunctionReturn(0); 394 } 395 396 #undef __FUNCT__ 397 #define __FUNCT__ "SNESDestroy_MS" 398 static PetscErrorCode SNESDestroy_MS(SNES snes) 399 { 400 PetscErrorCode ierr; 401 402 PetscFunctionBegin; 403 ierr = PetscFree(snes->data);CHKERRQ(ierr); 404 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 405 PetscFunctionReturn(0); 406 } 407 408 #undef __FUNCT__ 409 #define __FUNCT__ "SNESView_MS" 410 static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer) 411 { 412 PetscBool iascii; 413 PetscErrorCode ierr; 414 SNES_MS *ms = (SNES_MS*)snes->data; 415 416 PetscFunctionBegin; 417 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 418 if (iascii) { 419 SNESMSTableau tab = ms->tableau; 420 ierr = PetscViewerASCIIPrintf(viewer," multi-stage method type: %s\n",tab?tab->name:"not yet set");CHKERRQ(ierr); 421 } 422 PetscFunctionReturn(0); 423 } 424 425 #undef __FUNCT__ 426 #define __FUNCT__ "SNESSetFromOptions_MS" 427 static PetscErrorCode SNESSetFromOptions_MS(SNES snes) 428 { 429 SNES_MS *ms = (SNES_MS*)snes->data; 430 PetscErrorCode ierr; 431 432 PetscFunctionBegin; 433 ierr = PetscOptionsHead("SNES MS options");CHKERRQ(ierr); 434 { 435 SNESMSTableauLink link; 436 PetscInt count,choice; 437 PetscBool flg; 438 const char **namelist; 439 char mstype[256]; 440 441 ierr = PetscStrncpy(mstype,SNESMSDefault,sizeof(mstype));CHKERRQ(ierr); 442 for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ; 443 ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 444 for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 445 ierr = PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg);CHKERRQ(ierr); 446 ierr = SNESMSSetType(snes,flg ? namelist[choice] : mstype);CHKERRQ(ierr); 447 ierr = PetscFree(namelist);CHKERRQ(ierr); 448 ierr = PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",ms->damping,&ms->damping,PETSC_NULL);CHKERRQ(ierr); 449 ierr = PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,PETSC_NULL);CHKERRQ(ierr); 450 } 451 ierr = PetscOptionsTail();CHKERRQ(ierr); 452 PetscFunctionReturn(0); 453 } 454 455 EXTERN_C_BEGIN 456 #undef __FUNCT__ 457 #define __FUNCT__ "SNESMSSetType_MS" 458 PetscErrorCode SNESMSSetType_MS(SNES snes,SNESMSType mstype) 459 { 460 PetscErrorCode ierr; 461 SNES_MS *ms = (SNES_MS*)snes->data; 462 SNESMSTableauLink link; 463 PetscBool match; 464 465 PetscFunctionBegin; 466 if (ms->tableau) { 467 ierr = PetscStrcmp(ms->tableau->name,mstype,&match);CHKERRQ(ierr); 468 if (match) PetscFunctionReturn(0); 469 } 470 for (link = SNESMSTableauList; link; link=link->next) { 471 ierr = PetscStrcmp(link->tab.name,mstype,&match);CHKERRQ(ierr); 472 if (match) { 473 ierr = SNESReset_MS(snes);CHKERRQ(ierr); 474 ms->tableau = &link->tab; 475 PetscFunctionReturn(0); 476 } 477 } 478 SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype); 479 PetscFunctionReturn(0); 480 } 481 EXTERN_C_END 482 483 484 #undef __FUNCT__ 485 #define __FUNCT__ "SNESMSSetType" 486 /*@C 487 SNESMSSetType - Set the type of multistage smoother 488 489 Logically collective 490 491 Input Parameter: 492 + snes - nonlinear solver context 493 - mstype - type of multistage method 494 495 Level: beginner 496 497 .seealso: SNESMSGetType(), SNESMS 498 @*/ 499 PetscErrorCode SNESMSSetType(SNES snes,const SNESMSType rostype) 500 { 501 PetscErrorCode ierr; 502 503 PetscFunctionBegin; 504 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 505 ierr = PetscTryMethod(snes,"SNESMSSetType_C",(SNES,const SNESMSType),(snes,rostype));CHKERRQ(ierr); 506 PetscFunctionReturn(0); 507 } 508 509 /* -------------------------------------------------------------------------- */ 510 /*MC 511 SNESMS - multi-stage smoothers 512 513 Options Database: 514 515 + -snes_ms_type - type of multi-stage smoother 516 - -snes_ms_damping - damping for multi-stage method 517 518 Notes: 519 These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 520 FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev). 521 522 Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 523 524 The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister(). 525 526 References: 527 528 Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 529 530 Jameson (1983) Solution of the Euler equations for two dimensional transonic flow by a multigrid method. 531 532 Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes. 533 534 Level: beginner 535 536 .seealso: SNESCreate(), SNES, SNESSetType(), SNESMS, SNESFAS, KSPCHEBYSHEV 537 538 M*/ 539 EXTERN_C_BEGIN 540 #undef __FUNCT__ 541 #define __FUNCT__ "SNESCreate_MS" 542 PetscErrorCode SNESCreate_MS(SNES snes) 543 { 544 PetscErrorCode ierr; 545 SNES_MS *ms; 546 547 PetscFunctionBegin; 548 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 549 ierr = SNESMSInitializePackage(PETSC_NULL);CHKERRQ(ierr); 550 #endif 551 552 snes->ops->setup = SNESSetUp_MS; 553 snes->ops->solve = SNESSolve_MS; 554 snes->ops->destroy = SNESDestroy_MS; 555 snes->ops->setfromoptions = SNESSetFromOptions_MS; 556 snes->ops->view = SNESView_MS; 557 snes->ops->reset = SNESReset_MS; 558 559 snes->usespc = PETSC_FALSE; 560 snes->usesksp = PETSC_TRUE; 561 562 ierr = PetscNewLog(snes,SNES_MS,&ms);CHKERRQ(ierr); 563 snes->data = (void*)ms; 564 ms->damping = 0.9; 565 ms->norms = PETSC_FALSE; 566 567 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESMSSetType_C","SNESMSSetType_MS",SNESMSSetType_MS);CHKERRQ(ierr); 568 PetscFunctionReturn(0); 569 } 570 EXTERN_C_END 571