xref: /petsc/src/ts/impls/implicit/irk/irk.c (revision 8af18dd833b18fc2df02377e75d76d273dc343a9)
1 /*
2   Code for timestepping with implicit Runge-Kutta method
3 
4   Notes:
5   The general system is written as
6 
7   F(t,U,Udot) = 0
8 
9 */
10 #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
11 #include <petscdm.h>
12 #include <petscdt.h>
13 
14 static TSIRKType         TSIRKDefault = TSIRKGAUSS;
15 static PetscBool         TSIRKRegisterAllCalled;
16 static PetscBool         TSIRKPackageInitialized;
17 static PetscFunctionList TSIRKList;
18 
19 struct _IRKTableau{
20   PetscReal   *A,*b,*c;
21   PetscScalar *A_inv,*A_inv_rowsum,*I_s;
22   PetscReal   *binterp;   /* Dense output formula */
23 };
24 
25 typedef struct _IRKTableau *IRKTableau;
26 
27 typedef struct {
28   char         *method_name;
29   PetscInt     order;            /* Classical approximation order of the method */
30   PetscInt     nstages;          /* Number of stages */
31   PetscBool    stiffly_accurate;
32   PetscInt     pinterp;          /* Interpolation order */
33   IRKTableau   tableau;
34   Vec          U0;               /* Backup vector */
35   Vec          Z;                /* Combined stage vector */
36   Vec          *Y;               /* States computed during the step */
37   Vec          Ydot;             /* Work vector holding time derivatives during residual evaluation */
38   Vec          U;                /* U is used to compute Ydot = shift(Y-U) */
39   Vec          *YdotI;           /* Work vectors to hold the residual evaluation */
40   Mat          TJ;               /* KAIJ matrix for the Jacobian of the combined system */
41   PetscScalar  *work;            /* Scalar work */
42   TSStepStatus status;
43   PetscBool    rebuild_completion;
44   PetscReal    ccfl;
45 } TS_IRK;
46 
47 /*@C
48    TSIRKTableauCreate - create the tableau for TSIRK and provide the entries
49 
50    Not Collective
51 
52    Input Parameters:
53 +  ts - timestepping context
54 .  nstages - number of stages, this is the dimension of the matrices below
55 .  A - stage coefficients (dimension nstages*nstages, row-major)
56 .  b - step completion table (dimension nstages)
57 .  c - abscissa (dimension nstages)
58 .  binterp - coefficients of the interpolation formula (dimension nstages)
59 .  A_inv - inverse of A (dimension nstages*nstages, row-major)
60 .  A_inv_rowsum - row sum of the inverse of A (dimension nstages)
61 -  I_s - identity matrix (dimension nstages*nstages)
62 
63    Level: advanced
64 
65 .seealso: TSIRK, TSIRKRegister()
66 @*/
67 PetscErrorCode TSIRKTableauCreate(TS ts,PetscInt nstages,const PetscReal *A,const PetscReal *b,const PetscReal *c,const PetscReal *binterp,const PetscScalar *A_inv,const PetscScalar *A_inv_rowsum,const PetscScalar *I_s)
68 {
69   TS_IRK         *irk = (TS_IRK*)ts->data;
70   IRKTableau     tab = irk->tableau;
71   PetscErrorCode ierr;
72 
73   PetscFunctionBegin;
74   irk->order = nstages;
75   ierr = PetscMalloc3(PetscSqr(nstages),&tab->A,PetscSqr(nstages),&tab->A_inv,PetscSqr(nstages),&tab->I_s);CHKERRQ(ierr);
76   ierr = PetscMalloc4(nstages,&tab->b,nstages,&tab->c,nstages,&tab->binterp,nstages,&tab->A_inv_rowsum);CHKERRQ(ierr);
77   ierr = PetscArraycpy(tab->A,A,PetscSqr(nstages));CHKERRQ(ierr);
78   ierr = PetscArraycpy(tab->b,b,nstages);CHKERRQ(ierr);
79   ierr = PetscArraycpy(tab->c,c,nstages);CHKERRQ(ierr);
80   /* optional coefficient arrays */
81   if (binterp) {
82     ierr = PetscArraycpy(tab->binterp,binterp,nstages);CHKERRQ(ierr);
83   }
84   if (A_inv) {
85     ierr = PetscArraycpy(tab->A_inv,A_inv,PetscSqr(nstages));CHKERRQ(ierr);
86   }
87   if (A_inv_rowsum) {
88     ierr = PetscArraycpy(tab->A_inv_rowsum,A_inv_rowsum,nstages);CHKERRQ(ierr);
89   }
90   if (I_s) {
91     ierr = PetscArraycpy(tab->I_s,I_s,PetscSqr(nstages));CHKERRQ(ierr);
92   }
93   PetscFunctionReturn(0);
94 }
95 
96 /* Arrays should be freed with PetscFree3(A,b,c) */
97 static PetscErrorCode TSIRKCreate_Gauss(TS ts)
98 {
99   PetscInt       nstages;
100   PetscReal      *gauss_A_real,*gauss_b,*b,*gauss_c;
101   PetscScalar    *gauss_A,*gauss_A_inv,*gauss_A_inv_rowsum,*I_s;
102   PetscScalar    *G0,*G1;
103   PetscInt       i,j;
104   Mat            G0mat,G1mat,Amat;
105   PetscErrorCode ierr;
106 
107   PetscFunctionBegin;
108   ierr = TSIRKGetNumStages(ts,&nstages);CHKERRQ(ierr);
109   ierr = PetscMalloc3(PetscSqr(nstages),&gauss_A_real,nstages,&gauss_b,nstages,&gauss_c);CHKERRQ(ierr);
110   ierr = PetscMalloc4(PetscSqr(nstages),&gauss_A,PetscSqr(nstages),&gauss_A_inv,nstages,&gauss_A_inv_rowsum,PetscSqr(nstages),&I_s);CHKERRQ(ierr);
111   ierr = PetscMalloc3(nstages,&b,PetscSqr(nstages),&G0,PetscSqr(nstages),&G1);CHKERRQ(ierr);
112   ierr = PetscDTGaussQuadrature(nstages,0.,1.,gauss_c,b);CHKERRQ(ierr);
113   for (i=0; i<nstages; i++) gauss_b[i] = b[i]; /* copy to possibly-complex array */
114 
115   /* A^T = G0^{-1} G1 */
116   for (i=0; i<nstages; i++) {
117     for (j=0; j<nstages; j++) {
118       G0[i*nstages+j] = PetscPowRealInt(gauss_c[i],j);
119       G1[i*nstages+j] = PetscPowRealInt(gauss_c[i],j+1)/(j+1);
120     }
121   }
122   /* The arrays above are row-aligned, but we create dense matrices as the transpose */
123   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,G0,&G0mat);CHKERRQ(ierr);
124   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,G1,&G1mat);CHKERRQ(ierr);
125   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,gauss_A,&Amat);CHKERRQ(ierr);
126   ierr = MatLUFactor(G0mat,NULL,NULL,NULL);CHKERRQ(ierr);
127   ierr = MatMatSolve(G0mat,G1mat,Amat);CHKERRQ(ierr);
128   ierr = MatTranspose(Amat,MAT_INPLACE_MATRIX,&Amat);CHKERRQ(ierr);
129   for (i=0; i<nstages; i++)
130     for (j=0; j<nstages; j++)
131       gauss_A_real[i*nstages+j] = PetscRealPart(gauss_A[i*nstages+j]);
132 
133   ierr = MatDestroy(&G0mat);CHKERRQ(ierr);
134   ierr = MatDestroy(&G1mat);CHKERRQ(ierr);
135   ierr = MatDestroy(&Amat);CHKERRQ(ierr);
136   ierr = PetscFree3(b,G0,G1);CHKERRQ(ierr);
137 
138   {/* Invert A */
139     /* PETSc does not provide a routine to calculate the inverse of a general matrix.
140      * To get the inverse of A, we form a sequential BAIJ matrix from it, consisting of a single block with block size
141      * equal to the dimension of A, and then use MatInvertBlockDiagonal(). */
142     Mat               A_baij;
143     PetscInt          idxm[1]={0},idxn[1]={0};
144     const PetscScalar *A_inv;
145 
146     ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,nstages,nstages,nstages,1,NULL,&A_baij);CHKERRQ(ierr);
147     ierr = MatSetOption(A_baij,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
148     ierr = MatSetValuesBlocked(A_baij,1,idxm,1,idxn,gauss_A,INSERT_VALUES);CHKERRQ(ierr);
149     ierr = MatAssemblyBegin(A_baij,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150     ierr = MatAssemblyEnd(A_baij,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
151     ierr = MatInvertBlockDiagonal(A_baij,&A_inv);CHKERRQ(ierr);
152     ierr = PetscMemcpy(gauss_A_inv,A_inv,nstages*nstages*sizeof(PetscScalar));CHKERRQ(ierr);
153     ierr = MatDestroy(&A_baij);CHKERRQ(ierr);
154   }
155 
156   /* Compute row sums A_inv_rowsum and identity I_s */
157   for (i=0; i<nstages; i++) {
158     gauss_A_inv_rowsum[i] = 0;
159     for (j=0; j<nstages; j++) {
160       gauss_A_inv_rowsum[i] += gauss_A_inv[i+nstages*j];
161       I_s[i+nstages*j] = 1.*(i == j);
162     }
163   }
164   ierr = TSIRKTableauCreate(ts,nstages,gauss_A_real,gauss_b,gauss_c,NULL,gauss_A_inv,gauss_A_inv_rowsum,I_s);CHKERRQ(ierr);
165   ierr = PetscFree3(gauss_A_real,gauss_b,gauss_c);CHKERRQ(ierr);
166   ierr = PetscFree4(gauss_A,gauss_A_inv,gauss_A_inv_rowsum,I_s);CHKERRQ(ierr);
167   PetscFunctionReturn(0);
168 }
169 
170 /*@C
171    TSIRKRegister -  adds a TSIRK implementation
172 
173    Not Collective
174 
175    Input Parameters:
176 +  sname - name of user-defined IRK scheme
177 -  function - function to create method context
178 
179    Notes:
180    TSIRKRegister() may be called multiple times to add several user-defined families.
181 
182    Sample usage:
183 .vb
184    TSIRKRegister("my_scheme",MySchemeCreate);
185 .ve
186 
187    Then, your scheme can be chosen with the procedural interface via
188 $     TSIRKSetType(ts,"my_scheme")
189    or at runtime via the option
190 $     -ts_irk_type my_scheme
191 
192    Level: advanced
193 
194 .seealso: TSIRKRegisterAll()
195 @*/
196 PetscErrorCode TSIRKRegister(const char sname[],PetscErrorCode (*function)(TS))
197 {
198   PetscErrorCode ierr;
199 
200   PetscFunctionBegin;
201   ierr = TSIRKInitializePackage();CHKERRQ(ierr);
202   ierr = PetscFunctionListAdd(&TSIRKList,sname,function);CHKERRQ(ierr);
203   PetscFunctionReturn(0);
204 }
205 
206 /*@C
207   TSIRKRegisterAll - Registers all of the implicit Runge-Kutta methods in TSIRK
208 
209   Not Collective, but should be called by all processes which will need the schemes to be registered
210 
211   Level: advanced
212 
213 .seealso:  TSIRKRegisterDestroy()
214 @*/
215 PetscErrorCode TSIRKRegisterAll(void)
216 {
217   PetscErrorCode ierr;
218 
219   PetscFunctionBegin;
220   if (TSIRKRegisterAllCalled) PetscFunctionReturn(0);
221   TSIRKRegisterAllCalled = PETSC_TRUE;
222 
223   ierr = TSIRKRegister(TSIRKGAUSS,TSIRKCreate_Gauss);CHKERRQ(ierr);
224   PetscFunctionReturn(0);
225 }
226 
227 /*@C
228    TSIRKRegisterDestroy - Frees the list of schemes that were registered by TSIRKRegister().
229 
230    Not Collective
231 
232    Level: advanced
233 
234 .seealso: TSIRKRegister(), TSIRKRegisterAll()
235 @*/
236 PetscErrorCode TSIRKRegisterDestroy(void)
237 {
238   PetscFunctionBegin;
239   TSIRKRegisterAllCalled = PETSC_FALSE;
240   PetscFunctionReturn(0);
241 }
242 
243 /*@C
244   TSIRKInitializePackage - This function initializes everything in the TSIRK package. It is called
245   from TSInitializePackage().
246 
247   Level: developer
248 
249 .seealso: PetscInitialize()
250 @*/
251 PetscErrorCode TSIRKInitializePackage(void)
252 {
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256   if (TSIRKPackageInitialized) PetscFunctionReturn(0);
257   TSIRKPackageInitialized = PETSC_TRUE;
258   ierr = TSIRKRegisterAll();CHKERRQ(ierr);
259   ierr = PetscRegisterFinalize(TSIRKFinalizePackage);CHKERRQ(ierr);
260   PetscFunctionReturn(0);
261 }
262 
263 /*@C
264   TSIRKFinalizePackage - This function destroys everything in the TSIRK package. It is
265   called from PetscFinalize().
266 
267   Level: developer
268 
269 .seealso: PetscFinalize()
270 @*/
271 PetscErrorCode TSIRKFinalizePackage(void)
272 {
273   PetscErrorCode ierr;
274 
275   PetscFunctionBegin;
276   ierr = PetscFunctionListDestroy(&TSIRKList);CHKERRQ(ierr);
277   TSIRKPackageInitialized = PETSC_FALSE;
278   PetscFunctionReturn(0);
279 }
280 
281 /*
282  This function can be called before or after ts->vec_sol has been updated.
283 */
284 static PetscErrorCode TSEvaluateStep_IRK(TS ts,PetscInt order,Vec U,PetscBool *done)
285 {
286   TS_IRK         *irk = (TS_IRK*)ts->data;
287   IRKTableau     tab = irk->tableau;
288   Vec            *YdotI = irk->YdotI;
289   PetscScalar    *w = irk->work;
290   PetscReal      h;
291   PetscInt       j;
292   PetscErrorCode ierr;
293 
294   PetscFunctionBegin;
295   switch (irk->status) {
296   case TS_STEP_INCOMPLETE:
297   case TS_STEP_PENDING:
298     h = ts->time_step; break;
299   case TS_STEP_COMPLETE:
300     h = ts->ptime - ts->ptime_prev; break;
301   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
302   }
303 
304   ierr = VecCopy(ts->vec_sol,U);CHKERRQ(ierr);
305   for (j=0; j<irk->nstages; j++) w[j] = h*tab->b[j];
306   ierr = VecMAXPY(U,irk->nstages,w,YdotI);CHKERRQ(ierr);
307   PetscFunctionReturn(0);
308 }
309 
310 static PetscErrorCode TSRollBack_IRK(TS ts)
311 {
312   TS_IRK         *irk = (TS_IRK*)ts->data;
313   PetscErrorCode ierr;
314 
315   PetscFunctionBegin;
316   ierr = VecCopy(irk->U0,ts->vec_sol);CHKERRQ(ierr);
317   PetscFunctionReturn(0);
318 }
319 
320 static PetscErrorCode TSStep_IRK(TS ts)
321 {
322   TS_IRK          *irk = (TS_IRK*)ts->data;
323   IRKTableau      tab = irk->tableau;
324   PetscScalar     *A_inv = tab->A_inv,*A_inv_rowsum = tab->A_inv_rowsum;
325   const PetscInt  nstages = irk->nstages;
326   SNES            snes;
327   PetscInt        i,j,its,lits,bs;
328   TSAdapt         adapt;
329   PetscInt        rejections = 0;
330   PetscBool       accept = PETSC_TRUE;
331   PetscReal       next_time_step = ts->time_step;
332   PetscErrorCode  ierr;
333 
334   PetscFunctionBegin;
335   if (!ts->steprollback) {
336     ierr = VecCopy(ts->vec_sol,irk->U0);CHKERRQ(ierr);
337   }
338   ierr = VecGetBlockSize(ts->vec_sol,&bs);CHKERRQ(ierr);
339   for (i=0; i<nstages; i++) {
340     ierr = VecStrideScatter(ts->vec_sol,i*bs,irk->Z,INSERT_VALUES);CHKERRQ(ierr);
341   }
342 
343   irk->status = TS_STEP_INCOMPLETE;
344   while (!ts->reason && irk->status != TS_STEP_COMPLETE) {
345     ierr = VecCopy(ts->vec_sol,irk->U);CHKERRQ(ierr);
346     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
347     ierr = SNESSolve(snes,NULL,irk->Z);CHKERRQ(ierr);
348     ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
349     ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
350     ts->snes_its += its; ts->ksp_its += lits;
351     ierr = VecStrideGatherAll(irk->Z,irk->Y,INSERT_VALUES);CHKERRQ(ierr);
352     for (i=0; i<nstages; i++) {
353       ierr = VecZeroEntries(irk->YdotI[i]);CHKERRQ(ierr);
354       for (j=0; j<nstages; j++) {
355         ierr = VecAXPY(irk->YdotI[i],A_inv[i+j*nstages]/ts->time_step,irk->Y[j]);CHKERRQ(ierr);
356       }
357       ierr = VecAXPY(irk->YdotI[i],-A_inv_rowsum[i]/ts->time_step,irk->U);CHKERRQ(ierr);
358     }
359     irk->status = TS_STEP_INCOMPLETE;
360     ierr = TSEvaluateStep_IRK(ts,irk->order,ts->vec_sol,NULL);CHKERRQ(ierr);
361     irk->status = TS_STEP_PENDING;
362     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
363     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
364     irk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
365     if (!accept) {
366       ierr = TSRollBack_IRK(ts);CHKERRQ(ierr);
367       ts->time_step = next_time_step;
368       goto reject_step;
369     }
370 
371     ts->ptime += ts->time_step;
372     ts->time_step = next_time_step;
373     break;
374   reject_step:
375     ts->reject++; accept = PETSC_FALSE;
376     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
377       ts->reason = TS_DIVERGED_STEP_REJECTED;
378       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
379     }
380   }
381   PetscFunctionReturn(0);
382 }
383 
384 static PetscErrorCode TSInterpolate_IRK(TS ts,PetscReal itime,Vec U)
385 {
386   TS_IRK          *irk = (TS_IRK*)ts->data;
387   PetscInt        nstages = irk->nstages,pinterp = irk->pinterp,i,j;
388   PetscReal       h;
389   PetscReal       tt,t;
390   PetscScalar     *bt;
391   const PetscReal *B = irk->tableau->binterp;
392   PetscErrorCode  ierr;
393 
394   PetscFunctionBegin;
395   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSIRK %s does not have an interpolation formula",irk->method_name);
396   switch (irk->status) {
397   case TS_STEP_INCOMPLETE:
398   case TS_STEP_PENDING:
399     h = ts->time_step;
400     t = (itime - ts->ptime)/h;
401     break;
402   case TS_STEP_COMPLETE:
403     h = ts->ptime - ts->ptime_prev;
404     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
405     break;
406   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
407   }
408   ierr = PetscMalloc1(nstages,&bt);CHKERRQ(ierr);
409   for (i=0; i<nstages; i++) bt[i] = 0;
410   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
411     for (i=0; i<nstages; i++) {
412       bt[i] += h * B[i*pinterp+j] * tt;
413     }
414   }
415   ierr = VecMAXPY(U,nstages,bt,irk->YdotI);CHKERRQ(ierr);
416   PetscFunctionReturn(0);
417 }
418 
419 static PetscErrorCode TSIRKTableauReset(TS ts)
420 {
421   TS_IRK         *irk = (TS_IRK*)ts->data;
422   IRKTableau     tab = irk->tableau;
423   PetscErrorCode ierr;
424 
425   PetscFunctionBegin;
426   if (!tab) PetscFunctionReturn(0);
427   ierr = PetscFree3(tab->A,tab->A_inv,tab->I_s);CHKERRQ(ierr);
428   ierr = PetscFree4(tab->b,tab->c,tab->binterp,tab->A_inv_rowsum);CHKERRQ(ierr);
429   PetscFunctionReturn(0);
430 }
431 
432 static PetscErrorCode TSReset_IRK(TS ts)
433 {
434   TS_IRK         *irk = (TS_IRK*)ts->data;
435   PetscErrorCode ierr;
436 
437   PetscFunctionBegin;
438   ierr = TSIRKTableauReset(ts);CHKERRQ(ierr);
439   if (irk->tableau) {
440     ierr = PetscFree(irk->tableau);CHKERRQ(ierr);
441   }
442   if (irk->method_name) {
443     ierr = PetscFree(irk->method_name);CHKERRQ(ierr);
444   }
445   if (irk->work) {
446     ierr = PetscFree(irk->work);CHKERRQ(ierr);
447   }
448   ierr = VecDestroyVecs(irk->nstages,&irk->Y);CHKERRQ(ierr);
449   ierr = VecDestroyVecs(irk->nstages,&irk->YdotI);CHKERRQ(ierr);
450   ierr = VecDestroy(&irk->Ydot);CHKERRQ(ierr);
451   ierr = VecDestroy(&irk->Z);CHKERRQ(ierr);
452   ierr = VecDestroy(&irk->U);CHKERRQ(ierr);
453   ierr = VecDestroy(&irk->U0);CHKERRQ(ierr);
454   ierr = MatDestroy(&irk->TJ);CHKERRQ(ierr);
455   PetscFunctionReturn(0);
456 }
457 
458 static PetscErrorCode TSIRKGetVecs(TS ts,DM dm,Vec *U)
459 {
460   TS_IRK         *irk = (TS_IRK*)ts->data;
461   PetscErrorCode ierr;
462 
463   PetscFunctionBegin;
464   if (U) {
465     if (dm && dm != ts->dm) {
466       ierr = DMGetNamedGlobalVector(dm,"TSIRK_U",U);CHKERRQ(ierr);
467     } else *U = irk->U;
468   }
469   PetscFunctionReturn(0);
470 }
471 
472 static PetscErrorCode TSIRKRestoreVecs(TS ts,DM dm,Vec *U)
473 {
474   PetscErrorCode ierr;
475 
476   PetscFunctionBegin;
477   if (U) {
478     if (dm && dm != ts->dm) {
479       ierr = DMRestoreNamedGlobalVector(dm,"TSIRK_U",U);CHKERRQ(ierr);
480     }
481   }
482   PetscFunctionReturn(0);
483 }
484 
485 /*
486   This defines the nonlinear equations that is to be solved with SNES
487     G[e\otimes t + C*dt, Z, Zdot] = 0
488     Zdot = (In \otimes S)*Z - (In \otimes Se) U
489   where S = 1/(dt*A)
490 */
491 static PetscErrorCode SNESTSFormFunction_IRK(SNES snes,Vec ZC,Vec FC,TS ts)
492 {
493   TS_IRK            *irk = (TS_IRK*)ts->data;
494   IRKTableau        tab  = irk->tableau;
495   const PetscInt    nstages = irk->nstages;
496   const PetscReal   *c = tab->c;
497   const PetscScalar *A_inv = tab->A_inv,*A_inv_rowsum = tab->A_inv_rowsum;
498   DM                dm,dmsave;
499   Vec               U,*YdotI = irk->YdotI,Ydot = irk->Ydot,*Y = irk->Y;
500   PetscReal         h = ts->time_step;
501   PetscInt          i,j;
502   PetscErrorCode    ierr;
503 
504   PetscFunctionBegin;
505   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
506   ierr = TSIRKGetVecs(ts,dm,&U);CHKERRQ(ierr);
507   ierr = VecStrideGatherAll(ZC,Y,INSERT_VALUES);CHKERRQ(ierr);
508   dmsave = ts->dm;
509   ts->dm = dm;
510   for (i=0; i<nstages; i++) {
511     ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
512     for (j=0; j<nstages; j++) {
513       ierr = VecAXPY(Ydot,A_inv[j*nstages+i]/h,Y[j]);CHKERRQ(ierr);
514     }
515     ierr = VecAXPY(Ydot,-A_inv_rowsum[i]/h,U);CHKERRQ(ierr); /* Ydot = (S \otimes In)*Z - (Se \otimes In) U */
516     ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step*c[i],Y[i],Ydot,YdotI[i],PETSC_FALSE);CHKERRQ(ierr);
517   }
518   ierr = VecStrideScatterAll(YdotI,FC,INSERT_VALUES);CHKERRQ(ierr);
519   ts->dm = dmsave;
520   ierr   = TSIRKRestoreVecs(ts,dm,&U);CHKERRQ(ierr);
521   PetscFunctionReturn(0);
522 }
523 
524 /*
525    For explicit ODE, the Jacobian is
526      JC = I_n \otimes S - J \otimes I_s
527    For DAE, the Jacobian is
528      JC = M_n \otimes S - J \otimes I_s
529 */
530 static PetscErrorCode SNESTSFormJacobian_IRK(SNES snes,Vec ZC,Mat JC,Mat JCpre,TS ts)
531 {
532   TS_IRK          *irk = (TS_IRK*)ts->data;
533   IRKTableau      tab  = irk->tableau;
534   const PetscInt  nstages = irk->nstages;
535   const PetscReal *c = tab->c;
536   DM              dm,dmsave;
537   Vec             *Y = irk->Y,Ydot = irk->Ydot;
538   Mat             J;
539   PetscScalar     *S;
540   PetscInt        i,j,bs;
541   PetscErrorCode  ierr;
542 
543   PetscFunctionBegin;
544   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
545   /* irk->Ydot has already been computed in SNESTSFormFunction_IRK (SNES guarantees this) */
546   dmsave = ts->dm;
547   ts->dm = dm;
548   ierr = VecGetBlockSize(Y[nstages-1],&bs);CHKERRQ(ierr);
549   if (ts->equation_type <= TS_EQ_ODE_EXPLICIT) { /* Support explicit formulas only */
550     ierr = VecStrideGather(ZC,(nstages-1)*bs,Y[nstages-1],INSERT_VALUES);CHKERRQ(ierr);
551     ierr = MatKAIJGetAIJ(JC,&J);CHKERRQ(ierr);
552     ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step*c[nstages-1],Y[nstages-1],Ydot,0,J,J,PETSC_FALSE);CHKERRQ(ierr);
553     ierr = MatKAIJGetS(JC,NULL,NULL,&S);CHKERRQ(ierr);
554     for (i=0; i<nstages; i++)
555       for (j=0; j<nstages; j++)
556         S[i+nstages*j] = tab->A_inv[i+nstages*j]/ts->time_step;
557     ierr = MatKAIJRestoreS(JC,&S);CHKERRQ(ierr);
558   } else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSIRK %s does not support implicit formula",irk->method_name); /* TODO: need the mass matrix for DAE  */
559   ts->dm = dmsave;
560   PetscFunctionReturn(0);
561 }
562 
563 static PetscErrorCode DMCoarsenHook_TSIRK(DM fine,DM coarse,void *ctx)
564 {
565   PetscFunctionBegin;
566   PetscFunctionReturn(0);
567 }
568 
569 static PetscErrorCode DMRestrictHook_TSIRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
570 {
571   TS             ts = (TS)ctx;
572   Vec            U,U_c;
573   PetscErrorCode ierr;
574 
575   PetscFunctionBegin;
576   ierr = TSIRKGetVecs(ts,fine,&U);CHKERRQ(ierr);
577   ierr = TSIRKGetVecs(ts,coarse,&U_c);CHKERRQ(ierr);
578   ierr = MatRestrict(restrct,U,U_c);CHKERRQ(ierr);
579   ierr = VecPointwiseMult(U_c,rscale,U_c);CHKERRQ(ierr);
580   ierr = TSIRKRestoreVecs(ts,fine,&U);CHKERRQ(ierr);
581   ierr = TSIRKRestoreVecs(ts,coarse,&U_c);CHKERRQ(ierr);
582   PetscFunctionReturn(0);
583 }
584 
585 static PetscErrorCode DMSubDomainHook_TSIRK(DM dm,DM subdm,void *ctx)
586 {
587   PetscFunctionBegin;
588   PetscFunctionReturn(0);
589 }
590 
591 static PetscErrorCode DMSubDomainRestrictHook_TSIRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
592 {
593   TS             ts = (TS)ctx;
594   Vec            U,U_c;
595   PetscErrorCode ierr;
596 
597   PetscFunctionBegin;
598   ierr = TSIRKGetVecs(ts,dm,&U);CHKERRQ(ierr);
599   ierr = TSIRKGetVecs(ts,subdm,&U_c);CHKERRQ(ierr);
600 
601   ierr = VecScatterBegin(gscat,U,U_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
602   ierr = VecScatterEnd(gscat,U,U_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
603 
604   ierr = TSIRKRestoreVecs(ts,dm,&U);CHKERRQ(ierr);
605   ierr = TSIRKRestoreVecs(ts,subdm,&U_c);CHKERRQ(ierr);
606   PetscFunctionReturn(0);
607 }
608 
609 static PetscErrorCode TSSetUp_IRK(TS ts)
610 {
611   TS_IRK         *irk = (TS_IRK*)ts->data;
612   IRKTableau     tab = irk->tableau;
613   DM             dm;
614   Mat            J;
615   Vec            R;
616   const PetscInt nstages = irk->nstages;
617   PetscInt       vsize,bs;
618   PetscErrorCode ierr;
619 
620   PetscFunctionBegin;
621   if (!irk->work) {
622     ierr = PetscMalloc1(irk->nstages,&irk->work);CHKERRQ(ierr);
623   }
624   if (!irk->Y) {
625     ierr = VecDuplicateVecs(ts->vec_sol,irk->nstages,&irk->Y);CHKERRQ(ierr);
626   }
627   if (!irk->YdotI) {
628     ierr = VecDuplicateVecs(ts->vec_sol,irk->nstages,&irk->YdotI);CHKERRQ(ierr);
629   }
630   if (!irk->Ydot) {
631     ierr = VecDuplicate(ts->vec_sol,&irk->Ydot);CHKERRQ(ierr);
632   }
633   if (!irk->U) {
634     ierr = VecDuplicate(ts->vec_sol,&irk->U);CHKERRQ(ierr);
635   }
636   if (!irk->U0) {
637     ierr = VecDuplicate(ts->vec_sol,&irk->U0);CHKERRQ(ierr);
638   }
639   if (!irk->Z) {
640     ierr = VecCreate(PetscObjectComm((PetscObject)ts->vec_sol),&irk->Z);CHKERRQ(ierr);
641     ierr = VecGetSize(ts->vec_sol,&vsize);CHKERRQ(ierr);
642     ierr = VecSetSizes(irk->Z,PETSC_DECIDE,vsize*irk->nstages);CHKERRQ(ierr);
643     ierr = VecGetBlockSize(ts->vec_sol,&bs);CHKERRQ(ierr);
644     ierr = VecSetBlockSize(irk->Z,irk->nstages*bs);CHKERRQ(ierr);
645     ierr = VecSetFromOptions(irk->Z);CHKERRQ(ierr);
646   }
647   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
648   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSIRK,DMRestrictHook_TSIRK,ts);CHKERRQ(ierr);
649   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSIRK,DMSubDomainRestrictHook_TSIRK,ts);CHKERRQ(ierr);
650 
651   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
652   ierr = VecDuplicate(irk->Z,&R);CHKERRQ(ierr);
653   ierr = SNESSetFunction(ts->snes,R,SNESTSFormFunction,ts);CHKERRQ(ierr);
654   ierr = SNESGetJacobian(ts->snes,&J,NULL,NULL,NULL);CHKERRQ(ierr);
655   if (!irk->TJ) {
656     /* Create the KAIJ matrix for solving the stages */
657     ierr = MatCreateKAIJ(J,nstages,nstages,tab->A_inv,tab->I_s,&irk->TJ);CHKERRQ(ierr);
658   }
659   ierr = SNESSetJacobian(ts->snes,irk->TJ,irk->TJ,SNESTSFormJacobian,ts);CHKERRQ(ierr);
660   ierr = VecDestroy(&R);CHKERRQ(ierr);
661   PetscFunctionReturn(0);
662 }
663 
664 static PetscErrorCode TSSetFromOptions_IRK(PetscOptionItems *PetscOptionsObject,TS ts)
665 {
666   TS_IRK         *irk = (TS_IRK*)ts->data;
667   char           tname[256] = TSIRKGAUSS;
668   PetscErrorCode ierr;
669 
670   PetscFunctionBegin;
671   ierr = PetscOptionsHead(PetscOptionsObject,"IRK ODE solver options");CHKERRQ(ierr);
672   {
673     PetscBool flg1,flg2;
674     ierr = PetscOptionsInt("-ts_irk_nstages","Stages of the IRK method","TSIRKSetNumStages",irk->nstages,&irk->nstages,&flg1);CHKERRQ(ierr);
675     ierr = PetscOptionsFList("-ts_irk_type","Type of IRK method","TSIRKSetType",TSIRKList,irk->method_name[0] ? irk->method_name : tname,tname,sizeof(tname),&flg2);CHKERRQ(ierr);
676     if (flg1 ||flg2 || !irk->method_name[0]) { /* Create the method tableau after nstages or method is set */
677       ierr = TSIRKSetType(ts,tname);CHKERRQ(ierr);
678     }
679   }
680   ierr = PetscOptionsTail();CHKERRQ(ierr);
681   PetscFunctionReturn(0);
682 }
683 
684 static PetscErrorCode TSView_IRK(TS ts,PetscViewer viewer)
685 {
686   TS_IRK         *irk = (TS_IRK*)ts->data;
687   PetscBool      iascii;
688   PetscErrorCode ierr;
689 
690   PetscFunctionBegin;
691   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
692   if (iascii) {
693     IRKTableau    tab = irk->tableau;
694     TSIRKType irktype;
695     char          buf[512];
696 
697     ierr = TSIRKGetType(ts,&irktype);CHKERRQ(ierr);
698     ierr = PetscViewerASCIIPrintf(viewer,"  IRK type %s\n",irktype);CHKERRQ(ierr);
699     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",irk->nstages,tab->c);CHKERRQ(ierr);
700     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa       c = %s\n",buf);CHKERRQ(ierr);
701     ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",irk->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr);
702     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",PetscSqr(irk->nstages),tab->A);CHKERRQ(ierr);
703     ierr = PetscViewerASCIIPrintf(viewer,"  A coefficients       A = %s\n",buf);CHKERRQ(ierr);
704   }
705   PetscFunctionReturn(0);
706 }
707 
708 static PetscErrorCode TSLoad_IRK(TS ts,PetscViewer viewer)
709 {
710   PetscErrorCode ierr;
711   SNES           snes;
712   TSAdapt        adapt;
713 
714   PetscFunctionBegin;
715   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
716   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
717   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
718   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
719   /* function and Jacobian context for SNES when used with TS is always ts object */
720   ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr);
721   ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr);
722   PetscFunctionReturn(0);
723 }
724 
725 /*@C
726   TSIRKSetType - Set the type of IRK scheme
727 
728   Logically collective
729 
730   Input Parameters:
731 +  ts - timestepping context
732 -  irktype - type of IRK scheme
733 
734   Options Database:
735 .  -ts_irk_type <gauss>
736 
737   Level: intermediate
738 
739 .seealso: TSIRKGetType(), TSIRK, TSIRKType, TSIRKGAUSS
740 @*/
741 PetscErrorCode TSIRKSetType(TS ts,TSIRKType irktype)
742 {
743   PetscErrorCode ierr;
744 
745   PetscFunctionBegin;
746   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
747   PetscValidCharPointer(irktype,2);
748   ierr = PetscTryMethod(ts,"TSIRKSetType_C",(TS,TSIRKType),(ts,irktype));CHKERRQ(ierr);
749   PetscFunctionReturn(0);
750 }
751 
752 /*@C
753   TSIRKGetType - Get the type of IRK IMEX scheme
754 
755   Logically collective
756 
757   Input Parameter:
758 .  ts - timestepping context
759 
760   Output Parameter:
761 .  irktype - type of IRK-IMEX scheme
762 
763   Level: intermediate
764 
765 .seealso: TSIRKGetType()
766 @*/
767 PetscErrorCode TSIRKGetType(TS ts,TSIRKType *irktype)
768 {
769   PetscErrorCode ierr;
770 
771   PetscFunctionBegin;
772   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
773   ierr = PetscUseMethod(ts,"TSIRKGetType_C",(TS,TSIRKType*),(ts,irktype));CHKERRQ(ierr);
774   PetscFunctionReturn(0);
775 }
776 
777 /*@C
778   TSIRKSetNumStages - Set the number of stages of IRK scheme
779 
780   Logically collective
781 
782   Input Parameters:
783 +  ts - timestepping context
784 -  nstages - number of stages of IRK scheme
785 
786   Options Database:
787 .  -ts_irk_nstages <int>: Number of stages
788 
789   Level: intermediate
790 
791 .seealso: TSIRKGetNumStages(), TSIRK
792 @*/
793 PetscErrorCode TSIRKSetNumStages(TS ts,PetscInt nstages)
794 {
795   PetscErrorCode ierr;
796 
797   PetscFunctionBegin;
798   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
799   ierr = PetscTryMethod(ts,"TSIRKSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr);
800   PetscFunctionReturn(0);
801 }
802 
803 /*@C
804   TSIRKGetNumStages - Get the number of stages of IRK scheme
805 
806   Logically collective
807 
808   Input Parameters:
809 +  ts - timestepping context
810 -  nstages - number of stages of IRK scheme
811 
812   Level: intermediate
813 
814 .seealso: TSIRKSetNumStages(), TSIRK
815 @*/
816 PetscErrorCode TSIRKGetNumStages(TS ts,PetscInt *nstages)
817 {
818   PetscErrorCode ierr;
819 
820   PetscFunctionBegin;
821   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
822   PetscValidIntPointer(nstages,2);
823   ierr = PetscTryMethod(ts,"TSIRKGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr);
824   PetscFunctionReturn(0);
825 }
826 
827 static PetscErrorCode TSIRKGetType_IRK(TS ts,TSIRKType *irktype)
828 {
829   TS_IRK *irk = (TS_IRK*)ts->data;
830 
831   PetscFunctionBegin;
832   *irktype = irk->method_name;
833   PetscFunctionReturn(0);
834 }
835 
836 static PetscErrorCode TSIRKSetType_IRK(TS ts,TSIRKType irktype)
837 {
838   TS_IRK         *irk = (TS_IRK*)ts->data;
839   PetscErrorCode ierr,(*irkcreate)(TS);
840 
841   PetscFunctionBegin;
842   if (irk->method_name) {
843     ierr = PetscFree(irk->method_name);CHKERRQ(ierr);
844     ierr = TSIRKTableauReset(ts);CHKERRQ(ierr);
845   }
846   ierr = PetscFunctionListFind(TSIRKList,irktype,&irkcreate);CHKERRQ(ierr);
847   if (!irkcreate) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSIRK type \"%s\" given",irktype);
848   ierr = (*irkcreate)(ts);CHKERRQ(ierr);
849   ierr = PetscStrallocpy(irktype,&irk->method_name);CHKERRQ(ierr);
850   PetscFunctionReturn(0);
851 }
852 
853 static PetscErrorCode TSIRKSetNumStages_IRK(TS ts,PetscInt nstages)
854 {
855   TS_IRK *irk = (TS_IRK*)ts->data;
856 
857   PetscFunctionBegin;
858   if (nstages<=0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"input argument, %d, out of range",nstages);
859   irk->nstages = nstages;
860   PetscFunctionReturn(0);
861 }
862 
863 static PetscErrorCode TSIRKGetNumStages_IRK(TS ts,PetscInt *nstages)
864 {
865   TS_IRK *irk = (TS_IRK*)ts->data;
866 
867   PetscFunctionBegin;
868   PetscValidIntPointer(nstages,2);
869   *nstages = irk->nstages;
870   PetscFunctionReturn(0);
871 }
872 
873 static PetscErrorCode TSDestroy_IRK(TS ts)
874 {
875   PetscErrorCode ierr;
876 
877   PetscFunctionBegin;
878   ierr = TSReset_IRK(ts);CHKERRQ(ierr);
879   if (ts->dm) {
880     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSIRK,DMRestrictHook_TSIRK,ts);CHKERRQ(ierr);
881     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSIRK,DMSubDomainRestrictHook_TSIRK,ts);CHKERRQ(ierr);
882   }
883   ierr = PetscFree(ts->data);CHKERRQ(ierr);
884   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetType_C",NULL);CHKERRQ(ierr);
885   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetType_C",NULL);CHKERRQ(ierr);
886   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetNumStages_C",NULL);CHKERRQ(ierr);
887   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetNumStages_C",NULL);CHKERRQ(ierr);
888   PetscFunctionReturn(0);
889 }
890 
891 /*MC
892       TSIRK - ODE and DAE solver using Implicit Runge-Kutta schemes
893 
894   Notes:
895 
896   TSIRK uses the sparse Kronecker product matrix implementation of MATKAIJ to achieve good arithmetic intensity.
897 
898   The default is TSIRK3, it can be changed with TSIRKSetType() or -ts_irk_type
899 
900   If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information.
901 
902   Consider trying TSROSW if the stiff part is linear or weakly nonlinear.
903 
904   Level: beginner
905 
906 .seealso:  TSCreate(), TS, TSSetType(), TSIRKSetType(), TSIRKGetType(),
907            TSIRK1BEE, TSIRK2C, TSIRK2D, TSIRK2E, TSIRK3, TSIRKL2, TSIRKA2, TSIRKARS122,
908            TSIRK4, TSIRK5, TSIRKPRSSP2, TSIRKARS443, TSIRKBPR3, TSIRKType, TSIRKRegister()
909 
910 M*/
911 PETSC_EXTERN PetscErrorCode TSCreate_IRK(TS ts)
912 {
913   TS_IRK         *irk;
914   PetscErrorCode ierr;
915 
916   PetscFunctionBegin;
917   ierr = TSIRKInitializePackage();CHKERRQ(ierr);
918 
919   ts->ops->reset          = TSReset_IRK;
920   ts->ops->destroy        = TSDestroy_IRK;
921   ts->ops->view           = TSView_IRK;
922   ts->ops->load           = TSLoad_IRK;
923   ts->ops->setup          = TSSetUp_IRK;
924   ts->ops->step           = TSStep_IRK;
925   ts->ops->interpolate    = TSInterpolate_IRK;
926   ts->ops->evaluatestep   = TSEvaluateStep_IRK;
927   ts->ops->rollback       = TSRollBack_IRK;
928   ts->ops->setfromoptions = TSSetFromOptions_IRK;
929   ts->ops->snesfunction   = SNESTSFormFunction_IRK;
930   ts->ops->snesjacobian   = SNESTSFormJacobian_IRK;
931 
932   ts->usessnes = PETSC_TRUE;
933 
934   ierr = PetscNewLog(ts,&irk);CHKERRQ(ierr);
935   ts->data = (void*)irk;
936 
937   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetType_C",TSIRKSetType_IRK);CHKERRQ(ierr);
938   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetType_C",TSIRKGetType_IRK);CHKERRQ(ierr);
939   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetNumStages_C",TSIRKSetNumStages_IRK);CHKERRQ(ierr);
940   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetNumStages_C",TSIRKGetNumStages_IRK);CHKERRQ(ierr);
941   /* 3-stage IRK_Gauss is the default */
942   ierr = PetscNew(&irk->tableau);CHKERRQ(ierr);
943   irk->nstages = 3;
944   ierr = TSIRKSetType(ts,TSIRKDefault);CHKERRQ(ierr);
945   PetscFunctionReturn(0);
946 }
947