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