xref: /petsc/src/ts/impls/arkimex/arkimex.c (revision e24355fe3772bda2877cb2f6faabfa49e8b33b37)
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 
148a381b04SJed Brown static const TSARKIMEXType TSARKIMEXDefault = TSARKIMEX2D;
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;
218a381b04SJed Brown   PetscInt order;
228a381b04SJed Brown   PetscInt s;
238a381b04SJed Brown   PetscReal *At,*bt,*ct;
248a381b04SJed Brown   PetscReal *A,*b,*c;           /* Non-stiff tableau */
258a381b04SJed Brown };
268a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink;
278a381b04SJed Brown struct _ARKTableauLink {
288a381b04SJed Brown   struct _ARKTableau tab;
298a381b04SJed Brown   ARKTableauLink next;
308a381b04SJed Brown };
318a381b04SJed Brown static ARKTableauLink ARKTableauList;
328a381b04SJed Brown 
338a381b04SJed Brown typedef struct {
348a381b04SJed Brown   ARKTableau  tableau;
358a381b04SJed Brown   Vec         *Y;               /* States computed during the step */
368a381b04SJed Brown   Vec         *YdotI;           /* Time derivatives for the stiff part */
378a381b04SJed Brown   Vec         *YdotRHS;         /* Function evaluations for the non-stiff part */
388a381b04SJed Brown   Vec         Ydot;             /* Work vector holding Ydot during residual evaluation */
398a381b04SJed Brown   Vec         Work;             /* Generic work vector */
408a381b04SJed Brown   Vec         Z;                /* Ydot = shift(Y-Z) */
418a381b04SJed Brown   PetscScalar *work;            /* Scalar work */
428a381b04SJed Brown   PetscReal   shift;
438a381b04SJed Brown   PetscReal   stage_time;
448a381b04SJed Brown } TS_ARKIMEX;
458a381b04SJed Brown 
468a381b04SJed Brown #undef __FUNCT__
478a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterAll"
488a381b04SJed Brown /*@C
498a381b04SJed Brown   TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX
508a381b04SJed Brown 
518a381b04SJed Brown   Not Collective
528a381b04SJed Brown 
538a381b04SJed Brown   Level: advanced
548a381b04SJed Brown 
558a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all
568a381b04SJed Brown 
578a381b04SJed Brown .seealso:  TSARKIMEXRegisterDestroy()
588a381b04SJed Brown @*/
598a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void)
608a381b04SJed Brown {
618a381b04SJed Brown   PetscErrorCode ierr;
628a381b04SJed Brown 
638a381b04SJed Brown   PetscFunctionBegin;
648a381b04SJed Brown   if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0);
658a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_TRUE;
668a381b04SJed Brown 
678a381b04SJed Brown   {
688a381b04SJed Brown     const PetscReal
698a381b04SJed Brown       A[3][3] = {{0,0,0},
708a381b04SJed Brown                  {0.41421356237309504880,0,0},
718a381b04SJed Brown                  {0.75,0.25,0}},
728a381b04SJed Brown       At[3][3] = {{0,0,0},
738a381b04SJed Brown                   {0.12132034355964257320,0.29289321881345247560,0},
748a381b04SJed Brown                   {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}};
758a381b04SJed Brown       ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
768a381b04SJed Brown   }
77a3a57f36SJed Brown   {
78a3a57f36SJed Brown     const PetscReal s2 = sqrt(2),
79a3a57f36SJed Brown       A[3][3] = {{0,0,0},
80a3a57f36SJed Brown                  {2-s2,0,0},
81a3a57f36SJed Brown                  {(3-2*s2)/6,(3+2*s2)/6,0}},
82a3a57f36SJed Brown       At[3][3] = {{0,0,0},
83a3a57f36SJed Brown                   {1-1/s2,1-1/s2,0},
84a3a57f36SJed Brown                   {1/(2*s2),1/(2*s2),1-1/s2}};
85a3a57f36SJed Brown       ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
86a3a57f36SJed Brown   }
87a3a57f36SJed Brown   {
88a3a57f36SJed Brown     const PetscReal
89a3a57f36SJed Brown       A[4][4] = {{0,0,0,0},
90a3a57f36SJed Brown                  {1767732205903./2027836641118,0,0,0},
91a3a57f36SJed Brown                  {5535828885825./10492691773637,788022342437./10882634858940,0,0},
92a3a57f36SJed Brown                  {6485989280629./16251701735622,-4246266847089./9704473918619,10755448449292./10357097424841,0}},
93a3a57f36SJed Brown       At[4][4] = {{0,0,0,0},
94a3a57f36SJed Brown                   {1767732205903./4055673282236,1767732205903./4055673282236,0,0},
95a3a57f36SJed Brown                   {2746238789719./10658868560708,-640167445237./6845629431997,1767732205903./4055673282236,0},
96a3a57f36SJed Brown                   {1471266399579./7840856788654,-4482444167858./7529755066697,11266239266428./11593286722821,1767732205903./4055673282236}};
97a3a57f36SJed Brown       ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
98a3a57f36SJed Brown   }
99a3a57f36SJed Brown   {
100a3a57f36SJed Brown     const PetscReal
101a3a57f36SJed Brown       A[6][6] = {{0,0,0,0,0,0},
102a3a57f36SJed Brown                  {1./2,0,0,0,0,0},
103a3a57f36SJed Brown                  {13861./62500,6889./62500,0,0,0,0},
104a3a57f36SJed Brown                  {-116923316275./2393684061468,-2731218467317./15368042101831,9408046702089./11113171139209,0,0,0},
105a3a57f36SJed Brown                  {-451086348788./2902428689909,-2682348792572./7519795681897,12662868775082./11960479115383,3355817975965./11060851509271,0,0},
106a3a57f36SJed Brown                  {647845179188./3216320057751,73281519250./8382639484533,552539513391./3454668386233,3354512671639./8306763924573,4040./17871,0}},
107a3a57f36SJed Brown       At[6][6] = {{0,0,0,0,0,0},
108a3a57f36SJed Brown                   {1./4,1./4,0,0,0,0},
109a3a57f36SJed Brown                   {8611./62500,-1743./31250,1./4,0,0,0},
110a3a57f36SJed Brown                   {5012029./34652500,-654441./2922500,174375./388108,1./4,0,0},
111a3a57f36SJed Brown                   {15267082809./155376265600,-71443401./120774400,730878875./902184768,2285395./8070912,1./4,0},
112a3a57f36SJed Brown                   {82889./524892,0,15625./83664,69875./102672,-2260./8211,1./4}};
113a3a57f36SJed Brown       ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
114a3a57f36SJed Brown   }
115a3a57f36SJed Brown   {
116a3a57f36SJed Brown     const PetscReal
117a3a57f36SJed Brown       A[8][8] = {{0,0,0,0,0,0,0,0},
118a3a57f36SJed Brown                  {41./100,0,0,0,0,0,0,0},
119a3a57f36SJed Brown                  {367902744464./2072280473677,677623207551./8224143866563,0,0,0,0,0,0},
120a3a57f36SJed Brown                  {1268023523408./10340822734521,0,1029933939417./13636558850479,0,0,0,0,0},
121a3a57f36SJed Brown                  {14463281900351./6315353703477,0,66114435211212./5879490589093,-54053170152839./4284798021562,0,0,0,0},
122a3a57f36SJed Brown                  {14090043504691./34967701212078,0,15191511035443./11219624916014,-18461159152457./12425892160975,-281667163811./9011619295870,0,0,0},
123a3a57f36SJed Brown                  {19230459214898./13134317526959,0,21275331358303./2942455364971,-38145345988419./4862620318723,-1./8,-1./8,0,0},
124a3a57f36SJed Brown                  {-19977161125411./11928030595625,0,-40795976796054./6384907823539,177454434618887./12078138498510,782672205425./8267701900261,-69563011059811./9646580694205,7356628210526./4942186776405,0}},
125a3a57f36SJed Brown       At[8][8] = {{0,0,0,0,0,0,0,0},
126a3a57f36SJed Brown                   {41./200,41./200,0,0,0,0,0,0},
127a3a57f36SJed Brown                   {41./400,-567603406766./11931857230679,41./200,0,0,0,0,0},
128a3a57f36SJed Brown                   {683785636431./9252920307686,0,-110385047103./1367015193373,41./200,0,0,0,0},
129a3a57f36SJed Brown                   {3016520224154./10081342136671,0,30586259806659./12414158314087,-22760509404356./11113319521817,41./200,0,0,0},
130a3a57f36SJed Brown                   {218866479029./1489978393911,0,638256894668./5436446318841,-1179710474555./5321154724896,-60928119172./8023461067671,41./200,0,0},
131a3a57f36SJed Brown                   {1020004230633./5715676835656,0,25762820946817./25263940353407,-2161375909145./9755907335909,-211217309593./5846859502534,-4269925059573./7827059040749,41./200,0},
132a3a57f36SJed Brown                   {-872700587467./9133579230613,0,0,22348218063261./9555858737531,-1143369518992./8141816002931,-39379526789629./19018526304540,32727382324388./42900044865799,41./200}};
133a3a57f36SJed Brown       ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
134a3a57f36SJed Brown   }
135a3a57f36SJed Brown 
1368a381b04SJed Brown   PetscFunctionReturn(0);
1378a381b04SJed Brown }
1388a381b04SJed Brown 
1398a381b04SJed Brown #undef __FUNCT__
1408a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy"
1418a381b04SJed Brown /*@C
1428a381b04SJed Brown    TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
1438a381b04SJed Brown 
1448a381b04SJed Brown    Not Collective
1458a381b04SJed Brown 
1468a381b04SJed Brown    Level: advanced
1478a381b04SJed Brown 
1488a381b04SJed Brown .keywords: TSARKIMEX, register, destroy
1498a381b04SJed Brown .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic()
1508a381b04SJed Brown @*/
1518a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void)
1528a381b04SJed Brown {
1538a381b04SJed Brown   PetscErrorCode ierr;
1548a381b04SJed Brown   ARKTableauLink link;
1558a381b04SJed Brown 
1568a381b04SJed Brown   PetscFunctionBegin;
1578a381b04SJed Brown   while ((link = ARKTableauList)) {
1588a381b04SJed Brown     ARKTableau t = &link->tab;
1598a381b04SJed Brown     ARKTableauList = link->next;
1608a381b04SJed Brown     ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr);
1618a381b04SJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
1628a381b04SJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
1638a381b04SJed Brown   }
1648a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_FALSE;
1658a381b04SJed Brown   PetscFunctionReturn(0);
1668a381b04SJed Brown }
1678a381b04SJed Brown 
1688a381b04SJed Brown #undef __FUNCT__
1698a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage"
1708a381b04SJed Brown /*@C
1718a381b04SJed Brown   TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
1728a381b04SJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
1738a381b04SJed Brown   when using static libraries.
1748a381b04SJed Brown 
1758a381b04SJed Brown   Input Parameter:
1768a381b04SJed Brown   path - The dynamic library path, or PETSC_NULL
1778a381b04SJed Brown 
1788a381b04SJed Brown   Level: developer
1798a381b04SJed Brown 
1808a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package
1818a381b04SJed Brown .seealso: PetscInitialize()
1828a381b04SJed Brown @*/
1838a381b04SJed Brown PetscErrorCode TSARKIMEXInitializePackage(const char path[])
1848a381b04SJed Brown {
1858a381b04SJed Brown   PetscErrorCode ierr;
1868a381b04SJed Brown 
1878a381b04SJed Brown   PetscFunctionBegin;
1888a381b04SJed Brown   if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0);
1898a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_TRUE;
1908a381b04SJed Brown   ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr);
1918a381b04SJed Brown   ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr);
1928a381b04SJed Brown   PetscFunctionReturn(0);
1938a381b04SJed Brown }
1948a381b04SJed Brown 
1958a381b04SJed Brown #undef __FUNCT__
1968a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage"
1978a381b04SJed Brown /*@C
1988a381b04SJed Brown   TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
1998a381b04SJed Brown   called from PetscFinalize().
2008a381b04SJed Brown 
2018a381b04SJed Brown   Level: developer
2028a381b04SJed Brown 
2038a381b04SJed Brown .keywords: Petsc, destroy, package
2048a381b04SJed Brown .seealso: PetscFinalize()
2058a381b04SJed Brown @*/
2068a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void)
2078a381b04SJed Brown {
2088a381b04SJed Brown   PetscErrorCode ierr;
2098a381b04SJed Brown 
2108a381b04SJed Brown   PetscFunctionBegin;
2118a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_FALSE;
2128a381b04SJed Brown   ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr);
2138a381b04SJed Brown   PetscFunctionReturn(0);
2148a381b04SJed Brown }
2158a381b04SJed Brown 
2168a381b04SJed Brown #undef __FUNCT__
2178a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister"
2188a381b04SJed Brown PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType name,PetscInt order,PetscInt s,
2198a381b04SJed Brown                                  const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
2208a381b04SJed Brown                                  const PetscReal A[],const PetscReal b[],const PetscReal c[])
2218a381b04SJed Brown {
2228a381b04SJed Brown   PetscErrorCode ierr;
2238a381b04SJed Brown   ARKTableauLink link;
2248a381b04SJed Brown   ARKTableau t;
2258a381b04SJed Brown   PetscInt i,j;
2268a381b04SJed Brown 
2278a381b04SJed Brown   PetscFunctionBegin;
2288a381b04SJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
2298a381b04SJed Brown   t = &link->tab;
2308a381b04SJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
2318a381b04SJed Brown   t->order = order;
2328a381b04SJed Brown   t->s = s;
2338a381b04SJed 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);
2348a381b04SJed Brown   ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr);
2358a381b04SJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
2368a381b04SJed Brown   if (bt) {ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr);}
2378a381b04SJed Brown   else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
2388a381b04SJed Brown   if (b) {ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);}
2398a381b04SJed Brown   else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i];
2408a381b04SJed Brown   if (ct) {ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr);}
2418a381b04SJed 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];
2428a381b04SJed Brown   if (c) {ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr);}
2438a381b04SJed 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];
2448a381b04SJed Brown   link->next = ARKTableauList;
2458a381b04SJed Brown   ARKTableauList = link;
2468a381b04SJed Brown   PetscFunctionReturn(0);
2478a381b04SJed Brown }
2488a381b04SJed Brown 
2498a381b04SJed Brown #undef __FUNCT__
2508a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX"
2518a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts)
2528a381b04SJed Brown {
2538a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
2548a381b04SJed Brown   ARKTableau      tab  = ark->tableau;
2558a381b04SJed Brown   const PetscInt  s    = tab->s;
2568a381b04SJed Brown   const PetscReal *At  = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c;
257406d0ec2SJed Brown   PetscScalar     *w   = ark->work;
2588a381b04SJed Brown   Vec             *Y   = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,W = ark->Work,Z = ark->Z;
2598a381b04SJed Brown   SNES            snes;
2608a381b04SJed Brown   PetscInt        i,j,its,lits;
2618a381b04SJed Brown   PetscReal       h,t;
2628a381b04SJed Brown   PetscErrorCode  ierr;
2638a381b04SJed Brown 
2648a381b04SJed Brown   PetscFunctionBegin;
2658a381b04SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2668a381b04SJed Brown   h = ts->time_step = ts->next_time_step;
2678a381b04SJed Brown   t = ts->ptime;
2688a381b04SJed Brown 
2698a381b04SJed Brown   for (i=0; i<s; i++) {
2708a381b04SJed Brown     if (At[i*s+i] == 0) {           /* This stage is explicit */
2718a381b04SJed Brown       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
2728a381b04SJed Brown       for (j=0; j<i; j++) w[j] = -h*At[i*s+j];
2738a381b04SJed Brown       ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr);
2748a381b04SJed Brown       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
2758a381b04SJed Brown       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
2768a381b04SJed Brown     } else {
2778a381b04SJed Brown       ark->stage_time = t + h*ct[i];
2788a381b04SJed Brown       ark->shift = 1./(h*At[i*s+i]);
2798a381b04SJed Brown       /* Affine part */
2808a381b04SJed Brown       ierr = VecZeroEntries(W);CHKERRQ(ierr);
2818a381b04SJed Brown       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
2828a381b04SJed Brown       ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr);
2838a381b04SJed Brown       /* Ydot = shift*(Y-Z) */
2848a381b04SJed Brown       ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr);
2858a381b04SJed Brown       for (j=0; j<i; j++) w[j] = h*At[i*s+j];
2868a381b04SJed Brown       ierr = VecMAXPY(Z,i,w,YdotRHS);CHKERRQ(ierr);
2878a381b04SJed Brown       /* Initial guess taken from last stage */
2888a381b04SJed Brown       ierr = VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);CHKERRQ(ierr);
2898a381b04SJed Brown       ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr);
2908a381b04SJed Brown       ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
2918a381b04SJed Brown       ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
2928a381b04SJed Brown       ts->nonlinear_its += its; ts->linear_its += lits;
2938a381b04SJed Brown     }
2948a381b04SJed Brown     ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
2958a381b04SJed Brown     ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],PETSC_TRUE);CHKERRQ(ierr);
2968a381b04SJed Brown     ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
2978a381b04SJed Brown   }
298a17b68daSJed Brown   for (j=0; j<s; j++) w[j] = -h*bt[j];
2998a381b04SJed Brown   ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr);
3008a381b04SJed Brown   for (j=0; j<s; j++) w[j] = h*b[j];
3018a381b04SJed Brown   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
3028a381b04SJed Brown 
3038a381b04SJed Brown   ts->ptime          += ts->time_step;
3048a381b04SJed Brown   ts->next_time_step  = ts->time_step;
3058a381b04SJed Brown   ts->steps++;
3068a381b04SJed Brown   PetscFunctionReturn(0);
3078a381b04SJed Brown }
3088a381b04SJed Brown 
3098a381b04SJed Brown /*------------------------------------------------------------*/
3108a381b04SJed Brown #undef __FUNCT__
3118a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX"
3128a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts)
3138a381b04SJed Brown {
3148a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
3158a381b04SJed Brown   PetscInt        s;
3168a381b04SJed Brown   PetscErrorCode  ierr;
3178a381b04SJed Brown 
3188a381b04SJed Brown   PetscFunctionBegin;
3198a381b04SJed Brown   if (!ark->tableau) PetscFunctionReturn(0);
3208a381b04SJed Brown    s = ark->tableau->s;
3218a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr);
3228a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr);
3238a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr);
3248a381b04SJed Brown   ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr);
3258a381b04SJed Brown   ierr = VecDestroy(&ark->Work);CHKERRQ(ierr);
3268a381b04SJed Brown   ierr = VecDestroy(&ark->Z);CHKERRQ(ierr);
3278a381b04SJed Brown   ierr = PetscFree(ark->work);CHKERRQ(ierr);
3288a381b04SJed Brown   PetscFunctionReturn(0);
3298a381b04SJed Brown }
3308a381b04SJed Brown 
3318a381b04SJed Brown #undef __FUNCT__
3328a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX"
3338a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
3348a381b04SJed Brown {
3358a381b04SJed Brown   PetscErrorCode  ierr;
3368a381b04SJed Brown 
3378a381b04SJed Brown   PetscFunctionBegin;
3388a381b04SJed Brown   ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
3398a381b04SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
3408a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);CHKERRQ(ierr);
3418a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);CHKERRQ(ierr);
3428a381b04SJed Brown   PetscFunctionReturn(0);
3438a381b04SJed Brown }
3448a381b04SJed Brown 
3458a381b04SJed Brown /*
3468a381b04SJed Brown   This defines the nonlinear equation that is to be solved with SNES
3478a381b04SJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
3488a381b04SJed Brown */
3498a381b04SJed Brown #undef __FUNCT__
3508a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX"
3518a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
3528a381b04SJed Brown {
3538a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
3548a381b04SJed Brown   PetscErrorCode ierr;
3558a381b04SJed Brown 
3568a381b04SJed Brown   PetscFunctionBegin;
3578a381b04SJed Brown   ierr = VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
35831f6fcc0SJed Brown   ierr = TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,F,PETSC_TRUE);CHKERRQ(ierr);
3598a381b04SJed Brown   PetscFunctionReturn(0);
3608a381b04SJed Brown }
3618a381b04SJed Brown 
3628a381b04SJed Brown #undef __FUNCT__
3638a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX"
3648a381b04SJed Brown static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
3658a381b04SJed Brown {
3668a381b04SJed Brown   TS_ARKIMEX       *ark = (TS_ARKIMEX*)ts->data;
3678a381b04SJed Brown   PetscErrorCode ierr;
3688a381b04SJed Brown 
3698a381b04SJed Brown   PetscFunctionBegin;
3708a381b04SJed Brown   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
3718a381b04SJed Brown   ierr = TSComputeIJacobian(ts,ark->stage_time,X,ark->Ydot,ark->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
3728a381b04SJed Brown   PetscFunctionReturn(0);
3738a381b04SJed Brown }
3748a381b04SJed Brown 
3758a381b04SJed Brown #undef __FUNCT__
3768a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX"
3778a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
3788a381b04SJed Brown {
3798a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
3808a381b04SJed Brown   ARKTableau     tab  = ark->tableau;
3818a381b04SJed Brown   PetscInt       s = tab->s;
3828a381b04SJed Brown   PetscErrorCode ierr;
3838a381b04SJed Brown 
3848a381b04SJed Brown   PetscFunctionBegin;
3858a381b04SJed Brown   if (!ark->tableau) {
386*e24355feSJed Brown     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
3878a381b04SJed Brown   }
3888a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr);
3898a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr);
3908a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr);
3918a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr);
3928a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr);
3938a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr);
3948a381b04SJed Brown   ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr);
3958a381b04SJed Brown   PetscFunctionReturn(0);
3968a381b04SJed Brown }
3978a381b04SJed Brown /*------------------------------------------------------------*/
3988a381b04SJed Brown 
3998a381b04SJed Brown #undef __FUNCT__
4008a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX"
4018a381b04SJed Brown static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts)
4028a381b04SJed Brown {
4038a381b04SJed Brown   PetscErrorCode ierr;
4048a381b04SJed Brown   char           arktype[256];
4058a381b04SJed Brown 
4068a381b04SJed Brown   PetscFunctionBegin;
4078a381b04SJed Brown   ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr);
4088a381b04SJed Brown   {
4098a381b04SJed Brown     ARKTableauLink link;
4108a381b04SJed Brown     PetscInt count,choice;
4118a381b04SJed Brown     PetscBool flg;
4128a381b04SJed Brown     const char **namelist;
4138a381b04SJed Brown     ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof arktype);CHKERRQ(ierr);
4148a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
4158a381b04SJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
4168a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
4178a381b04SJed Brown     ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr);
4188a381b04SJed Brown     ierr = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr);
4198a381b04SJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
4208a381b04SJed Brown   }
4218a381b04SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
4228a381b04SJed Brown   PetscFunctionReturn(0);
4238a381b04SJed Brown }
4248a381b04SJed Brown 
4258a381b04SJed Brown #undef __FUNCT__
4268a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray"
4278a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
4288a381b04SJed Brown {
4298a381b04SJed Brown   int i,left,count;
4308a381b04SJed Brown   char *p;
4318a381b04SJed Brown 
4328a381b04SJed Brown   PetscFunctionBegin;
4338a381b04SJed Brown   for (i=0,p=buf,left=(int)len; i<n; i++) {
4348a381b04SJed Brown     count = snprintf(p,left,fmt,x[i]);
4358a381b04SJed Brown     if (count < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"snprintf()");
4368a381b04SJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
4378a381b04SJed Brown     left -= count;
4388a381b04SJed Brown     p += count;
4398a381b04SJed Brown     *p++ = ' ';
4408a381b04SJed Brown   }
4418a381b04SJed Brown   p[i ? 0 : -1] = 0;
4428a381b04SJed Brown   PetscFunctionReturn(0);
4438a381b04SJed Brown }
4448a381b04SJed Brown 
4458a381b04SJed Brown #undef __FUNCT__
4468a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX"
4478a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
4488a381b04SJed Brown {
4498a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
4508a381b04SJed Brown   ARKTableau     tab = ark->tableau;
4518a381b04SJed Brown   PetscBool      iascii;
4528a381b04SJed Brown   PetscErrorCode ierr;
4538a381b04SJed Brown 
4548a381b04SJed Brown   PetscFunctionBegin;
4558a381b04SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4568a381b04SJed Brown   if (iascii) {
4578a381b04SJed Brown     const TSARKIMEXType arktype;
4588a381b04SJed Brown     char buf[512];
4598a381b04SJed Brown     ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr);
4608a381b04SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype);CHKERRQ(ierr);
4618a381b04SJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr);
46231f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf);CHKERRQ(ierr);
46331f6fcc0SJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
46431f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf);CHKERRQ(ierr);
4658a381b04SJed Brown   }
4668a381b04SJed Brown   PetscFunctionReturn(0);
4678a381b04SJed Brown }
4688a381b04SJed Brown 
4698a381b04SJed Brown #undef __FUNCT__
4708a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType"
4718a381b04SJed Brown /*@C
4728a381b04SJed Brown   TSARKIMEXSetType - Set the type of ARK IMEX scheme
4738a381b04SJed Brown 
4748a381b04SJed Brown   Logically collective
4758a381b04SJed Brown 
4768a381b04SJed Brown   Input Parameter:
4778a381b04SJed Brown +  ts - timestepping context
4788a381b04SJed Brown -  arktype - type of ARK-IMEX scheme
4798a381b04SJed Brown 
4808a381b04SJed Brown   Level: intermediate
4818a381b04SJed Brown 
4828a381b04SJed Brown .seealso: TSARKIMEXGetType()
4838a381b04SJed Brown @*/
4848a381b04SJed Brown PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType arktype)
4858a381b04SJed Brown {
4868a381b04SJed Brown   PetscErrorCode ierr;
4878a381b04SJed Brown 
4888a381b04SJed Brown   PetscFunctionBegin;
4898a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4908a381b04SJed Brown   ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,const TSARKIMEXType),(ts,arktype));CHKERRQ(ierr);
4918a381b04SJed Brown   PetscFunctionReturn(0);
4928a381b04SJed Brown }
4938a381b04SJed Brown 
4948a381b04SJed Brown #undef __FUNCT__
4958a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType"
4968a381b04SJed Brown /*@C
4978a381b04SJed Brown   TSARKIMEXGetType - Get the type of ARK IMEX scheme
4988a381b04SJed Brown 
4998a381b04SJed Brown   Logically collective
5008a381b04SJed Brown 
5018a381b04SJed Brown   Input Parameter:
5028a381b04SJed Brown .  ts - timestepping context
5038a381b04SJed Brown 
5048a381b04SJed Brown   Output Parameter:
5058a381b04SJed Brown .  arktype - type of ARK-IMEX scheme
5068a381b04SJed Brown 
5078a381b04SJed Brown   Level: intermediate
5088a381b04SJed Brown 
5098a381b04SJed Brown .seealso: TSARKIMEXGetType()
5108a381b04SJed Brown @*/
5118a381b04SJed Brown PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType *arktype)
5128a381b04SJed Brown {
5138a381b04SJed Brown   PetscErrorCode ierr;
5148a381b04SJed Brown 
5158a381b04SJed Brown   PetscFunctionBegin;
5168a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5178a381b04SJed Brown   ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,const TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr);
5188a381b04SJed Brown   PetscFunctionReturn(0);
5198a381b04SJed Brown }
5208a381b04SJed Brown 
5218a381b04SJed Brown EXTERN_C_BEGIN
5228a381b04SJed Brown #undef __FUNCT__
5238a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX"
5248a381b04SJed Brown PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,const TSARKIMEXType *arktype)
5258a381b04SJed Brown {
5268a381b04SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
5278a381b04SJed Brown   PetscErrorCode ierr;
5288a381b04SJed Brown 
5298a381b04SJed Brown   PetscFunctionBegin;
5308a381b04SJed Brown   if (!ark->tableau) {ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);}
5318a381b04SJed Brown   *arktype = ark->tableau->name;
5328a381b04SJed Brown   PetscFunctionReturn(0);
5338a381b04SJed Brown }
5348a381b04SJed Brown #undef __FUNCT__
5358a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX"
5368a381b04SJed Brown PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,const TSARKIMEXType arktype)
5378a381b04SJed Brown {
5388a381b04SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
5398a381b04SJed Brown   PetscErrorCode ierr;
5408a381b04SJed Brown   PetscBool match;
5418a381b04SJed Brown   ARKTableauLink link;
5428a381b04SJed Brown 
5438a381b04SJed Brown   PetscFunctionBegin;
5448a381b04SJed Brown   if (ark->tableau) {
5458a381b04SJed Brown     ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr);
5468a381b04SJed Brown     if (match) PetscFunctionReturn(0);
5478a381b04SJed Brown   }
5488a381b04SJed Brown   for (link = ARKTableauList; link; link=link->next) {
5498a381b04SJed Brown     ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr);
5508a381b04SJed Brown     if (match) {
5518a381b04SJed Brown       ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
5528a381b04SJed Brown       ark->tableau = &link->tab;
5538a381b04SJed Brown       PetscFunctionReturn(0);
5548a381b04SJed Brown     }
5558a381b04SJed Brown   }
5568a381b04SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
5578a381b04SJed Brown   PetscFunctionReturn(0);
5588a381b04SJed Brown }
5598a381b04SJed Brown EXTERN_C_END
5608a381b04SJed Brown 
5618a381b04SJed Brown /* ------------------------------------------------------------ */
5628a381b04SJed Brown /*MC
5638a381b04SJed Brown       TSARKIMEX - ODE solver using Additive Runge-Kutta IMEX schemes
5648a381b04SJed Brown 
5658a381b04SJed Brown   Level: beginner
5668a381b04SJed Brown 
5678a381b04SJed Brown .seealso:  TSCreate(), TS, TSSetType()
5688a381b04SJed Brown 
5698a381b04SJed Brown M*/
5708a381b04SJed Brown EXTERN_C_BEGIN
5718a381b04SJed Brown #undef __FUNCT__
5728a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX"
5738a381b04SJed Brown PetscErrorCode  TSCreate_ARKIMEX(TS ts)
5748a381b04SJed Brown {
5758a381b04SJed Brown   TS_ARKIMEX       *th;
5768a381b04SJed Brown   PetscErrorCode ierr;
5778a381b04SJed Brown 
5788a381b04SJed Brown   PetscFunctionBegin;
5798a381b04SJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
5808a381b04SJed Brown   ierr = TSARKIMEXInitializePackage(PETSC_NULL);CHKERRQ(ierr);
5818a381b04SJed Brown #endif
5828a381b04SJed Brown 
5838a381b04SJed Brown   ts->ops->reset          = TSReset_ARKIMEX;
5848a381b04SJed Brown   ts->ops->destroy        = TSDestroy_ARKIMEX;
5858a381b04SJed Brown   ts->ops->view           = TSView_ARKIMEX;
5868a381b04SJed Brown   ts->ops->setup          = TSSetUp_ARKIMEX;
5878a381b04SJed Brown   ts->ops->step           = TSStep_ARKIMEX;
5888a381b04SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
5898a381b04SJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
5908a381b04SJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;
5918a381b04SJed Brown 
5928a381b04SJed Brown   ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr);
5938a381b04SJed Brown   ts->data = (void*)th;
5948a381b04SJed Brown 
5958a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr);
5968a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr);
5978a381b04SJed Brown   PetscFunctionReturn(0);
5988a381b04SJed Brown }
5998a381b04SJed Brown EXTERN_C_END
600