xref: /petsc/src/ts/impls/arkimex/arkimex.c (revision 4f38528163c37a794bf829ba4a4ecb99bdc1de85)
18a381b04SJed Brown /*
28a381b04SJed Brown   Code for timestepping with additive Runge-Kutta IMEX method
38a381b04SJed Brown 
48a381b04SJed Brown   Notes:
58a381b04SJed Brown   The general system is written as
68a381b04SJed Brown 
78a381b04SJed Brown   F(t,X,Xdot) = G(t,X)
88a381b04SJed Brown 
98a381b04SJed Brown   where F represents the stiff part of the physics and G represents the non-stiff part.
108a381b04SJed Brown 
118a381b04SJed Brown */
128a381b04SJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
138a381b04SJed Brown 
14*4f385281SJed Brown static const TSARKIMEXType TSARKIMEXDefault = TSARKIMEX2E;
158a381b04SJed Brown static PetscBool TSARKIMEXRegisterAllCalled;
168a381b04SJed Brown static PetscBool TSARKIMEXPackageInitialized;
178a381b04SJed Brown 
188a381b04SJed Brown typedef struct _ARKTableau *ARKTableau;
198a381b04SJed Brown struct _ARKTableau {
208a381b04SJed Brown   char *name;
21*4f385281SJed Brown   PetscInt order;               /* Classical approximation order of the method */
22*4f385281SJed Brown   PetscInt s;                   /* Number of stages */
23*4f385281SJed Brown   PetscInt pinterp;             /* Interpolation order */
24*4f385281SJed Brown   PetscReal *At,*bt,*ct;        /* Stiff tableau */
258a381b04SJed Brown   PetscReal *A,*b,*c;           /* Non-stiff tableau */
26cd652676SJed Brown   PetscReal *binterpt,*binterp; /* Dense output formula */
278a381b04SJed Brown };
288a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink;
298a381b04SJed Brown struct _ARKTableauLink {
308a381b04SJed Brown   struct _ARKTableau tab;
318a381b04SJed Brown   ARKTableauLink next;
328a381b04SJed Brown };
338a381b04SJed Brown static ARKTableauLink ARKTableauList;
348a381b04SJed Brown 
358a381b04SJed Brown typedef struct {
368a381b04SJed Brown   ARKTableau  tableau;
378a381b04SJed Brown   Vec         *Y;               /* States computed during the step */
388a381b04SJed Brown   Vec         *YdotI;           /* Time derivatives for the stiff part */
398a381b04SJed Brown   Vec         *YdotRHS;         /* Function evaluations for the non-stiff part */
408a381b04SJed Brown   Vec         Ydot;             /* Work vector holding Ydot during residual evaluation */
418a381b04SJed Brown   Vec         Work;             /* Generic work vector */
428a381b04SJed Brown   Vec         Z;                /* Ydot = shift(Y-Z) */
438a381b04SJed Brown   PetscScalar *work;            /* Scalar work */
448a381b04SJed Brown   PetscReal   shift;
458a381b04SJed Brown   PetscReal   stage_time;
468a381b04SJed Brown } TS_ARKIMEX;
478a381b04SJed Brown 
488a381b04SJed Brown #undef __FUNCT__
498a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterAll"
508a381b04SJed Brown /*@C
518a381b04SJed Brown   TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX
528a381b04SJed Brown 
53fca742c7SJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
548a381b04SJed Brown 
558a381b04SJed Brown   Level: advanced
568a381b04SJed Brown 
578a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all
588a381b04SJed Brown 
598a381b04SJed Brown .seealso:  TSARKIMEXRegisterDestroy()
608a381b04SJed Brown @*/
618a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void)
628a381b04SJed Brown {
638a381b04SJed Brown   PetscErrorCode ierr;
648a381b04SJed Brown 
658a381b04SJed Brown   PetscFunctionBegin;
668a381b04SJed Brown   if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0);
678a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_TRUE;
688a381b04SJed Brown 
698a381b04SJed Brown   {
708a381b04SJed Brown     const PetscReal
718a381b04SJed Brown       A[3][3] = {{0,0,0},
728a381b04SJed Brown                  {0.41421356237309504880,0,0},
738a381b04SJed Brown                  {0.75,0.25,0}},
748a381b04SJed Brown       At[3][3] = {{0,0,0},
758a381b04SJed Brown                   {0.12132034355964257320,0.29289321881345247560,0},
768a381b04SJed Brown                   {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}};
77cd652676SJed Brown     ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,0,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
788a381b04SJed Brown   }
79a3a57f36SJed Brown   {
80a3a57f36SJed Brown     const PetscReal s2 = sqrt(2),
81a3a57f36SJed Brown       A[3][3] = {{0,0,0},
82a3a57f36SJed Brown                  {2-s2,0,0},
83a3a57f36SJed Brown                  {(3-2*s2)/6,(3+2*s2)/6,0}},
84a3a57f36SJed Brown       At[3][3] = {{0,0,0},
85a3a57f36SJed Brown                   {1-1/s2,1-1/s2,0},
86cd652676SJed Brown                   {1/(2*s2),1/(2*s2),1-1/s2}},
87cd652676SJed Brown       binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}};
88cd652676SJed Brown     ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
89a3a57f36SJed Brown   }
90a3a57f36SJed Brown   {
91a3a57f36SJed Brown     const PetscReal
92a3a57f36SJed Brown       A[4][4] = {{0,0,0,0},
934040e9f2SJed Brown                  {1767732205903./2027836641118.,0,0,0},
944040e9f2SJed Brown                  {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
954040e9f2SJed Brown                  {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
96a3a57f36SJed Brown       At[4][4] = {{0,0,0,0},
974040e9f2SJed Brown                   {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
984040e9f2SJed Brown                   {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
994040e9f2SJed Brown                   {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
1004040e9f2SJed Brown       binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
1014040e9f2SJed Brown                         {-18682724506714./9892148508045.,17870216137069./13817060693119.},
1024040e9f2SJed Brown                         {34259539580243./13192909600954.,-28141676662227./17317692491321.},
1034040e9f2SJed Brown                         {584795268549./6622622206610.,   2508943948391./7218656332882.}};
104cd652676SJed Brown     ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
105a3a57f36SJed Brown   }
106a3a57f36SJed Brown   {
107a3a57f36SJed Brown     const PetscReal
108a3a57f36SJed Brown       A[6][6] = {{0,0,0,0,0,0},
109a3a57f36SJed Brown                  {1./2,0,0,0,0,0},
1104040e9f2SJed Brown                  {13861./62500.,6889./62500.,0,0,0,0},
1114040e9f2SJed Brown                  {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
1124040e9f2SJed Brown                  {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
1134040e9f2SJed Brown                  {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
114a3a57f36SJed Brown       At[6][6] = {{0,0,0,0,0,0},
115a3a57f36SJed Brown                   {1./4,1./4,0,0,0,0},
1164040e9f2SJed Brown                   {8611./62500.,-1743./31250.,1./4,0,0,0},
1174040e9f2SJed Brown                   {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
1184040e9f2SJed Brown                   {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
1194040e9f2SJed Brown                   {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
1204040e9f2SJed Brown       binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
121cd652676SJed Brown                         {0,0,0},
1224040e9f2SJed Brown                         {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
1234040e9f2SJed Brown                         {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
1244040e9f2SJed Brown                         {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
1254040e9f2SJed Brown                         {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
126cd652676SJed Brown     ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
127a3a57f36SJed Brown   }
128a3a57f36SJed Brown   {
129a3a57f36SJed Brown     const PetscReal
130a3a57f36SJed Brown       A[8][8] = {{0,0,0,0,0,0,0,0},
131a3a57f36SJed Brown                  {41./100,0,0,0,0,0,0,0},
1324040e9f2SJed Brown                  {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
1334040e9f2SJed Brown                  {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
1344040e9f2SJed Brown                  {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
1354040e9f2SJed Brown                  {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
1364040e9f2SJed Brown                  {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
1374040e9f2SJed Brown                  {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
138a3a57f36SJed Brown       At[8][8] = {{0,0,0,0,0,0,0,0},
1394040e9f2SJed Brown                   {41./200.,41./200.,0,0,0,0,0,0},
1404040e9f2SJed Brown                   {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
1414040e9f2SJed Brown                   {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
1424040e9f2SJed Brown                   {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
1434040e9f2SJed Brown                   {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
1444040e9f2SJed Brown                   {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
1454040e9f2SJed Brown                   {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
1464040e9f2SJed Brown       binterpt[8][3] = {{-17674230611817./10670229744614. ,  43486358583215./12773830924787. , -9257016797708./5021505065439.},
147cd652676SJed Brown                         {0                               ,  0                              , 0                            },
148cd652676SJed Brown                         {0                               ,  0                              , 0                            },
1494040e9f2SJed Brown                         {65168852399939./7868540260826.  ,  -91478233927265./11067650958493., 26096422576131./11239449250142.},
1504040e9f2SJed Brown                         {15494834004392./5936557850923.  ,  -79368583304911./10890268929626., 92396832856987./20362823103730.},
1514040e9f2SJed Brown                         {-99329723586156./26959484932159.,  -12239297817655./9152339842473. , 30029262896817./10175596800299.},
1524040e9f2SJed Brown                         {-19024464361622./5461577185407. ,  115839755401235./10719374521269., -26136350496073./3983972220547.},
1534040e9f2SJed Brown                         {-6511271360970./6095937251113.  ,  5843115559534./2180450260947.   , -5289405421727./3760307252460. }};
154cd652676SJed Brown     ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
155a3a57f36SJed Brown   }
156a3a57f36SJed Brown 
1578a381b04SJed Brown   PetscFunctionReturn(0);
1588a381b04SJed Brown }
1598a381b04SJed Brown 
1608a381b04SJed Brown #undef __FUNCT__
1618a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy"
1628a381b04SJed Brown /*@C
1638a381b04SJed Brown    TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
1648a381b04SJed Brown 
1658a381b04SJed Brown    Not Collective
1668a381b04SJed Brown 
1678a381b04SJed Brown    Level: advanced
1688a381b04SJed Brown 
1698a381b04SJed Brown .keywords: TSARKIMEX, register, destroy
1708a381b04SJed Brown .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic()
1718a381b04SJed Brown @*/
1728a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void)
1738a381b04SJed Brown {
1748a381b04SJed Brown   PetscErrorCode ierr;
1758a381b04SJed Brown   ARKTableauLink link;
1768a381b04SJed Brown 
1778a381b04SJed Brown   PetscFunctionBegin;
1788a381b04SJed Brown   while ((link = ARKTableauList)) {
1798a381b04SJed Brown     ARKTableau t = &link->tab;
1808a381b04SJed Brown     ARKTableauList = link->next;
1818a381b04SJed Brown     ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr);
182cd652676SJed Brown     ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr);
1838a381b04SJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
1848a381b04SJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
1858a381b04SJed Brown   }
1868a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_FALSE;
1878a381b04SJed Brown   PetscFunctionReturn(0);
1888a381b04SJed Brown }
1898a381b04SJed Brown 
1908a381b04SJed Brown #undef __FUNCT__
1918a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage"
1928a381b04SJed Brown /*@C
1938a381b04SJed Brown   TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
1948a381b04SJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
1958a381b04SJed Brown   when using static libraries.
1968a381b04SJed Brown 
1978a381b04SJed Brown   Input Parameter:
1988a381b04SJed Brown   path - The dynamic library path, or PETSC_NULL
1998a381b04SJed Brown 
2008a381b04SJed Brown   Level: developer
2018a381b04SJed Brown 
2028a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package
2038a381b04SJed Brown .seealso: PetscInitialize()
2048a381b04SJed Brown @*/
2058a381b04SJed Brown PetscErrorCode TSARKIMEXInitializePackage(const char path[])
2068a381b04SJed Brown {
2078a381b04SJed Brown   PetscErrorCode ierr;
2088a381b04SJed Brown 
2098a381b04SJed Brown   PetscFunctionBegin;
2108a381b04SJed Brown   if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0);
2118a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_TRUE;
2128a381b04SJed Brown   ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr);
2138a381b04SJed Brown   ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr);
2148a381b04SJed Brown   PetscFunctionReturn(0);
2158a381b04SJed Brown }
2168a381b04SJed Brown 
2178a381b04SJed Brown #undef __FUNCT__
2188a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage"
2198a381b04SJed Brown /*@C
2208a381b04SJed Brown   TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
2218a381b04SJed Brown   called from PetscFinalize().
2228a381b04SJed Brown 
2238a381b04SJed Brown   Level: developer
2248a381b04SJed Brown 
2258a381b04SJed Brown .keywords: Petsc, destroy, package
2268a381b04SJed Brown .seealso: PetscFinalize()
2278a381b04SJed Brown @*/
2288a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void)
2298a381b04SJed Brown {
2308a381b04SJed Brown   PetscErrorCode ierr;
2318a381b04SJed Brown 
2328a381b04SJed Brown   PetscFunctionBegin;
2338a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_FALSE;
2348a381b04SJed Brown   ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr);
2358a381b04SJed Brown   PetscFunctionReturn(0);
2368a381b04SJed Brown }
2378a381b04SJed Brown 
2388a381b04SJed Brown #undef __FUNCT__
2398a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister"
240cd652676SJed Brown /*@C
241cd652676SJed Brown    TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
242cd652676SJed Brown 
243cd652676SJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
244cd652676SJed Brown 
245cd652676SJed Brown    Input Parameters:
246cd652676SJed Brown +  name - identifier for method
247cd652676SJed Brown .  order - approximation order of method
248cd652676SJed Brown .  s - number of stages, this is the dimension of the matrices below
249cd652676SJed Brown .  At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
250cd652676SJed Brown .  bt - Butcher table for completing the stiff part of the step (dimension s; PETSC_NULL to use the last row of At)
251cd652676SJed Brown .  ct - Abscissa of each stiff stage (dimension s, PETSC_NULL to use row sums of At)
252cd652676SJed Brown .  A - Non-stiff stage coefficients (dimension s*s, row-major)
253cd652676SJed Brown .  b - Non-stiff step completion table (dimension s; PETSC_NULL to use last row of At)
254cd652676SJed Brown .  c - Non-stiff abscissa (dimension s; PETSC_NULL to use row sums of A)
255cd652676SJed Brown .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
256cd652676SJed Brown .  binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
257cd652676SJed Brown -  binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; PETSC_NULL to reuse binterpt)
258cd652676SJed Brown 
259cd652676SJed Brown    Notes:
260cd652676SJed Brown    Several ARK IMEX methods are provided, this function is only needed to create new methods.
261cd652676SJed Brown 
262cd652676SJed Brown    Level: advanced
263cd652676SJed Brown 
264cd652676SJed Brown .keywords: TS, register
265cd652676SJed Brown 
266cd652676SJed Brown .seealso: TSARKIMEX
267cd652676SJed Brown @*/
2688a381b04SJed Brown PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType name,PetscInt order,PetscInt s,
2698a381b04SJed Brown                                  const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
270cd652676SJed Brown                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
271cd652676SJed Brown                                  PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
2728a381b04SJed Brown {
2738a381b04SJed Brown   PetscErrorCode ierr;
2748a381b04SJed Brown   ARKTableauLink link;
2758a381b04SJed Brown   ARKTableau t;
2768a381b04SJed Brown   PetscInt i,j;
2778a381b04SJed Brown 
2788a381b04SJed Brown   PetscFunctionBegin;
2798a381b04SJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
280cd652676SJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
2818a381b04SJed Brown   t = &link->tab;
2828a381b04SJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
2838a381b04SJed Brown   t->order = order;
2848a381b04SJed Brown   t->s = s;
2858a381b04SJed Brown   ierr = PetscMalloc6(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s,PetscReal,&t->ct,s*s,PetscReal,&t->A,s,PetscReal,&t->b,s,PetscReal,&t->c);CHKERRQ(ierr);
2868a381b04SJed Brown   ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr);
2878a381b04SJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
2888a381b04SJed Brown   if (bt) {ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr);}
2898a381b04SJed Brown   else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
2908a381b04SJed Brown   if (b) {ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);}
2918a381b04SJed Brown   else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i];
2928a381b04SJed Brown   if (ct) {ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr);}
2938a381b04SJed Brown   else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j];
2948a381b04SJed Brown   if (c) {ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr);}
2958a381b04SJed Brown   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
296*4f385281SJed Brown   t->pinterp = pinterp;
297cd652676SJed Brown   ierr = PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);CHKERRQ(ierr);
298cd652676SJed Brown   ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
299cd652676SJed Brown   ierr = PetscMemcpy(t->binterp,binterp?binterp:binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
3008a381b04SJed Brown   link->next = ARKTableauList;
3018a381b04SJed Brown   ARKTableauList = link;
3028a381b04SJed Brown   PetscFunctionReturn(0);
3038a381b04SJed Brown }
3048a381b04SJed Brown 
3058a381b04SJed Brown #undef __FUNCT__
3068a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX"
3078a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts)
3088a381b04SJed Brown {
3098a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
3108a381b04SJed Brown   ARKTableau      tab  = ark->tableau;
3118a381b04SJed Brown   const PetscInt  s    = tab->s;
3128a381b04SJed Brown   const PetscReal *At  = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c;
313406d0ec2SJed Brown   PetscScalar     *w   = ark->work;
3148a381b04SJed Brown   Vec             *Y   = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,W = ark->Work,Z = ark->Z;
3158a381b04SJed Brown   SNES            snes;
3168a381b04SJed Brown   PetscInt        i,j,its,lits;
3178a381b04SJed Brown   PetscReal       h,t;
3188a381b04SJed Brown   PetscErrorCode  ierr;
3198a381b04SJed Brown 
3208a381b04SJed Brown   PetscFunctionBegin;
3218a381b04SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3228a381b04SJed Brown   h = ts->time_step = ts->next_time_step;
3238a381b04SJed Brown   t = ts->ptime;
3248a381b04SJed Brown 
3258a381b04SJed Brown   for (i=0; i<s; i++) {
3268a381b04SJed Brown     if (At[i*s+i] == 0) {           /* This stage is explicit */
3278a381b04SJed Brown       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
3288a381b04SJed Brown       for (j=0; j<i; j++) w[j] = -h*At[i*s+j];
3298a381b04SJed Brown       ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr);
3308a381b04SJed Brown       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
3318a381b04SJed Brown       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
3328a381b04SJed Brown     } else {
3338a381b04SJed Brown       ark->stage_time = t + h*ct[i];
3348a381b04SJed Brown       ark->shift = 1./(h*At[i*s+i]);
3358a381b04SJed Brown       /* Affine part */
3368a381b04SJed Brown       ierr = VecZeroEntries(W);CHKERRQ(ierr);
3378a381b04SJed Brown       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
3388a381b04SJed Brown       ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr);
3398a381b04SJed Brown       /* Ydot = shift*(Y-Z) */
3408a381b04SJed Brown       ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr);
341*4f385281SJed Brown       for (j=0; j<i; j++) w[j] = -h*At[i*s+j];
342*4f385281SJed Brown       ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr);
3438a381b04SJed Brown       /* Initial guess taken from last stage */
3448a381b04SJed Brown       ierr = VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);CHKERRQ(ierr);
3458a381b04SJed Brown       ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr);
3468a381b04SJed Brown       ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
3478a381b04SJed Brown       ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
3488a381b04SJed Brown       ts->nonlinear_its += its; ts->linear_its += lits;
3498a381b04SJed Brown     }
3508a381b04SJed Brown     ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
3518a381b04SJed Brown     ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],PETSC_TRUE);CHKERRQ(ierr);
3528a381b04SJed Brown     ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
3538a381b04SJed Brown   }
354a17b68daSJed Brown   for (j=0; j<s; j++) w[j] = -h*bt[j];
3558a381b04SJed Brown   ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr);
3568a381b04SJed Brown   for (j=0; j<s; j++) w[j] = h*b[j];
3578a381b04SJed Brown   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
3588a381b04SJed Brown 
3598a381b04SJed Brown   ts->ptime          += ts->time_step;
3608a381b04SJed Brown   ts->next_time_step  = ts->time_step;
3618a381b04SJed Brown   ts->steps++;
3628a381b04SJed Brown   PetscFunctionReturn(0);
3638a381b04SJed Brown }
3648a381b04SJed Brown 
365cd652676SJed Brown #undef __FUNCT__
366cd652676SJed Brown #define __FUNCT__ "TSInterpolate_ARKIMEX"
367cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
368cd652676SJed Brown {
369cd652676SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
370*4f385281SJed Brown   PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
371*4f385281SJed Brown   PetscReal tt,t = (itime - ts->ptime)/ts->time_step + 1; /* In the interval [0,1] */
372cd652676SJed Brown   PetscScalar *bt,*b;
373cd652676SJed Brown   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
374cd652676SJed Brown   PetscErrorCode ierr;
375cd652676SJed Brown 
376cd652676SJed Brown   PetscFunctionBegin;
377cd652676SJed Brown   if (!Bt || !B) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
378cd652676SJed Brown   ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr);
379cd652676SJed Brown   for (i=0; i<s; i++) bt[i] = b[i] = 0;
380*4f385281SJed Brown   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
381cd652676SJed Brown     for (i=0; i<s; i++) {
382*4f385281SJed Brown       bt[i] += ts->time_step * Bt[i*pinterp+j] * tt * -1.0;
383*4f385281SJed Brown       b[i]  += ts->time_step * B[i*pinterp+j] * tt;
384cd652676SJed Brown     }
385cd652676SJed Brown   }
386cd652676SJed Brown   if (ark->tableau->At[0*s+0] != 0.0) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_SUP,"First stage not explicit so starting stage not saved");
387cd652676SJed Brown   ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr);
388cd652676SJed Brown   ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr);
389cd652676SJed Brown   ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr);
390cd652676SJed Brown   ierr = PetscFree2(bt,b);CHKERRQ(ierr);
391cd652676SJed Brown   PetscFunctionReturn(0);
392cd652676SJed Brown }
393cd652676SJed Brown 
3948a381b04SJed Brown /*------------------------------------------------------------*/
3958a381b04SJed Brown #undef __FUNCT__
3968a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX"
3978a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts)
3988a381b04SJed Brown {
3998a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
4008a381b04SJed Brown   PetscInt        s;
4018a381b04SJed Brown   PetscErrorCode  ierr;
4028a381b04SJed Brown 
4038a381b04SJed Brown   PetscFunctionBegin;
4048a381b04SJed Brown   if (!ark->tableau) PetscFunctionReturn(0);
4058a381b04SJed Brown    s = ark->tableau->s;
4068a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr);
4078a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr);
4088a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr);
4098a381b04SJed Brown   ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr);
4108a381b04SJed Brown   ierr = VecDestroy(&ark->Work);CHKERRQ(ierr);
4118a381b04SJed Brown   ierr = VecDestroy(&ark->Z);CHKERRQ(ierr);
4128a381b04SJed Brown   ierr = PetscFree(ark->work);CHKERRQ(ierr);
4138a381b04SJed Brown   PetscFunctionReturn(0);
4148a381b04SJed Brown }
4158a381b04SJed Brown 
4168a381b04SJed Brown #undef __FUNCT__
4178a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX"
4188a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
4198a381b04SJed Brown {
4208a381b04SJed Brown   PetscErrorCode  ierr;
4218a381b04SJed Brown 
4228a381b04SJed Brown   PetscFunctionBegin;
4238a381b04SJed Brown   ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
4248a381b04SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
4258a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);CHKERRQ(ierr);
4268a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);CHKERRQ(ierr);
4278a381b04SJed Brown   PetscFunctionReturn(0);
4288a381b04SJed Brown }
4298a381b04SJed Brown 
4308a381b04SJed Brown /*
4318a381b04SJed Brown   This defines the nonlinear equation that is to be solved with SNES
4328a381b04SJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
4338a381b04SJed Brown */
4348a381b04SJed Brown #undef __FUNCT__
4358a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX"
4368a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
4378a381b04SJed Brown {
4388a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
4398a381b04SJed Brown   PetscErrorCode ierr;
4408a381b04SJed Brown 
4418a381b04SJed Brown   PetscFunctionBegin;
4428a381b04SJed Brown   ierr = VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
44331f6fcc0SJed Brown   ierr = TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,F,PETSC_TRUE);CHKERRQ(ierr);
4448a381b04SJed Brown   PetscFunctionReturn(0);
4458a381b04SJed Brown }
4468a381b04SJed Brown 
4478a381b04SJed Brown #undef __FUNCT__
4488a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX"
4498a381b04SJed Brown static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
4508a381b04SJed Brown {
4518a381b04SJed Brown   TS_ARKIMEX       *ark = (TS_ARKIMEX*)ts->data;
4528a381b04SJed Brown   PetscErrorCode ierr;
4538a381b04SJed Brown 
4548a381b04SJed Brown   PetscFunctionBegin;
4558a381b04SJed Brown   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
4568a381b04SJed Brown   ierr = TSComputeIJacobian(ts,ark->stage_time,X,ark->Ydot,ark->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
4578a381b04SJed Brown   PetscFunctionReturn(0);
4588a381b04SJed Brown }
4598a381b04SJed Brown 
4608a381b04SJed Brown #undef __FUNCT__
4618a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX"
4628a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
4638a381b04SJed Brown {
4648a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
4658a381b04SJed Brown   ARKTableau     tab  = ark->tableau;
4668a381b04SJed Brown   PetscInt       s = tab->s;
4678a381b04SJed Brown   PetscErrorCode ierr;
4688a381b04SJed Brown 
4698a381b04SJed Brown   PetscFunctionBegin;
4708a381b04SJed Brown   if (!ark->tableau) {
471e24355feSJed Brown     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
4728a381b04SJed Brown   }
4738a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr);
4748a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr);
4758a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr);
4768a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr);
4778a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr);
4788a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr);
4798a381b04SJed Brown   ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr);
4808a381b04SJed Brown   PetscFunctionReturn(0);
4818a381b04SJed Brown }
4828a381b04SJed Brown /*------------------------------------------------------------*/
4838a381b04SJed Brown 
4848a381b04SJed Brown #undef __FUNCT__
4858a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX"
4868a381b04SJed Brown static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts)
4878a381b04SJed Brown {
4888a381b04SJed Brown   PetscErrorCode ierr;
4898a381b04SJed Brown   char           arktype[256];
4908a381b04SJed Brown 
4918a381b04SJed Brown   PetscFunctionBegin;
4928a381b04SJed Brown   ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr);
4938a381b04SJed Brown   {
4948a381b04SJed Brown     ARKTableauLink link;
4958a381b04SJed Brown     PetscInt count,choice;
4968a381b04SJed Brown     PetscBool flg;
4978a381b04SJed Brown     const char **namelist;
4988a381b04SJed Brown     ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof arktype);CHKERRQ(ierr);
4998a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
5008a381b04SJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
5018a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
5028a381b04SJed Brown     ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr);
5038a381b04SJed Brown     ierr = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr);
5048a381b04SJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
5058a381b04SJed Brown   }
5068a381b04SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
5078a381b04SJed Brown   PetscFunctionReturn(0);
5088a381b04SJed Brown }
5098a381b04SJed Brown 
5108a381b04SJed Brown #undef __FUNCT__
5118a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray"
5128a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
5138a381b04SJed Brown {
5148a381b04SJed Brown   int i,left,count;
5158a381b04SJed Brown   char *p;
5168a381b04SJed Brown 
5178a381b04SJed Brown   PetscFunctionBegin;
5188a381b04SJed Brown   for (i=0,p=buf,left=(int)len; i<n; i++) {
5198a381b04SJed Brown     count = snprintf(p,left,fmt,x[i]);
5208a381b04SJed Brown     if (count < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"snprintf()");
5218a381b04SJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
5228a381b04SJed Brown     left -= count;
5238a381b04SJed Brown     p += count;
5248a381b04SJed Brown     *p++ = ' ';
5258a381b04SJed Brown   }
5268a381b04SJed Brown   p[i ? 0 : -1] = 0;
5278a381b04SJed Brown   PetscFunctionReturn(0);
5288a381b04SJed Brown }
5298a381b04SJed Brown 
5308a381b04SJed Brown #undef __FUNCT__
5318a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX"
5328a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
5338a381b04SJed Brown {
5348a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
5358a381b04SJed Brown   ARKTableau     tab = ark->tableau;
5368a381b04SJed Brown   PetscBool      iascii;
5378a381b04SJed Brown   PetscErrorCode ierr;
5388a381b04SJed Brown 
5398a381b04SJed Brown   PetscFunctionBegin;
5408a381b04SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5418a381b04SJed Brown   if (iascii) {
5428a381b04SJed Brown     const TSARKIMEXType arktype;
5438a381b04SJed Brown     char buf[512];
5448a381b04SJed Brown     ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr);
5458a381b04SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype);CHKERRQ(ierr);
5468a381b04SJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr);
54731f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf);CHKERRQ(ierr);
54831f6fcc0SJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
54931f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf);CHKERRQ(ierr);
5508a381b04SJed Brown   }
5518a381b04SJed Brown   PetscFunctionReturn(0);
5528a381b04SJed Brown }
5538a381b04SJed Brown 
5548a381b04SJed Brown #undef __FUNCT__
5558a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType"
5568a381b04SJed Brown /*@C
5578a381b04SJed Brown   TSARKIMEXSetType - Set the type of ARK IMEX scheme
5588a381b04SJed Brown 
5598a381b04SJed Brown   Logically collective
5608a381b04SJed Brown 
5618a381b04SJed Brown   Input Parameter:
5628a381b04SJed Brown +  ts - timestepping context
5638a381b04SJed Brown -  arktype - type of ARK-IMEX scheme
5648a381b04SJed Brown 
5658a381b04SJed Brown   Level: intermediate
5668a381b04SJed Brown 
5678a381b04SJed Brown .seealso: TSARKIMEXGetType()
5688a381b04SJed Brown @*/
5698a381b04SJed Brown PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType arktype)
5708a381b04SJed Brown {
5718a381b04SJed Brown   PetscErrorCode ierr;
5728a381b04SJed Brown 
5738a381b04SJed Brown   PetscFunctionBegin;
5748a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5758a381b04SJed Brown   ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,const TSARKIMEXType),(ts,arktype));CHKERRQ(ierr);
5768a381b04SJed Brown   PetscFunctionReturn(0);
5778a381b04SJed Brown }
5788a381b04SJed Brown 
5798a381b04SJed Brown #undef __FUNCT__
5808a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType"
5818a381b04SJed Brown /*@C
5828a381b04SJed Brown   TSARKIMEXGetType - Get the type of ARK IMEX scheme
5838a381b04SJed Brown 
5848a381b04SJed Brown   Logically collective
5858a381b04SJed Brown 
5868a381b04SJed Brown   Input Parameter:
5878a381b04SJed Brown .  ts - timestepping context
5888a381b04SJed Brown 
5898a381b04SJed Brown   Output Parameter:
5908a381b04SJed Brown .  arktype - type of ARK-IMEX scheme
5918a381b04SJed Brown 
5928a381b04SJed Brown   Level: intermediate
5938a381b04SJed Brown 
5948a381b04SJed Brown .seealso: TSARKIMEXGetType()
5958a381b04SJed Brown @*/
5968a381b04SJed Brown PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType *arktype)
5978a381b04SJed Brown {
5988a381b04SJed Brown   PetscErrorCode ierr;
5998a381b04SJed Brown 
6008a381b04SJed Brown   PetscFunctionBegin;
6018a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6028a381b04SJed Brown   ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,const TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr);
6038a381b04SJed Brown   PetscFunctionReturn(0);
6048a381b04SJed Brown }
6058a381b04SJed Brown 
6068a381b04SJed Brown EXTERN_C_BEGIN
6078a381b04SJed Brown #undef __FUNCT__
6088a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX"
6098a381b04SJed Brown PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,const TSARKIMEXType *arktype)
6108a381b04SJed Brown {
6118a381b04SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
6128a381b04SJed Brown   PetscErrorCode ierr;
6138a381b04SJed Brown 
6148a381b04SJed Brown   PetscFunctionBegin;
6158a381b04SJed Brown   if (!ark->tableau) {ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);}
6168a381b04SJed Brown   *arktype = ark->tableau->name;
6178a381b04SJed Brown   PetscFunctionReturn(0);
6188a381b04SJed Brown }
6198a381b04SJed Brown #undef __FUNCT__
6208a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX"
6218a381b04SJed Brown PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,const TSARKIMEXType arktype)
6228a381b04SJed Brown {
6238a381b04SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
6248a381b04SJed Brown   PetscErrorCode ierr;
6258a381b04SJed Brown   PetscBool match;
6268a381b04SJed Brown   ARKTableauLink link;
6278a381b04SJed Brown 
6288a381b04SJed Brown   PetscFunctionBegin;
6298a381b04SJed Brown   if (ark->tableau) {
6308a381b04SJed Brown     ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr);
6318a381b04SJed Brown     if (match) PetscFunctionReturn(0);
6328a381b04SJed Brown   }
6338a381b04SJed Brown   for (link = ARKTableauList; link; link=link->next) {
6348a381b04SJed Brown     ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr);
6358a381b04SJed Brown     if (match) {
6368a381b04SJed Brown       ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
6378a381b04SJed Brown       ark->tableau = &link->tab;
6388a381b04SJed Brown       PetscFunctionReturn(0);
6398a381b04SJed Brown     }
6408a381b04SJed Brown   }
6418a381b04SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
6428a381b04SJed Brown   PetscFunctionReturn(0);
6438a381b04SJed Brown }
6448a381b04SJed Brown EXTERN_C_END
6458a381b04SJed Brown 
6468a381b04SJed Brown /* ------------------------------------------------------------ */
6478a381b04SJed Brown /*MC
6488a381b04SJed Brown       TSARKIMEX - ODE solver using Additive Runge-Kutta IMEX schemes
6498a381b04SJed Brown 
650fca742c7SJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
651fca742c7SJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
652fca742c7SJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
653fca742c7SJed Brown 
654fca742c7SJed Brown   Notes:
655fca742c7SJed Brown   This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X).
656fca742c7SJed Brown 
6578a381b04SJed Brown   Level: beginner
6588a381b04SJed Brown 
659fca742c7SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSARKIMEXRegister()
6608a381b04SJed Brown 
6618a381b04SJed Brown M*/
6628a381b04SJed Brown EXTERN_C_BEGIN
6638a381b04SJed Brown #undef __FUNCT__
6648a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX"
6658a381b04SJed Brown PetscErrorCode  TSCreate_ARKIMEX(TS ts)
6668a381b04SJed Brown {
6678a381b04SJed Brown   TS_ARKIMEX       *th;
6688a381b04SJed Brown   PetscErrorCode ierr;
6698a381b04SJed Brown 
6708a381b04SJed Brown   PetscFunctionBegin;
6718a381b04SJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
6728a381b04SJed Brown   ierr = TSARKIMEXInitializePackage(PETSC_NULL);CHKERRQ(ierr);
6738a381b04SJed Brown #endif
6748a381b04SJed Brown 
6758a381b04SJed Brown   ts->ops->reset          = TSReset_ARKIMEX;
6768a381b04SJed Brown   ts->ops->destroy        = TSDestroy_ARKIMEX;
6778a381b04SJed Brown   ts->ops->view           = TSView_ARKIMEX;
6788a381b04SJed Brown   ts->ops->setup          = TSSetUp_ARKIMEX;
6798a381b04SJed Brown   ts->ops->step           = TSStep_ARKIMEX;
680cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_ARKIMEX;
6818a381b04SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
6828a381b04SJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
6838a381b04SJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;
6848a381b04SJed Brown 
6858a381b04SJed Brown   ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr);
6868a381b04SJed Brown   ts->data = (void*)th;
6878a381b04SJed Brown 
6888a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr);
6898a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr);
6908a381b04SJed Brown   PetscFunctionReturn(0);
6918a381b04SJed Brown }
6928a381b04SJed Brown EXTERN_C_END
693