xref: /petsc/src/ts/impls/arkimex/arkimex.c (revision 26885f590fc0384dd77b642cf8e35a4e520c99b3)
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 
144f385281SJed 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;
214f385281SJed Brown   PetscInt order;               /* Classical approximation order of the method */
224f385281SJed Brown   PetscInt s;                   /* Number of stages */
234f385281SJed Brown   PetscInt pinterp;             /* Interpolation order */
244f385281SJed 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;
464cc180ffSJed Brown   PetscBool   imex;
478a381b04SJed Brown } TS_ARKIMEX;
488a381b04SJed Brown 
4964f491ddSJed Brown /*MC
5064f491ddSJed Brown      TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part.
5164f491ddSJed Brown 
5264f491ddSJed Brown      This method has one explicit stage and two implicit stages.
5364f491ddSJed Brown 
5464f491ddSJed Brown .seealso: TSARKIMEX
5564f491ddSJed Brown M*/
5664f491ddSJed Brown /*MC
5764f491ddSJed Brown      TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part.
5864f491ddSJed Brown 
5964f491ddSJed Brown      This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu.
6064f491ddSJed Brown 
6164f491ddSJed Brown .seealso: TSARKIMEX
6264f491ddSJed Brown M*/
6364f491ddSJed Brown /*MC
6464f491ddSJed Brown      TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part.
6564f491ddSJed Brown 
6664f491ddSJed Brown      This method has one explicit stage and three implicit stages.
6764f491ddSJed Brown 
6864f491ddSJed Brown      References:
6964f491ddSJed Brown      Kennedy and Carpenter 2003.
7064f491ddSJed Brown 
7164f491ddSJed Brown .seealso: TSARKIMEX
7264f491ddSJed Brown M*/
7364f491ddSJed Brown /*MC
7464f491ddSJed Brown      TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part.
7564f491ddSJed Brown 
7664f491ddSJed Brown      This method has one explicit stage and four implicit stages.
7764f491ddSJed Brown 
7864f491ddSJed Brown      References:
7964f491ddSJed Brown      Kennedy and Carpenter 2003.
8064f491ddSJed Brown 
8164f491ddSJed Brown .seealso: TSARKIMEX
8264f491ddSJed Brown M*/
8364f491ddSJed Brown /*MC
8464f491ddSJed Brown      TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part.
8564f491ddSJed Brown 
8664f491ddSJed Brown      This method has one explicit stage and five implicit stages.
8764f491ddSJed Brown 
8864f491ddSJed Brown      References:
8964f491ddSJed Brown      Kennedy and Carpenter 2003.
9064f491ddSJed Brown 
9164f491ddSJed Brown .seealso: TSARKIMEX
9264f491ddSJed Brown M*/
9364f491ddSJed Brown 
948a381b04SJed Brown #undef __FUNCT__
958a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterAll"
968a381b04SJed Brown /*@C
978a381b04SJed Brown   TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX
988a381b04SJed Brown 
99fca742c7SJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
1008a381b04SJed Brown 
1018a381b04SJed Brown   Level: advanced
1028a381b04SJed Brown 
1038a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all
1048a381b04SJed Brown 
1058a381b04SJed Brown .seealso:  TSARKIMEXRegisterDestroy()
1068a381b04SJed Brown @*/
1078a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void)
1088a381b04SJed Brown {
1098a381b04SJed Brown   PetscErrorCode ierr;
1108a381b04SJed Brown 
1118a381b04SJed Brown   PetscFunctionBegin;
1128a381b04SJed Brown   if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0);
1138a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_TRUE;
1148a381b04SJed Brown 
1158a381b04SJed Brown   {
1168a381b04SJed Brown     const PetscReal
1178a381b04SJed Brown       A[3][3] = {{0,0,0},
1188a381b04SJed Brown                  {0.41421356237309504880,0,0},
1198a381b04SJed Brown                  {0.75,0.25,0}},
1208a381b04SJed Brown       At[3][3] = {{0,0,0},
1218a381b04SJed Brown                   {0.12132034355964257320,0.29289321881345247560,0},
12206db7b1cSJed Brown                   {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}},
12306db7b1cSJed Brown       binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}};
12406db7b1cSJed Brown     ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
1258a381b04SJed Brown   }
12606db7b1cSJed Brown   {                             /* Optimal for linear implicit part */
12703403c7fSJed Brown     const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),
128a3a57f36SJed Brown       A[3][3] = {{0,0,0},
129a3a57f36SJed Brown                  {2-s2,0,0},
130a3a57f36SJed Brown                  {(3-2*s2)/6,(3+2*s2)/6,0}},
131a3a57f36SJed Brown       At[3][3] = {{0,0,0},
132a3a57f36SJed Brown                   {1-1/s2,1-1/s2,0},
133cd652676SJed Brown                   {1/(2*s2),1/(2*s2),1-1/s2}},
134cd652676SJed Brown       binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}};
135cd652676SJed 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);
136a3a57f36SJed Brown   }
137a3a57f36SJed Brown   {
138a3a57f36SJed Brown     const PetscReal
139a3a57f36SJed Brown       A[4][4] = {{0,0,0,0},
1404040e9f2SJed Brown                  {1767732205903./2027836641118.,0,0,0},
1414040e9f2SJed Brown                  {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
1424040e9f2SJed Brown                  {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
143a3a57f36SJed Brown       At[4][4] = {{0,0,0,0},
1444040e9f2SJed Brown                   {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
1454040e9f2SJed Brown                   {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
1464040e9f2SJed Brown                   {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
1474040e9f2SJed Brown       binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
1484040e9f2SJed Brown                         {-18682724506714./9892148508045.,17870216137069./13817060693119.},
1494040e9f2SJed Brown                         {34259539580243./13192909600954.,-28141676662227./17317692491321.},
1504040e9f2SJed Brown                         {584795268549./6622622206610.,   2508943948391./7218656332882.}};
151cd652676SJed 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);
152a3a57f36SJed Brown   }
153a3a57f36SJed Brown   {
154a3a57f36SJed Brown     const PetscReal
155a3a57f36SJed Brown       A[6][6] = {{0,0,0,0,0,0},
156a3a57f36SJed Brown                  {1./2,0,0,0,0,0},
1574040e9f2SJed Brown                  {13861./62500.,6889./62500.,0,0,0,0},
1584040e9f2SJed Brown                  {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
1594040e9f2SJed Brown                  {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
1604040e9f2SJed Brown                  {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
161a3a57f36SJed Brown       At[6][6] = {{0,0,0,0,0,0},
162a3a57f36SJed Brown                   {1./4,1./4,0,0,0,0},
1634040e9f2SJed Brown                   {8611./62500.,-1743./31250.,1./4,0,0,0},
1644040e9f2SJed Brown                   {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
1654040e9f2SJed Brown                   {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
1664040e9f2SJed Brown                   {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
1674040e9f2SJed Brown       binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
168cd652676SJed Brown                         {0,0,0},
1694040e9f2SJed Brown                         {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
1704040e9f2SJed Brown                         {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
1714040e9f2SJed Brown                         {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
1724040e9f2SJed Brown                         {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
173cd652676SJed 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);
174a3a57f36SJed Brown   }
175a3a57f36SJed Brown   {
176a3a57f36SJed Brown     const PetscReal
177a3a57f36SJed Brown       A[8][8] = {{0,0,0,0,0,0,0,0},
178a3a57f36SJed Brown                  {41./100,0,0,0,0,0,0,0},
1794040e9f2SJed Brown                  {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
1804040e9f2SJed Brown                  {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
1814040e9f2SJed Brown                  {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
1824040e9f2SJed Brown                  {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
1834040e9f2SJed Brown                  {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
1844040e9f2SJed Brown                  {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
185a3a57f36SJed Brown       At[8][8] = {{0,0,0,0,0,0,0,0},
1864040e9f2SJed Brown                   {41./200.,41./200.,0,0,0,0,0,0},
1874040e9f2SJed Brown                   {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
1884040e9f2SJed Brown                   {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
1894040e9f2SJed Brown                   {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
1904040e9f2SJed Brown                   {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
1914040e9f2SJed Brown                   {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
1924040e9f2SJed Brown                   {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
1934040e9f2SJed Brown       binterpt[8][3] = {{-17674230611817./10670229744614. ,  43486358583215./12773830924787. , -9257016797708./5021505065439.},
194cd652676SJed Brown                         {0                               ,  0                              , 0                            },
195cd652676SJed Brown                         {0                               ,  0                              , 0                            },
1964040e9f2SJed Brown                         {65168852399939./7868540260826.  ,  -91478233927265./11067650958493., 26096422576131./11239449250142.},
1974040e9f2SJed Brown                         {15494834004392./5936557850923.  ,  -79368583304911./10890268929626., 92396832856987./20362823103730.},
1984040e9f2SJed Brown                         {-99329723586156./26959484932159.,  -12239297817655./9152339842473. , 30029262896817./10175596800299.},
1994040e9f2SJed Brown                         {-19024464361622./5461577185407. ,  115839755401235./10719374521269., -26136350496073./3983972220547.},
2004040e9f2SJed Brown                         {-6511271360970./6095937251113.  ,  5843115559534./2180450260947.   , -5289405421727./3760307252460. }};
201cd652676SJed 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);
202a3a57f36SJed Brown   }
203a3a57f36SJed Brown 
2048a381b04SJed Brown   PetscFunctionReturn(0);
2058a381b04SJed Brown }
2068a381b04SJed Brown 
2078a381b04SJed Brown #undef __FUNCT__
2088a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy"
2098a381b04SJed Brown /*@C
2108a381b04SJed Brown    TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
2118a381b04SJed Brown 
2128a381b04SJed Brown    Not Collective
2138a381b04SJed Brown 
2148a381b04SJed Brown    Level: advanced
2158a381b04SJed Brown 
2168a381b04SJed Brown .keywords: TSARKIMEX, register, destroy
2178a381b04SJed Brown .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic()
2188a381b04SJed Brown @*/
2198a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void)
2208a381b04SJed Brown {
2218a381b04SJed Brown   PetscErrorCode ierr;
2228a381b04SJed Brown   ARKTableauLink link;
2238a381b04SJed Brown 
2248a381b04SJed Brown   PetscFunctionBegin;
2258a381b04SJed Brown   while ((link = ARKTableauList)) {
2268a381b04SJed Brown     ARKTableau t = &link->tab;
2278a381b04SJed Brown     ARKTableauList = link->next;
2288a381b04SJed Brown     ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr);
229cd652676SJed Brown     ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr);
2308a381b04SJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
2318a381b04SJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
2328a381b04SJed Brown   }
2338a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_FALSE;
2348a381b04SJed Brown   PetscFunctionReturn(0);
2358a381b04SJed Brown }
2368a381b04SJed Brown 
2378a381b04SJed Brown #undef __FUNCT__
2388a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage"
2398a381b04SJed Brown /*@C
2408a381b04SJed Brown   TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
2418a381b04SJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
2428a381b04SJed Brown   when using static libraries.
2438a381b04SJed Brown 
2448a381b04SJed Brown   Input Parameter:
2458a381b04SJed Brown   path - The dynamic library path, or PETSC_NULL
2468a381b04SJed Brown 
2478a381b04SJed Brown   Level: developer
2488a381b04SJed Brown 
2498a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package
2508a381b04SJed Brown .seealso: PetscInitialize()
2518a381b04SJed Brown @*/
2528a381b04SJed Brown PetscErrorCode TSARKIMEXInitializePackage(const char path[])
2538a381b04SJed Brown {
2548a381b04SJed Brown   PetscErrorCode ierr;
2558a381b04SJed Brown 
2568a381b04SJed Brown   PetscFunctionBegin;
2578a381b04SJed Brown   if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0);
2588a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_TRUE;
2598a381b04SJed Brown   ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr);
2608a381b04SJed Brown   ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr);
2618a381b04SJed Brown   PetscFunctionReturn(0);
2628a381b04SJed Brown }
2638a381b04SJed Brown 
2648a381b04SJed Brown #undef __FUNCT__
2658a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage"
2668a381b04SJed Brown /*@C
2678a381b04SJed Brown   TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
2688a381b04SJed Brown   called from PetscFinalize().
2698a381b04SJed Brown 
2708a381b04SJed Brown   Level: developer
2718a381b04SJed Brown 
2728a381b04SJed Brown .keywords: Petsc, destroy, package
2738a381b04SJed Brown .seealso: PetscFinalize()
2748a381b04SJed Brown @*/
2758a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void)
2768a381b04SJed Brown {
2778a381b04SJed Brown   PetscErrorCode ierr;
2788a381b04SJed Brown 
2798a381b04SJed Brown   PetscFunctionBegin;
2808a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_FALSE;
2818a381b04SJed Brown   ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr);
2828a381b04SJed Brown   PetscFunctionReturn(0);
2838a381b04SJed Brown }
2848a381b04SJed Brown 
2858a381b04SJed Brown #undef __FUNCT__
2868a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister"
287cd652676SJed Brown /*@C
288cd652676SJed Brown    TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
289cd652676SJed Brown 
290cd652676SJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
291cd652676SJed Brown 
292cd652676SJed Brown    Input Parameters:
293cd652676SJed Brown +  name - identifier for method
294cd652676SJed Brown .  order - approximation order of method
295cd652676SJed Brown .  s - number of stages, this is the dimension of the matrices below
296cd652676SJed Brown .  At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
297cd652676SJed Brown .  bt - Butcher table for completing the stiff part of the step (dimension s; PETSC_NULL to use the last row of At)
298cd652676SJed Brown .  ct - Abscissa of each stiff stage (dimension s, PETSC_NULL to use row sums of At)
299cd652676SJed Brown .  A - Non-stiff stage coefficients (dimension s*s, row-major)
300cd652676SJed Brown .  b - Non-stiff step completion table (dimension s; PETSC_NULL to use last row of At)
301cd652676SJed Brown .  c - Non-stiff abscissa (dimension s; PETSC_NULL to use row sums of A)
302cd652676SJed Brown .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
303cd652676SJed Brown .  binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
304cd652676SJed Brown -  binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; PETSC_NULL to reuse binterpt)
305cd652676SJed Brown 
306cd652676SJed Brown    Notes:
307cd652676SJed Brown    Several ARK IMEX methods are provided, this function is only needed to create new methods.
308cd652676SJed Brown 
309cd652676SJed Brown    Level: advanced
310cd652676SJed Brown 
311cd652676SJed Brown .keywords: TS, register
312cd652676SJed Brown 
313cd652676SJed Brown .seealso: TSARKIMEX
314cd652676SJed Brown @*/
3158a381b04SJed Brown PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType name,PetscInt order,PetscInt s,
3168a381b04SJed Brown                                  const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
317cd652676SJed Brown                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
318cd652676SJed Brown                                  PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
3198a381b04SJed Brown {
3208a381b04SJed Brown   PetscErrorCode ierr;
3218a381b04SJed Brown   ARKTableauLink link;
3228a381b04SJed Brown   ARKTableau t;
3238a381b04SJed Brown   PetscInt i,j;
3248a381b04SJed Brown 
3258a381b04SJed Brown   PetscFunctionBegin;
3268a381b04SJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
327cd652676SJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
3288a381b04SJed Brown   t = &link->tab;
3298a381b04SJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
3308a381b04SJed Brown   t->order = order;
3318a381b04SJed Brown   t->s = s;
3328a381b04SJed 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);
3338a381b04SJed Brown   ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr);
3348a381b04SJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
3358a381b04SJed Brown   if (bt) {ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr);}
3368a381b04SJed Brown   else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
3378a381b04SJed Brown   if (b) {ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);}
3388a381b04SJed Brown   else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i];
3398a381b04SJed Brown   if (ct) {ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr);}
3408a381b04SJed 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];
3418a381b04SJed Brown   if (c) {ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr);}
3428a381b04SJed 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];
3434f385281SJed Brown   t->pinterp = pinterp;
344cd652676SJed Brown   ierr = PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);CHKERRQ(ierr);
345cd652676SJed Brown   ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
346cd652676SJed Brown   ierr = PetscMemcpy(t->binterp,binterp?binterp:binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
3478a381b04SJed Brown   link->next = ARKTableauList;
3488a381b04SJed Brown   ARKTableauList = link;
3498a381b04SJed Brown   PetscFunctionReturn(0);
3508a381b04SJed Brown }
3518a381b04SJed Brown 
3528a381b04SJed Brown #undef __FUNCT__
3538a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX"
3548a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts)
3558a381b04SJed Brown {
3568a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
3578a381b04SJed Brown   ARKTableau      tab  = ark->tableau;
3588a381b04SJed Brown   const PetscInt  s    = tab->s;
3598a381b04SJed Brown   const PetscReal *At  = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c;
360406d0ec2SJed Brown   PetscScalar     *w   = ark->work;
3618a381b04SJed Brown   Vec             *Y   = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,W = ark->Work,Z = ark->Z;
3628a381b04SJed Brown   SNES            snes;
3638a381b04SJed Brown   PetscInt        i,j,its,lits;
364cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
3658a381b04SJed Brown   PetscReal       h,t;
3668a381b04SJed Brown   PetscErrorCode  ierr;
3678a381b04SJed Brown 
3688a381b04SJed Brown   PetscFunctionBegin;
3698a381b04SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
370cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
371cdbf8f93SLisandro Dalcin   h = ts->time_step;
3728a381b04SJed Brown   t = ts->ptime;
3738a381b04SJed Brown 
3748a381b04SJed Brown   for (i=0; i<s; i++) {
3758a381b04SJed Brown     if (At[i*s+i] == 0) {           /* This stage is explicit */
3768a381b04SJed Brown       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
3778a381b04SJed Brown       for (j=0; j<i; j++) w[j] = -h*At[i*s+j];
3788a381b04SJed Brown       ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr);
3798a381b04SJed Brown       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
3808a381b04SJed Brown       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
3818a381b04SJed Brown     } else {
3828a381b04SJed Brown       ark->stage_time = t + h*ct[i];
3838a381b04SJed Brown       ark->shift = 1./(h*At[i*s+i]);
3848a381b04SJed Brown       /* Affine part */
3858a381b04SJed Brown       ierr = VecZeroEntries(W);CHKERRQ(ierr);
3868a381b04SJed Brown       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
3878a381b04SJed Brown       ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr);
3888a381b04SJed Brown       /* Ydot = shift*(Y-Z) */
3898a381b04SJed Brown       ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr);
3904f385281SJed Brown       for (j=0; j<i; j++) w[j] = -h*At[i*s+j];
3914f385281SJed Brown       ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr);
3928a381b04SJed Brown       /* Initial guess taken from last stage */
3938a381b04SJed Brown       ierr = VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);CHKERRQ(ierr);
3948a381b04SJed Brown       ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr);
3958a381b04SJed Brown       ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
3968a381b04SJed Brown       ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
3978a381b04SJed Brown       ts->nonlinear_its += its; ts->linear_its += lits;
3988a381b04SJed Brown     }
3998a381b04SJed Brown     ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
4004cc180ffSJed Brown     ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);
4014cc180ffSJed Brown     if (ark->imex) {
4028a381b04SJed Brown       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
4034cc180ffSJed Brown     } else {
4044cc180ffSJed Brown       ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr);
4054cc180ffSJed Brown     }
4068a381b04SJed Brown   }
407a17b68daSJed Brown   for (j=0; j<s; j++) w[j] = -h*bt[j];
4088a381b04SJed Brown   ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr);
4098a381b04SJed Brown   for (j=0; j<s; j++) w[j] = h*b[j];
4108a381b04SJed Brown   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
4118a381b04SJed Brown 
4128a381b04SJed Brown   ts->ptime += ts->time_step;
413cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
4148a381b04SJed Brown   ts->steps++;
4158a381b04SJed Brown   PetscFunctionReturn(0);
4168a381b04SJed Brown }
4178a381b04SJed Brown 
418cd652676SJed Brown #undef __FUNCT__
419cd652676SJed Brown #define __FUNCT__ "TSInterpolate_ARKIMEX"
420cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
421cd652676SJed Brown {
422cd652676SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
4234f385281SJed Brown   PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
4244f385281SJed Brown   PetscReal tt,t = (itime - ts->ptime)/ts->time_step + 1; /* In the interval [0,1] */
425cd652676SJed Brown   PetscScalar *bt,*b;
426cd652676SJed Brown   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
427cd652676SJed Brown   PetscErrorCode ierr;
428cd652676SJed Brown 
429cd652676SJed Brown   PetscFunctionBegin;
430cd652676SJed Brown   if (!Bt || !B) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
431cd652676SJed Brown   ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr);
432cd652676SJed Brown   for (i=0; i<s; i++) bt[i] = b[i] = 0;
4334f385281SJed Brown   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
434cd652676SJed Brown     for (i=0; i<s; i++) {
4354f385281SJed Brown       bt[i] += ts->time_step * Bt[i*pinterp+j] * tt * -1.0;
4364f385281SJed Brown       b[i]  += ts->time_step * B[i*pinterp+j] * tt;
437cd652676SJed Brown     }
438cd652676SJed Brown   }
439cd652676SJed 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");
440cd652676SJed Brown   ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr);
441cd652676SJed Brown   ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr);
442cd652676SJed Brown   ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr);
443cd652676SJed Brown   ierr = PetscFree2(bt,b);CHKERRQ(ierr);
444cd652676SJed Brown   PetscFunctionReturn(0);
445cd652676SJed Brown }
446cd652676SJed Brown 
4478a381b04SJed Brown /*------------------------------------------------------------*/
4488a381b04SJed Brown #undef __FUNCT__
4498a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX"
4508a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts)
4518a381b04SJed Brown {
4528a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
4538a381b04SJed Brown   PetscInt        s;
4548a381b04SJed Brown   PetscErrorCode  ierr;
4558a381b04SJed Brown 
4568a381b04SJed Brown   PetscFunctionBegin;
4578a381b04SJed Brown   if (!ark->tableau) PetscFunctionReturn(0);
4588a381b04SJed Brown    s = ark->tableau->s;
4598a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr);
4608a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr);
4618a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr);
4628a381b04SJed Brown   ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr);
4638a381b04SJed Brown   ierr = VecDestroy(&ark->Work);CHKERRQ(ierr);
4648a381b04SJed Brown   ierr = VecDestroy(&ark->Z);CHKERRQ(ierr);
4658a381b04SJed Brown   ierr = PetscFree(ark->work);CHKERRQ(ierr);
4668a381b04SJed Brown   PetscFunctionReturn(0);
4678a381b04SJed Brown }
4688a381b04SJed Brown 
4698a381b04SJed Brown #undef __FUNCT__
4708a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX"
4718a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
4728a381b04SJed Brown {
4738a381b04SJed Brown   PetscErrorCode  ierr;
4748a381b04SJed Brown 
4758a381b04SJed Brown   PetscFunctionBegin;
4768a381b04SJed Brown   ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
4778a381b04SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
4788a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);CHKERRQ(ierr);
4798a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);CHKERRQ(ierr);
4808a381b04SJed Brown   PetscFunctionReturn(0);
4818a381b04SJed Brown }
4828a381b04SJed Brown 
4838a381b04SJed Brown /*
4848a381b04SJed Brown   This defines the nonlinear equation that is to be solved with SNES
4858a381b04SJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
4868a381b04SJed Brown */
4878a381b04SJed Brown #undef __FUNCT__
4888a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX"
4898a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
4908a381b04SJed Brown {
4918a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
4928a381b04SJed Brown   PetscErrorCode ierr;
4938a381b04SJed Brown 
4948a381b04SJed Brown   PetscFunctionBegin;
4958a381b04SJed Brown   ierr = VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
4964cc180ffSJed Brown   ierr = TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,F,ark->imex);CHKERRQ(ierr);
4978a381b04SJed Brown   PetscFunctionReturn(0);
4988a381b04SJed Brown }
4998a381b04SJed Brown 
5008a381b04SJed Brown #undef __FUNCT__
5018a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX"
5028a381b04SJed Brown static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
5038a381b04SJed Brown {
5048a381b04SJed Brown   TS_ARKIMEX       *ark = (TS_ARKIMEX*)ts->data;
5058a381b04SJed Brown   PetscErrorCode ierr;
5068a381b04SJed Brown 
5078a381b04SJed Brown   PetscFunctionBegin;
5088a381b04SJed Brown   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
5098a381b04SJed Brown   ierr = TSComputeIJacobian(ts,ark->stage_time,X,ark->Ydot,ark->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
5108a381b04SJed Brown   PetscFunctionReturn(0);
5118a381b04SJed Brown }
5128a381b04SJed Brown 
5138a381b04SJed Brown #undef __FUNCT__
5148a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX"
5158a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
5168a381b04SJed Brown {
5178a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
5188a381b04SJed Brown   ARKTableau     tab  = ark->tableau;
5198a381b04SJed Brown   PetscInt       s = tab->s;
5208a381b04SJed Brown   PetscErrorCode ierr;
5218a381b04SJed Brown 
5228a381b04SJed Brown   PetscFunctionBegin;
5238a381b04SJed Brown   if (!ark->tableau) {
524e24355feSJed Brown     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
5258a381b04SJed Brown   }
5268a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr);
5278a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr);
5288a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr);
5298a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr);
5308a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr);
5318a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr);
5328a381b04SJed Brown   ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr);
5338a381b04SJed Brown   PetscFunctionReturn(0);
5348a381b04SJed Brown }
5358a381b04SJed Brown /*------------------------------------------------------------*/
5368a381b04SJed Brown 
5378a381b04SJed Brown #undef __FUNCT__
5388a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX"
5398a381b04SJed Brown static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts)
5408a381b04SJed Brown {
5414cc180ffSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
5428a381b04SJed Brown   PetscErrorCode ierr;
5438a381b04SJed Brown   char           arktype[256];
5448a381b04SJed Brown 
5458a381b04SJed Brown   PetscFunctionBegin;
5468a381b04SJed Brown   ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr);
5478a381b04SJed Brown   {
5488a381b04SJed Brown     ARKTableauLink link;
5498a381b04SJed Brown     PetscInt count,choice;
5508a381b04SJed Brown     PetscBool flg;
5518a381b04SJed Brown     const char **namelist;
5528a381b04SJed Brown     ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof arktype);CHKERRQ(ierr);
5538a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
5548a381b04SJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
5558a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
5568a381b04SJed Brown     ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr);
5578a381b04SJed Brown     ierr = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr);
5588a381b04SJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
5594cc180ffSJed Brown     flg = (PetscBool)!ark->imex;
5604cc180ffSJed Brown     ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
5614cc180ffSJed Brown     ark->imex = (PetscBool)!flg;
562d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
5638a381b04SJed Brown   }
5648a381b04SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
5658a381b04SJed Brown   PetscFunctionReturn(0);
5668a381b04SJed Brown }
5678a381b04SJed Brown 
5688a381b04SJed Brown #undef __FUNCT__
5698a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray"
5708a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
5718a381b04SJed Brown {
572257d2499SJed Brown   PetscErrorCode ierr;
5738a381b04SJed Brown   int i,left,count;
5748a381b04SJed Brown   char *p;
5758a381b04SJed Brown 
5768a381b04SJed Brown   PetscFunctionBegin;
5778a381b04SJed Brown   for (i=0,p=buf,left=(int)len; i<n; i++) {
578257d2499SJed Brown     ierr = PetscSNPrintf(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
5798a381b04SJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
5808a381b04SJed Brown     left -= count;
5818a381b04SJed Brown     p += count;
5828a381b04SJed Brown     *p++ = ' ';
5838a381b04SJed Brown   }
5848a381b04SJed Brown   p[i ? 0 : -1] = 0;
5858a381b04SJed Brown   PetscFunctionReturn(0);
5868a381b04SJed Brown }
5878a381b04SJed Brown 
5888a381b04SJed Brown #undef __FUNCT__
5898a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX"
5908a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
5918a381b04SJed Brown {
5928a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
5938a381b04SJed Brown   ARKTableau     tab = ark->tableau;
5948a381b04SJed Brown   PetscBool      iascii;
5958a381b04SJed Brown   PetscErrorCode ierr;
5968a381b04SJed Brown 
5978a381b04SJed Brown   PetscFunctionBegin;
5988a381b04SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5998a381b04SJed Brown   if (iascii) {
6008a381b04SJed Brown     const TSARKIMEXType arktype;
6018a381b04SJed Brown     char buf[512];
6028a381b04SJed Brown     ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr);
6038a381b04SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype);CHKERRQ(ierr);
6048a381b04SJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr);
60531f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf);CHKERRQ(ierr);
60631f6fcc0SJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
60731f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf);CHKERRQ(ierr);
6088a381b04SJed Brown   }
609d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
6108a381b04SJed Brown   PetscFunctionReturn(0);
6118a381b04SJed Brown }
6128a381b04SJed Brown 
6138a381b04SJed Brown #undef __FUNCT__
6148a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType"
6158a381b04SJed Brown /*@C
6168a381b04SJed Brown   TSARKIMEXSetType - Set the type of ARK IMEX scheme
6178a381b04SJed Brown 
6188a381b04SJed Brown   Logically collective
6198a381b04SJed Brown 
6208a381b04SJed Brown   Input Parameter:
6218a381b04SJed Brown +  ts - timestepping context
6228a381b04SJed Brown -  arktype - type of ARK-IMEX scheme
6238a381b04SJed Brown 
6248a381b04SJed Brown   Level: intermediate
6258a381b04SJed Brown 
6268a381b04SJed Brown .seealso: TSARKIMEXGetType()
6278a381b04SJed Brown @*/
6288a381b04SJed Brown PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType arktype)
6298a381b04SJed Brown {
6308a381b04SJed Brown   PetscErrorCode ierr;
6318a381b04SJed Brown 
6328a381b04SJed Brown   PetscFunctionBegin;
6338a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6348a381b04SJed Brown   ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,const TSARKIMEXType),(ts,arktype));CHKERRQ(ierr);
6358a381b04SJed Brown   PetscFunctionReturn(0);
6368a381b04SJed Brown }
6378a381b04SJed Brown 
6388a381b04SJed Brown #undef __FUNCT__
6398a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType"
6408a381b04SJed Brown /*@C
6418a381b04SJed Brown   TSARKIMEXGetType - Get the type of ARK IMEX scheme
6428a381b04SJed Brown 
6438a381b04SJed Brown   Logically collective
6448a381b04SJed Brown 
6458a381b04SJed Brown   Input Parameter:
6468a381b04SJed Brown .  ts - timestepping context
6478a381b04SJed Brown 
6488a381b04SJed Brown   Output Parameter:
6498a381b04SJed Brown .  arktype - type of ARK-IMEX scheme
6508a381b04SJed Brown 
6518a381b04SJed Brown   Level: intermediate
6528a381b04SJed Brown 
6538a381b04SJed Brown .seealso: TSARKIMEXGetType()
6548a381b04SJed Brown @*/
6558a381b04SJed Brown PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType *arktype)
6568a381b04SJed Brown {
6578a381b04SJed Brown   PetscErrorCode ierr;
6588a381b04SJed Brown 
6598a381b04SJed Brown   PetscFunctionBegin;
6608a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6618a381b04SJed Brown   ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,const TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr);
6628a381b04SJed Brown   PetscFunctionReturn(0);
6638a381b04SJed Brown }
6648a381b04SJed Brown 
6654cc180ffSJed Brown #undef __FUNCT__
6664cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit"
6674cc180ffSJed Brown /*@C
6684cc180ffSJed Brown   TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly
6694cc180ffSJed Brown 
6704cc180ffSJed Brown   Logically collective
6714cc180ffSJed Brown 
6724cc180ffSJed Brown   Input Parameter:
6734cc180ffSJed Brown +  ts - timestepping context
6744cc180ffSJed Brown -  flg - PETSC_TRUE for fully implicit
6754cc180ffSJed Brown 
6764cc180ffSJed Brown   Level: intermediate
6774cc180ffSJed Brown 
6784cc180ffSJed Brown .seealso: TSARKIMEXGetType()
6794cc180ffSJed Brown @*/
6804cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
6814cc180ffSJed Brown {
6824cc180ffSJed Brown   PetscErrorCode ierr;
6834cc180ffSJed Brown 
6844cc180ffSJed Brown   PetscFunctionBegin;
6854cc180ffSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6864cc180ffSJed Brown   ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
6874cc180ffSJed Brown   PetscFunctionReturn(0);
6884cc180ffSJed Brown }
6894cc180ffSJed Brown 
6908a381b04SJed Brown EXTERN_C_BEGIN
6918a381b04SJed Brown #undef __FUNCT__
6928a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX"
6938a381b04SJed Brown PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,const TSARKIMEXType *arktype)
6948a381b04SJed Brown {
6958a381b04SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
6968a381b04SJed Brown   PetscErrorCode ierr;
6978a381b04SJed Brown 
6988a381b04SJed Brown   PetscFunctionBegin;
6998a381b04SJed Brown   if (!ark->tableau) {ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);}
7008a381b04SJed Brown   *arktype = ark->tableau->name;
7018a381b04SJed Brown   PetscFunctionReturn(0);
7028a381b04SJed Brown }
7038a381b04SJed Brown #undef __FUNCT__
7048a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX"
7058a381b04SJed Brown PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,const TSARKIMEXType arktype)
7068a381b04SJed Brown {
7078a381b04SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
7088a381b04SJed Brown   PetscErrorCode ierr;
7098a381b04SJed Brown   PetscBool match;
7108a381b04SJed Brown   ARKTableauLink link;
7118a381b04SJed Brown 
7128a381b04SJed Brown   PetscFunctionBegin;
7138a381b04SJed Brown   if (ark->tableau) {
7148a381b04SJed Brown     ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr);
7158a381b04SJed Brown     if (match) PetscFunctionReturn(0);
7168a381b04SJed Brown   }
7178a381b04SJed Brown   for (link = ARKTableauList; link; link=link->next) {
7188a381b04SJed Brown     ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr);
7198a381b04SJed Brown     if (match) {
7208a381b04SJed Brown       ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
7218a381b04SJed Brown       ark->tableau = &link->tab;
7228a381b04SJed Brown       PetscFunctionReturn(0);
7238a381b04SJed Brown     }
7248a381b04SJed Brown   }
7258a381b04SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
7268a381b04SJed Brown   PetscFunctionReturn(0);
7278a381b04SJed Brown }
7284cc180ffSJed Brown #undef __FUNCT__
7294cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX"
7304cc180ffSJed Brown PetscErrorCode  TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
7314cc180ffSJed Brown {
7324cc180ffSJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
7334cc180ffSJed Brown 
7344cc180ffSJed Brown   PetscFunctionBegin;
7354cc180ffSJed Brown   ark->imex = (PetscBool)!flg;
7364cc180ffSJed Brown   PetscFunctionReturn(0);
7374cc180ffSJed Brown }
7388a381b04SJed Brown EXTERN_C_END
7398a381b04SJed Brown 
7408a381b04SJed Brown /* ------------------------------------------------------------ */
7418a381b04SJed Brown /*MC
7428a381b04SJed Brown       TSARKIMEX - ODE solver using Additive Runge-Kutta IMEX schemes
7438a381b04SJed Brown 
744fca742c7SJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
745fca742c7SJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
746fca742c7SJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
747fca742c7SJed Brown 
748fca742c7SJed Brown   Notes:
749c8058688SBarry Smith   The default is TSARKIMEX2E, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type
750c8058688SBarry Smith 
751fca742c7SJed Brown   This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X).
752fca742c7SJed Brown 
7538a381b04SJed Brown   Level: beginner
7548a381b04SJed Brown 
755c8058688SBarry Smith .seealso:  TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3,
756*26885f59SBarry Smith            TSARKIMEX4, TSARKIMEX5, TSARKIMEXType, TSARKIMEXRegister()
7578a381b04SJed Brown 
7588a381b04SJed Brown M*/
7598a381b04SJed Brown EXTERN_C_BEGIN
7608a381b04SJed Brown #undef __FUNCT__
7618a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX"
7628a381b04SJed Brown PetscErrorCode  TSCreate_ARKIMEX(TS ts)
7638a381b04SJed Brown {
7648a381b04SJed Brown   TS_ARKIMEX       *th;
7658a381b04SJed Brown   PetscErrorCode ierr;
7668a381b04SJed Brown 
7678a381b04SJed Brown   PetscFunctionBegin;
7688a381b04SJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
7698a381b04SJed Brown   ierr = TSARKIMEXInitializePackage(PETSC_NULL);CHKERRQ(ierr);
7708a381b04SJed Brown #endif
7718a381b04SJed Brown 
7728a381b04SJed Brown   ts->ops->reset          = TSReset_ARKIMEX;
7738a381b04SJed Brown   ts->ops->destroy        = TSDestroy_ARKIMEX;
7748a381b04SJed Brown   ts->ops->view           = TSView_ARKIMEX;
7758a381b04SJed Brown   ts->ops->setup          = TSSetUp_ARKIMEX;
7768a381b04SJed Brown   ts->ops->step           = TSStep_ARKIMEX;
777cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_ARKIMEX;
7788a381b04SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
7798a381b04SJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
7808a381b04SJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;
7818a381b04SJed Brown 
7828a381b04SJed Brown   ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr);
7838a381b04SJed Brown   ts->data = (void*)th;
7844cc180ffSJed Brown   th->imex = PETSC_TRUE;
7858a381b04SJed Brown 
7868a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr);
7878a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr);
7884cc180ffSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","TSARKIMEXSetFullyImplicit_ARKIMEX",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr);
7898a381b04SJed Brown   PetscFunctionReturn(0);
7908a381b04SJed Brown }
7918a381b04SJed Brown EXTERN_C_END
792