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