xref: /petsc/src/ts/impls/explicit/ssp/ssp.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 /*
2        Code for Timestepping with explicit SSP.
3 */
4 #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
5 
6 PetscFunctionList TSSSPList = NULL;
7 static PetscBool  TSSSPPackageInitialized;
8 
9 typedef struct {
10   PetscErrorCode (*onestep)(TS, PetscReal, PetscReal, Vec);
11   char     *type_name;
12   PetscInt  nstages;
13   Vec      *work;
14   PetscInt  nwork;
15   PetscBool workout;
16 } TS_SSP;
17 
18 static PetscErrorCode TSSSPGetWorkVectors(TS ts, PetscInt n, Vec **work) {
19   TS_SSP *ssp = (TS_SSP *)ts->data;
20 
21   PetscFunctionBegin;
22   PetscCheck(!ssp->workout, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Work vectors already gotten");
23   if (ssp->nwork < n) {
24     if (ssp->nwork > 0) PetscCall(VecDestroyVecs(ssp->nwork, &ssp->work));
25     PetscCall(VecDuplicateVecs(ts->vec_sol, n, &ssp->work));
26     ssp->nwork = n;
27   }
28   *work        = ssp->work;
29   ssp->workout = PETSC_TRUE;
30   PetscFunctionReturn(0);
31 }
32 
33 static PetscErrorCode TSSSPRestoreWorkVectors(TS ts, PetscInt n, Vec **work) {
34   TS_SSP *ssp = (TS_SSP *)ts->data;
35 
36   PetscFunctionBegin;
37   PetscCheck(ssp->workout, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Work vectors have not been gotten");
38   PetscCheck(*work == ssp->work, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong work vectors checked out");
39   ssp->workout = PETSC_FALSE;
40   *work        = NULL;
41   PetscFunctionReturn(0);
42 }
43 
44 /*MC
45    TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s
46 
47    Pseudocode 2 of Ketcheson 2008
48 
49    Level: beginner
50 
51 .seealso: `TSSSP`, `TSSSPSetType()`, `TSSSPSetNumStages()`
52 M*/
53 static PetscErrorCode TSSSPStep_RK_2(TS ts, PetscReal t0, PetscReal dt, Vec sol) {
54   TS_SSP  *ssp = (TS_SSP *)ts->data;
55   Vec     *work, F;
56   PetscInt i, s;
57 
58   PetscFunctionBegin;
59   s = ssp->nstages;
60   PetscCall(TSSSPGetWorkVectors(ts, 2, &work));
61   F = work[1];
62   PetscCall(VecCopy(sol, work[0]));
63   for (i = 0; i < s - 1; i++) {
64     PetscReal stage_time = t0 + dt * (i / (s - 1.));
65     PetscCall(TSPreStage(ts, stage_time));
66     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
67     PetscCall(VecAXPY(work[0], dt / (s - 1.), F));
68   }
69   PetscCall(TSComputeRHSFunction(ts, t0 + dt, work[0], F));
70   PetscCall(VecAXPBYPCZ(sol, (s - 1.) / s, dt / s, 1. / s, work[0], F));
71   PetscCall(TSSSPRestoreWorkVectors(ts, 2, &work));
72   PetscFunctionReturn(0);
73 }
74 
75 /*MC
76    TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(PetscSqrtReal(s)-1)/PetscSqrtReal(s), where PetscSqrtReal(s) is an integer
77 
78    Pseudocode 2 of Ketcheson 2008
79 
80    Level: beginner
81 
82 .seealso: `TSSSP`, `TSSSPSetType()`, `TSSSPSetNumStages()`
83 M*/
84 static PetscErrorCode TSSSPStep_RK_3(TS ts, PetscReal t0, PetscReal dt, Vec sol) {
85   TS_SSP   *ssp = (TS_SSP *)ts->data;
86   Vec      *work, F;
87   PetscInt  i, s, n, r;
88   PetscReal c, stage_time;
89 
90   PetscFunctionBegin;
91   s = ssp->nstages;
92   n = (PetscInt)(PetscSqrtReal((PetscReal)s) + 0.001);
93   r = s - n;
94   PetscCheck(n * n == s, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for optimal third order schemes with %" PetscInt_FMT " stages, must be a square number at least 4", s);
95   PetscCall(TSSSPGetWorkVectors(ts, 3, &work));
96   F = work[2];
97   PetscCall(VecCopy(sol, work[0]));
98   for (i = 0; i < (n - 1) * (n - 2) / 2; i++) {
99     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
100     stage_time = t0 + c * dt;
101     PetscCall(TSPreStage(ts, stage_time));
102     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
103     PetscCall(VecAXPY(work[0], dt / r, F));
104   }
105   PetscCall(VecCopy(work[0], work[1]));
106   for (; i < n * (n + 1) / 2 - 1; i++) {
107     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
108     stage_time = t0 + c * dt;
109     PetscCall(TSPreStage(ts, stage_time));
110     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
111     PetscCall(VecAXPY(work[0], dt / r, F));
112   }
113   {
114     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
115     stage_time = t0 + c * dt;
116     PetscCall(TSPreStage(ts, stage_time));
117     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
118     PetscCall(VecAXPBYPCZ(work[0], 1. * n / (2 * n - 1.), (n - 1.) * dt / (r * (2 * n - 1)), (n - 1.) / (2 * n - 1.), work[1], F));
119     i++;
120   }
121   for (; i < s; i++) {
122     c          = (i < n * (n + 1) / 2) ? 1. * i / (s - n) : (1. * i - n) / (s - n);
123     stage_time = t0 + c * dt;
124     PetscCall(TSPreStage(ts, stage_time));
125     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
126     PetscCall(VecAXPY(work[0], dt / r, F));
127   }
128   PetscCall(VecCopy(work[0], sol));
129   PetscCall(TSSSPRestoreWorkVectors(ts, 3, &work));
130   PetscFunctionReturn(0);
131 }
132 
133 /*MC
134    TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6
135 
136    SSPRK(10,4), Pseudocode 3 of Ketcheson 2008
137 
138    Level: beginner
139 
140 .seealso: `TSSSP`, `TSSSPSetType()`
141 M*/
142 static PetscErrorCode TSSSPStep_RK_10_4(TS ts, PetscReal t0, PetscReal dt, Vec sol) {
143   const PetscReal c[10] = {0, 1. / 6, 2. / 6, 3. / 6, 4. / 6, 2. / 6, 3. / 6, 4. / 6, 5. / 6, 1};
144   Vec            *work, F;
145   PetscInt        i;
146   PetscReal       stage_time;
147 
148   PetscFunctionBegin;
149   PetscCall(TSSSPGetWorkVectors(ts, 3, &work));
150   F = work[2];
151   PetscCall(VecCopy(sol, work[0]));
152   for (i = 0; i < 5; i++) {
153     stage_time = t0 + c[i] * dt;
154     PetscCall(TSPreStage(ts, stage_time));
155     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
156     PetscCall(VecAXPY(work[0], dt / 6, F));
157   }
158   PetscCall(VecAXPBYPCZ(work[1], 1. / 25, 9. / 25, 0, sol, work[0]));
159   PetscCall(VecAXPBY(work[0], 15, -5, work[1]));
160   for (; i < 9; i++) {
161     stage_time = t0 + c[i] * dt;
162     PetscCall(TSPreStage(ts, stage_time));
163     PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
164     PetscCall(VecAXPY(work[0], dt / 6, F));
165   }
166   stage_time = t0 + dt;
167   PetscCall(TSPreStage(ts, stage_time));
168   PetscCall(TSComputeRHSFunction(ts, stage_time, work[0], F));
169   PetscCall(VecAXPBYPCZ(work[1], 3. / 5, dt / 10, 1, work[0], F));
170   PetscCall(VecCopy(work[1], sol));
171   PetscCall(TSSSPRestoreWorkVectors(ts, 3, &work));
172   PetscFunctionReturn(0);
173 }
174 
175 static PetscErrorCode TSSetUp_SSP(TS ts) {
176   PetscFunctionBegin;
177   PetscCall(TSCheckImplicitTerm(ts));
178   PetscCall(TSGetAdapt(ts, &ts->adapt));
179   PetscCall(TSAdaptCandidatesClear(ts->adapt));
180   PetscFunctionReturn(0);
181 }
182 
183 static PetscErrorCode TSStep_SSP(TS ts) {
184   TS_SSP   *ssp = (TS_SSP *)ts->data;
185   Vec       sol = ts->vec_sol;
186   PetscBool stageok, accept = PETSC_TRUE;
187   PetscReal next_time_step = ts->time_step;
188 
189   PetscFunctionBegin;
190   PetscCall((*ssp->onestep)(ts, ts->ptime, ts->time_step, sol));
191   PetscCall(TSPostStage(ts, ts->ptime, 0, &sol));
192   PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, sol, &stageok));
193   if (!stageok) {
194     ts->reason = TS_DIVERGED_STEP_REJECTED;
195     PetscFunctionReturn(0);
196   }
197 
198   PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
199   if (!accept) {
200     ts->reason = TS_DIVERGED_STEP_REJECTED;
201     PetscFunctionReturn(0);
202   }
203 
204   ts->ptime += ts->time_step;
205   ts->time_step = next_time_step;
206   PetscFunctionReturn(0);
207 }
208 /*------------------------------------------------------------*/
209 
210 static PetscErrorCode TSReset_SSP(TS ts) {
211   TS_SSP *ssp = (TS_SSP *)ts->data;
212 
213   PetscFunctionBegin;
214   if (ssp->work) PetscCall(VecDestroyVecs(ssp->nwork, &ssp->work));
215   ssp->nwork   = 0;
216   ssp->workout = PETSC_FALSE;
217   PetscFunctionReturn(0);
218 }
219 
220 static PetscErrorCode TSDestroy_SSP(TS ts) {
221   TS_SSP *ssp = (TS_SSP *)ts->data;
222 
223   PetscFunctionBegin;
224   PetscCall(TSReset_SSP(ts));
225   PetscCall(PetscFree(ssp->type_name));
226   PetscCall(PetscFree(ts->data));
227   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetType_C", NULL));
228   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetType_C", NULL));
229   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetNumStages_C", NULL));
230   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetNumStages_C", NULL));
231   PetscFunctionReturn(0);
232 }
233 /*------------------------------------------------------------*/
234 
235 /*@C
236    TSSSPSetType - set the SSP time integration scheme to use
237 
238    Logically Collective
239 
240    Input Parameters:
241 +  ts - time stepping object
242 -  ssptype - type of scheme to use
243 
244    Options Database Keys:
245    -ts_ssp_type <rks2>               : Type of SSP method (one of) rks2 rks3 rk104
246    -ts_ssp_nstages<rks2: 5, rks3: 9> : Number of stages
247 
248    Level: beginner
249 
250 .seealso: `TSSSP`, `TSSSPGetType()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
251 @*/
252 PetscErrorCode TSSSPSetType(TS ts, TSSSPType ssptype) {
253   PetscFunctionBegin;
254   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
255   PetscValidCharPointer(ssptype, 2);
256   PetscTryMethod(ts, "TSSSPSetType_C", (TS, TSSSPType), (ts, ssptype));
257   PetscFunctionReturn(0);
258 }
259 
260 /*@C
261    TSSSPGetType - get the SSP time integration scheme
262 
263    Logically Collective
264 
265    Input Parameter:
266 .  ts - time stepping object
267 
268    Output Parameter:
269 .  type - type of scheme being used
270 
271    Level: beginner
272 
273 .seealso: `TSSSP`, `TSSSPSettype()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
274 @*/
275 PetscErrorCode TSSSPGetType(TS ts, TSSSPType *type) {
276   PetscFunctionBegin;
277   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
278   PetscUseMethod(ts, "TSSSPGetType_C", (TS, TSSSPType *), (ts, type));
279   PetscFunctionReturn(0);
280 }
281 
282 /*@
283    TSSSPSetNumStages - set the number of stages to use with the SSP method. Must be called after
284    TSSSPSetType().
285 
286    Logically Collective
287 
288    Input Parameters:
289 +  ts - time stepping object
290 -  nstages - number of stages
291 
292    Options Database Keys:
293    -ts_ssp_type <rks2>               : Type of SSP method (one of) rks2 rks3 rk104
294    -ts_ssp_nstages<rks2: 5, rks3: 9> : Number of stages
295 
296    Level: beginner
297 
298 .seealso: `TSSSP`, `TSSSPGetNumStages()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
299 @*/
300 PetscErrorCode TSSSPSetNumStages(TS ts, PetscInt nstages) {
301   PetscFunctionBegin;
302   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
303   PetscTryMethod(ts, "TSSSPSetNumStages_C", (TS, PetscInt), (ts, nstages));
304   PetscFunctionReturn(0);
305 }
306 
307 /*@
308    TSSSPGetNumStages - get the number of stages in the SSP time integration scheme
309 
310    Logically Collective
311 
312    Input Parameter:
313 .  ts - time stepping object
314 
315    Output Parameter:
316 .  nstages - number of stages
317 
318    Level: beginner
319 
320 .seealso: `TSSSP`, `TSSSPGetType()`, `TSSSPSetNumStages()`, `TSSSPRKS2`, `TSSSPRKS3`, `TSSSPRK104`
321 @*/
322 PetscErrorCode TSSSPGetNumStages(TS ts, PetscInt *nstages) {
323   PetscFunctionBegin;
324   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
325   PetscUseMethod(ts, "TSSSPGetNumStages_C", (TS, PetscInt *), (ts, nstages));
326   PetscFunctionReturn(0);
327 }
328 
329 static PetscErrorCode TSSSPSetType_SSP(TS ts, TSSSPType type) {
330   TS_SSP *ssp = (TS_SSP *)ts->data;
331   PetscErrorCode (*r)(TS, PetscReal, PetscReal, Vec);
332   PetscBool flag;
333 
334   PetscFunctionBegin;
335   PetscCall(PetscFunctionListFind(TSSSPList, type, &r));
336   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TS_SSP type %s given", type);
337   ssp->onestep = r;
338   PetscCall(PetscFree(ssp->type_name));
339   PetscCall(PetscStrallocpy(type, &ssp->type_name));
340   ts->default_adapt_type = TSADAPTNONE;
341   PetscCall(PetscStrcmp(type, TSSSPRKS2, &flag));
342   if (flag) {
343     ssp->nstages = 5;
344   } else {
345     PetscCall(PetscStrcmp(type, TSSSPRKS3, &flag));
346     if (flag) {
347       ssp->nstages = 9;
348     } else {
349       ssp->nstages = 5;
350     }
351   }
352   PetscFunctionReturn(0);
353 }
354 static PetscErrorCode TSSSPGetType_SSP(TS ts, TSSSPType *type) {
355   TS_SSP *ssp = (TS_SSP *)ts->data;
356 
357   PetscFunctionBegin;
358   *type = ssp->type_name;
359   PetscFunctionReturn(0);
360 }
361 static PetscErrorCode TSSSPSetNumStages_SSP(TS ts, PetscInt nstages) {
362   TS_SSP *ssp = (TS_SSP *)ts->data;
363 
364   PetscFunctionBegin;
365   ssp->nstages = nstages;
366   PetscFunctionReturn(0);
367 }
368 static PetscErrorCode TSSSPGetNumStages_SSP(TS ts, PetscInt *nstages) {
369   TS_SSP *ssp = (TS_SSP *)ts->data;
370 
371   PetscFunctionBegin;
372   *nstages = ssp->nstages;
373   PetscFunctionReturn(0);
374 }
375 
376 static PetscErrorCode TSSetFromOptions_SSP(TS ts, PetscOptionItems *PetscOptionsObject) {
377   char      tname[256] = TSSSPRKS2;
378   TS_SSP   *ssp        = (TS_SSP *)ts->data;
379   PetscBool flg;
380 
381   PetscFunctionBegin;
382   PetscOptionsHeadBegin(PetscOptionsObject, "SSP ODE solver options");
383   {
384     PetscCall(PetscOptionsFList("-ts_ssp_type", "Type of SSP method", "TSSSPSetType", TSSSPList, tname, tname, sizeof(tname), &flg));
385     if (flg) PetscCall(TSSSPSetType(ts, tname));
386     PetscCall(PetscOptionsInt("-ts_ssp_nstages", "Number of stages", "TSSSPSetNumStages", ssp->nstages, &ssp->nstages, NULL));
387   }
388   PetscOptionsHeadEnd();
389   PetscFunctionReturn(0);
390 }
391 
392 static PetscErrorCode TSView_SSP(TS ts, PetscViewer viewer) {
393   TS_SSP   *ssp = (TS_SSP *)ts->data;
394   PetscBool ascii;
395 
396   PetscFunctionBegin;
397   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &ascii));
398   if (ascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  Scheme: %s\n", ssp->type_name));
399   PetscFunctionReturn(0);
400 }
401 
402 /* ------------------------------------------------------------ */
403 
404 /*MC
405       TSSSP - Explicit strong stability preserving ODE solver
406 
407   Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation
408   bounded (TVB) although these solutions often contain discontinuities.  Spatial discretizations such as Godunov's
409   scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties,
410   but they are usually formulated using a forward Euler time discretization or by coupling the space and time
411   discretization as in the classical Lax-Wendroff scheme.  When the space and time discretization is coupled, it is very
412   difficult to produce schemes with high temporal accuracy while preserving TVD properties.  An alternative is the
413   semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a
414   time discretization that preserves the TVD property.  Such integrators are called strong stability preserving (SSP).
415 
416   Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while
417   still being SSP.  Some theoretical bounds
418 
419   1. There are no explicit methods with c_eff > 1.
420 
421   2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0.
422 
423   3. There are no implicit methods with order greater than 1 and c_eff > 2.
424 
425   This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff.  More stages allows
426   for larger values of c_eff which improves efficiency.  These implementations are low-memory and only use 2 or 3 work
427   vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice.
428 
429   Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104}
430 
431   rks2: Second order methods with any number s>1 of stages.  c_eff = (s-1)/s
432 
433   rks3: Third order methods with s=n^2 stages, n>1.  c_eff = (s-n)/s
434 
435   rk104: A 10-stage fourth order method.  c_eff = 0.6
436 
437   Level: beginner
438 
439   References:
440 +  * - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008.
441 -  * - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009.
442 
443 .seealso: `TSCreate()`, `TS`, `TSSetType()`
444 
445 M*/
446 PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts) {
447   TS_SSP *ssp;
448 
449   PetscFunctionBegin;
450   PetscCall(TSSSPInitializePackage());
451 
452   ts->ops->setup          = TSSetUp_SSP;
453   ts->ops->step           = TSStep_SSP;
454   ts->ops->reset          = TSReset_SSP;
455   ts->ops->destroy        = TSDestroy_SSP;
456   ts->ops->setfromoptions = TSSetFromOptions_SSP;
457   ts->ops->view           = TSView_SSP;
458 
459   PetscCall(PetscNewLog(ts, &ssp));
460   ts->data = (void *)ssp;
461 
462   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetType_C", TSSSPGetType_SSP));
463   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetType_C", TSSSPSetType_SSP));
464   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPGetNumStages_C", TSSSPGetNumStages_SSP));
465   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSSPSetNumStages_C", TSSSPSetNumStages_SSP));
466 
467   PetscCall(TSSSPSetType(ts, TSSSPRKS2));
468   PetscFunctionReturn(0);
469 }
470 
471 /*@C
472   TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called
473   from TSInitializePackage().
474 
475   Level: developer
476 
477 .seealso: `PetscInitialize()`
478 @*/
479 PetscErrorCode TSSSPInitializePackage(void) {
480   PetscFunctionBegin;
481   if (TSSSPPackageInitialized) PetscFunctionReturn(0);
482   TSSSPPackageInitialized = PETSC_TRUE;
483   PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRKS2, TSSSPStep_RK_2));
484   PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRKS3, TSSSPStep_RK_3));
485   PetscCall(PetscFunctionListAdd(&TSSSPList, TSSSPRK104, TSSSPStep_RK_10_4));
486   PetscCall(PetscRegisterFinalize(TSSSPFinalizePackage));
487   PetscFunctionReturn(0);
488 }
489 
490 /*@C
491   TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is
492   called from PetscFinalize().
493 
494   Level: developer
495 
496 .seealso: `PetscFinalize()`
497 @*/
498 PetscErrorCode TSSSPFinalizePackage(void) {
499   PetscFunctionBegin;
500   TSSSPPackageInitialized = PETSC_FALSE;
501   PetscCall(PetscFunctionListDestroy(&TSSSPList));
502   PetscFunctionReturn(0);
503 }
504