xref: /petsc/src/ts/tests/ex29.c (revision 2e16c0ce58b3a4ec287cbc0a0807bfb0a0fa5ac9)
1 static char help[] ="Tests TS time span \n\n";
2 
3 #include <petscts.h>
4 
5 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
6 {
7   PetscInt          i,n;
8   const PetscScalar *xx;
9   PetscScalar       *ff;
10 
11   PetscFunctionBegin;
12   PetscCall(VecGetLocalSize(X,&n));
13   PetscCall(VecGetArrayRead(X,&xx));
14   PetscCall(VecGetArray(F,&ff));
15   if (n >= 1) ff[0] = 1;
16   for (i = 1; i < n; i++) ff[i] = (i+1)*(xx[i-1]+PetscPowReal(t,i))/2;
17   PetscCall(VecRestoreArrayRead(X,&xx));
18   PetscCall(VecRestoreArray(F,&ff));
19   PetscFunctionReturn(0);
20 }
21 
22 int main(int argc, char *argv[])
23 {
24   TS              ts;
25   Vec             X,*Xs;
26   PetscInt        i,n,N = 9;
27   PetscReal       tspan[8] = {16.0, 16.1, 16.2, 16.3, 16.4, 16.5, 16.6, 16.7};
28   const PetscReal *tspan2;
29 
30   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
31   PetscCall(TSCreate(PETSC_COMM_SELF,&ts));
32   PetscCall(TSSetType(ts,TSRK));
33   PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,NULL));
34   PetscCall(VecCreateSeq(PETSC_COMM_SELF,N,&X));
35   PetscCall(VecZeroEntries(X));
36   PetscCall(TSSetTimeStep(ts,0.001));
37   PetscCall(TSSetTimeSpan(ts,8,tspan));
38   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
39   PetscCall(TSSetFromOptions(ts));
40   PetscCall(TSSolve(ts,X));
41   PetscCall(TSGetTimeSpanSolutions(ts,&n,&Xs));
42   PetscCall(TSGetTimeSpan(ts,&n,&tspan2));
43   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Time Span: "));
44   for (i=0; i<n; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %g",(double)tspan2[i]));
45   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
46   PetscCall(TSDestroy(&ts));
47   PetscCall(VecDestroy(&X));
48   PetscCall(PetscFinalize());
49   return 0;
50 }
51 
52 /*TEST
53 
54 testset:
55   test:
56     suffix: 1
57     args: -ts_monitor
58   test:
59     suffix: 2
60     requires: !single
61     args: -ts_monitor -ts_adapt_type none
62 TEST*/
63