xref: /petsc/src/ts/tests/ex29.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
14a658b32SHong Zhang static char help[] ="Tests TS time span \n\n";
24a658b32SHong Zhang 
34a658b32SHong Zhang #include <petscts.h>
44a658b32SHong Zhang 
54a658b32SHong Zhang static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
64a658b32SHong Zhang {
74a658b32SHong Zhang   PetscInt          i,n;
84a658b32SHong Zhang   const PetscScalar *xx;
9e1db57b0SHong Zhang   PetscScalar       *ff;
104a658b32SHong Zhang 
114a658b32SHong Zhang   PetscFunctionBegin;
124a658b32SHong Zhang   PetscCall(VecGetLocalSize(X,&n));
134a658b32SHong Zhang   PetscCall(VecGetArrayRead(X,&xx));
144a658b32SHong Zhang   PetscCall(VecGetArray(F,&ff));
154a658b32SHong Zhang   if (n >= 1) ff[0] = 1;
164a658b32SHong Zhang   for (i = 1; i < n; i++) ff[i] = (i+1)*(xx[i-1]+PetscPowReal(t,i))/2;
174a658b32SHong Zhang   PetscCall(VecRestoreArrayRead(X,&xx));
184a658b32SHong Zhang   PetscCall(VecRestoreArray(F,&ff));
194a658b32SHong Zhang   PetscFunctionReturn(0);
204a658b32SHong Zhang }
214a658b32SHong Zhang 
224a658b32SHong Zhang int main(int argc, char *argv[])
234a658b32SHong Zhang {
244a658b32SHong Zhang   TS              ts;
254a658b32SHong Zhang   Vec             X,*Xs;
264a658b32SHong Zhang   PetscInt        i,n,N = 9;
27e1db57b0SHong Zhang   PetscReal       tspan[8] = {16.0, 16.1, 16.2, 16.3, 16.4, 16.5, 16.6, 16.7};
284a658b32SHong Zhang   const PetscReal *tspan2;
294a658b32SHong Zhang 
30*327415f7SBarry Smith   PetscFunctionBeginUser;
314a658b32SHong Zhang   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
324a658b32SHong Zhang   PetscCall(TSCreate(PETSC_COMM_SELF,&ts));
334a658b32SHong Zhang   PetscCall(TSSetType(ts,TSRK));
344a658b32SHong Zhang   PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,NULL));
354a658b32SHong Zhang   PetscCall(VecCreateSeq(PETSC_COMM_SELF,N,&X));
364a658b32SHong Zhang   PetscCall(VecZeroEntries(X));
37e1db57b0SHong Zhang   PetscCall(TSSetTimeStep(ts,0.001));
38e1db57b0SHong Zhang   PetscCall(TSSetTimeSpan(ts,8,tspan));
394a658b32SHong Zhang   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
404a658b32SHong Zhang   PetscCall(TSSetFromOptions(ts));
414a658b32SHong Zhang   PetscCall(TSSolve(ts,X));
424a658b32SHong Zhang   PetscCall(TSGetTimeSpanSolutions(ts,&n,&Xs));
434a658b32SHong Zhang   PetscCall(TSGetTimeSpan(ts,&n,&tspan2));
444a658b32SHong Zhang   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Time Span: "));
454a658b32SHong Zhang   for (i=0; i<n; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %g",(double)tspan2[i]));
464a658b32SHong Zhang   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
474a658b32SHong Zhang   PetscCall(TSDestroy(&ts));
484a658b32SHong Zhang   PetscCall(VecDestroy(&X));
494a658b32SHong Zhang   PetscCall(PetscFinalize());
504a658b32SHong Zhang   return 0;
514a658b32SHong Zhang }
524a658b32SHong Zhang 
534a658b32SHong Zhang /*TEST
544a658b32SHong Zhang 
554a658b32SHong Zhang testset:
564a658b32SHong Zhang   test:
574a658b32SHong Zhang     suffix: 1
584a658b32SHong Zhang     args: -ts_monitor
594a658b32SHong Zhang   test:
604a658b32SHong Zhang     suffix: 2
61e1db57b0SHong Zhang     requires: !single
62e1db57b0SHong Zhang     args: -ts_monitor -ts_adapt_type none
634a658b32SHong Zhang TEST*/
64