Actual source code: arkimex.c

petsc-3.4.5 2014-06-29
  1: /*
  2:   Code for timestepping with additive Runge-Kutta IMEX method

  4:   Notes:
  5:   The general system is written as

  7:   F(t,U,Udot) = G(t,U)

  9:   where F represents the stiff part of the physics and G represents the non-stiff part.

 11: */
 12: #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
 13: #include <petscdm.h>

 15: static TSARKIMEXType TSARKIMEXDefault = TSARKIMEX3;
 16: static PetscBool     TSARKIMEXRegisterAllCalled;
 17: static PetscBool     TSARKIMEXPackageInitialized;
 18: static PetscInt      explicit_stage_time_id;

 20: typedef struct _ARKTableau *ARKTableau;
 21: struct _ARKTableau {
 22:   char      *name;
 23:   PetscInt  order;                /* Classical approximation order of the method */
 24:   PetscInt  s;                    /* Number of stages */
 25:   PetscBool stiffly_accurate;     /* The implicit part is stiffly accurate*/
 26:   PetscBool FSAL_implicit;        /* The implicit part is FSAL*/
 27:   PetscBool explicit_first_stage; /* The implicit part has an explicit first stage*/
 28:   PetscInt  pinterp;              /* Interpolation order */
 29:   PetscReal *At,*bt,*ct;          /* Stiff tableau */
 30:   PetscReal *A,*b,*c;             /* Non-stiff tableau */
 31:   PetscReal *bembedt,*bembed;     /* Embedded formula of order one less (order-1) */
 32:   PetscReal *binterpt,*binterp;   /* Dense output formula */
 33:   PetscReal ccfl;                 /* Placeholder for CFL coefficient relative to forward Euler */
 34: };
 35: