xref: /petsc/src/ts/tests/ex29.c (revision a336c15037c72f93cd561f5a5e11e93175f2efd9)
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, PetscCtx ctx)
6 {
7   PetscInt           i, n;
8   const PetscScalar *xx;
9   PetscScalar       *ff;
10 
11   PetscFunctionBeginUser;
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(PETSC_SUCCESS);
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, *sol_times;
29 
30   PetscFunctionBeginUser;
31   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
32   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
33   PetscCall(TSSetType(ts, TSRK));
34   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL));
35   PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &X));
36   PetscCall(VecZeroEntries(X));
37   PetscCall(TSSetTimeStep(ts, 0.001));
38   PetscCall(TSSetTimeSpan(ts, 8, tspan));
39   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
40   PetscCall(TSSetFromOptions(ts));
41   PetscCall(TSSolve(ts, X));
42   PetscCall(TSGetEvaluationSolutions(ts, &n, &sol_times, &Xs));
43   PetscCall(TSGetEvaluationTimes(ts, &n, &tspan2));
44   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time Span: "));
45   for (i = 0; i < n; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %g", (double)tspan2[i]));
46   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
47   for (PetscInt i = 0; i < n; i++) {
48     PetscCheck(PetscIsCloseAtTol(tspan2[i], sol_times[i], 1e-6, 1e2 * PETSC_MACHINE_EPSILON), PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Requested solution at time %g, but received time at %g", (double)tspan[i], (double)sol_times[i]);
49   }
50   PetscCall(TSDestroy(&ts));
51   PetscCall(VecDestroy(&X));
52   PetscCall(PetscFinalize());
53   return 0;
54 }
55 
56 /*TEST
57 
58 testset:
59   test:
60     suffix: 1
61     args: -ts_monitor
62   test:
63     suffix: 2
64     requires: !single
65     args: -ts_monitor -ts_adapt_type none
66 TEST*/
67