xref: /petsc/src/ts/tests/ex2.c (revision eaf6e66d5871d10bd4e10f2c5f553a50279d20f9)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown        Formatted test for TS routines.
3c4762a1bSJed Brown 
4c4762a1bSJed Brown           Solves U_t=F(t,u)
5c4762a1bSJed Brown           Where:
6c4762a1bSJed Brown 
7c4762a1bSJed Brown                   [2*u1+u2
8c4762a1bSJed Brown           F(t,u)= [u1+2*u2+u3
9c4762a1bSJed Brown                   [   u2+2*u3
10c4762a1bSJed Brown        We can compare the solutions from euler, beuler and SUNDIALS to
11c4762a1bSJed Brown        see what is the difference.
12c4762a1bSJed Brown 
13c4762a1bSJed Brown */
14c4762a1bSJed Brown 
15c4762a1bSJed Brown static char help[] = "Solves a linear ODE. \n\n";
16c4762a1bSJed Brown 
17c4762a1bSJed Brown #include <petscts.h>
18c4762a1bSJed Brown #include <petscpc.h>
19c4762a1bSJed Brown 
20c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
21c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
22c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
23c4762a1bSJed Brown extern PetscErrorCode Initial(Vec, void *);
24c4762a1bSJed Brown extern PetscErrorCode MyMatMult(Mat, Vec, Vec);
25c4762a1bSJed Brown 
26c4762a1bSJed Brown extern PetscReal solx(PetscReal);
27c4762a1bSJed Brown extern PetscReal soly(PetscReal);
28c4762a1bSJed Brown extern PetscReal solz(PetscReal);
29c4762a1bSJed Brown 
30d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
31d71ae5a4SJacob Faibussowitsch {
32c4762a1bSJed Brown   PetscInt  time_steps = 100, steps;
33c4762a1bSJed Brown   Vec       global;
34c4762a1bSJed Brown   PetscReal dt, ftime;
35c4762a1bSJed Brown   TS        ts;
36c4762a1bSJed Brown   Mat       A = 0, S;
37*eaf6e66dSStefano Zampini   PetscBool nest;
38c4762a1bSJed Brown 
39327415f7SBarry Smith   PetscFunctionBeginUser;
409566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-time", &time_steps, NULL));
42*eaf6e66dSStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-nest", &nest, NULL));
43c4762a1bSJed Brown 
44*eaf6e66dSStefano Zampini   /* create vector to hold state */
45*eaf6e66dSStefano Zampini   if (nest) {
46*eaf6e66dSStefano Zampini     Vec g[3];
47*eaf6e66dSStefano Zampini 
48*eaf6e66dSStefano Zampini     PetscCall(VecCreate(PETSC_COMM_WORLD, &g[0]));
49*eaf6e66dSStefano Zampini     PetscCall(VecSetSizes(g[0], PETSC_DECIDE, 1));
50*eaf6e66dSStefano Zampini     PetscCall(VecSetFromOptions(g[0]));
51*eaf6e66dSStefano Zampini     PetscCall(VecDuplicate(g[0], &g[1]));
52*eaf6e66dSStefano Zampini     PetscCall(VecDuplicate(g[0], &g[2]));
53*eaf6e66dSStefano Zampini     PetscCall(VecCreateNest(PETSC_COMM_WORLD, 3, NULL, g, &global));
54*eaf6e66dSStefano Zampini     PetscCall(VecDestroy(&g[0]));
55*eaf6e66dSStefano Zampini     PetscCall(VecDestroy(&g[1]));
56*eaf6e66dSStefano Zampini     PetscCall(VecDestroy(&g[2]));
57*eaf6e66dSStefano Zampini   } else {
589566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD, &global));
599566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(global, PETSC_DECIDE, 3));
609566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(global));
61*eaf6e66dSStefano Zampini   }
62*eaf6e66dSStefano Zampini 
63*eaf6e66dSStefano Zampini   /* set initial conditions */
649566063dSJacob Faibussowitsch   PetscCall(Initial(global, NULL));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* make timestep context */
679566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
689566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
699566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, Monitor, NULL, NULL));
70c4762a1bSJed Brown   dt = 0.001;
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /*
73c4762a1bSJed Brown     The user provides the RHS and Jacobian
74c4762a1bSJed Brown   */
759566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL));
769566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
779566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 3, 3));
789566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
799566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
809566063dSJacob Faibussowitsch   PetscCall(RHSJacobian(ts, 0.0, global, A, A, NULL));
819566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL));
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, 3, 3, 3, 3, NULL, &S));
849566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S, MATOP_MULT, (void (*)(void))MyMatMult));
859566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, S, A, RHSJacobian, NULL));
86c4762a1bSJed Brown 
879566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
889566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
919566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps));
929566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 1));
939566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, global));
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, global));
969566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
979566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* free the memories */
100c4762a1bSJed Brown 
1019566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
1039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&S));
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
107b122ec5aSJacob Faibussowitsch   return 0;
108c4762a1bSJed Brown }
109c4762a1bSJed Brown 
110d71ae5a4SJacob Faibussowitsch PetscErrorCode MyMatMult(Mat S, Vec x, Vec y)
111d71ae5a4SJacob Faibussowitsch {
112c4762a1bSJed Brown   const PetscScalar *inptr;
113c4762a1bSJed Brown   PetscScalar       *outptr;
114c4762a1bSJed Brown 
1157510d9b0SBarry Smith   PetscFunctionBeginUser;
1169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &inptr));
1179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y, &outptr));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   outptr[0] = 2.0 * inptr[0] + inptr[1];
120c4762a1bSJed Brown   outptr[1] = inptr[0] + 2.0 * inptr[1] + inptr[2];
121c4762a1bSJed Brown   outptr[2] = inptr[1] + 2.0 * inptr[2];
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &inptr));
1249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y, &outptr));
1253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126c4762a1bSJed Brown }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown /* this test problem has initial values (1,1,1).                      */
129d71ae5a4SJacob Faibussowitsch PetscErrorCode Initial(Vec global, void *ctx)
130d71ae5a4SJacob Faibussowitsch {
131c4762a1bSJed Brown   PetscScalar *localptr;
132c4762a1bSJed Brown   PetscInt     i, mybase, myend, locsize;
133c4762a1bSJed Brown 
1343ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
135c4762a1bSJed Brown   /* determine starting point of each processor */
1369566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(global, &mybase, &myend));
1379566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global, &locsize));
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* Initialize the array */
1409566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(global, &localptr));
141c4762a1bSJed Brown   for (i = 0; i < locsize; i++) localptr[i] = 1.0;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   if (mybase == 0) localptr[0] = 1.0;
144c4762a1bSJed Brown 
1459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(global, &localptr));
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
147c4762a1bSJed Brown }
148c4762a1bSJed Brown 
149d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec global, void *ctx)
150d71ae5a4SJacob Faibussowitsch {
151c4762a1bSJed Brown   VecScatter         scatter;
152c4762a1bSJed Brown   IS                 from, to;
153c4762a1bSJed Brown   PetscInt           i, n, *idx;
154c4762a1bSJed Brown   Vec                tmp_vec;
155c4762a1bSJed Brown   const PetscScalar *tmp;
156c4762a1bSJed Brown 
1573ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
158c4762a1bSJed Brown   /* Get the size of the vector */
1599566063dSJacob Faibussowitsch   PetscCall(VecGetSize(global, &n));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* Set the index sets */
1629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idx));
163c4762a1bSJed Brown   for (i = 0; i < n; i++) idx[i] = i;
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /* Create local sequential vectors */
1669566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &tmp_vec));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* Create scatter context */
1699566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &from));
1709566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &to));
1719566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global, from, tmp_vec, to, &scatter));
1729566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD));
1739566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD));
174c4762a1bSJed Brown 
1759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_vec, &tmp));
1769371c9d4SSatish Balay   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "At t =%14.6e u = %14.6e  %14.6e  %14.6e \n", (double)time, (double)PetscRealPart(tmp[0]), (double)PetscRealPart(tmp[1]), (double)PetscRealPart(tmp[2])));
1779371c9d4SSatish Balay   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "At t =%14.6e errors = %14.6e  %14.6e  %14.6e \n", (double)time, (double)PetscRealPart(tmp[0] - solx(time)), (double)PetscRealPart(tmp[1] - soly(time)), (double)PetscRealPart(tmp[2] - solz(time))));
1789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_vec, &tmp));
1799566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
1809566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
1819566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
1829566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
1839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_vec));
1843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx)
188d71ae5a4SJacob Faibussowitsch {
189c4762a1bSJed Brown   PetscScalar       *outptr;
190c4762a1bSJed Brown   const PetscScalar *inptr;
191c4762a1bSJed Brown   PetscInt           i, n, *idx;
192c4762a1bSJed Brown   IS                 from, to;
193c4762a1bSJed Brown   VecScatter         scatter;
194c4762a1bSJed Brown   Vec                tmp_in, tmp_out;
195c4762a1bSJed Brown 
1963ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
197c4762a1bSJed Brown   /* Get the length of parallel vector */
1989566063dSJacob Faibussowitsch   PetscCall(VecGetSize(globalin, &n));
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* Set the index sets */
2019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idx));
202c4762a1bSJed Brown   for (i = 0; i < n; i++) idx[i] = i;
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* Create local sequential vectors */
2059566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &tmp_in));
2069566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tmp_in, &tmp_out));
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   /* Create scatter context */
2099566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &from));
2109566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &to));
2119566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(globalin, from, tmp_in, to, &scatter));
2129566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD));
2139566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD));
2149566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   /*Extract income array */
2179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_in, &inptr));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* Extract outcome array*/
2209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(tmp_out, &outptr));
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   outptr[0] = 2.0 * inptr[0] + inptr[1];
223c4762a1bSJed Brown   outptr[1] = inptr[0] + 2.0 * inptr[1] + inptr[2];
224c4762a1bSJed Brown   outptr[2] = inptr[1] + 2.0 * inptr[2];
225c4762a1bSJed Brown 
2269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_in, &inptr));
2279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(tmp_out, &outptr));
228c4762a1bSJed Brown 
2299566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(tmp_out, from, globalout, to, &scatter));
2309566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD));
2319566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* Destroy idx aand scatter */
2349566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
2359566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
2369566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
2379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_in));
2389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_out));
2399566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
2403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
241c4762a1bSJed Brown }
242c4762a1bSJed Brown 
243d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat BB, void *ctx)
244d71ae5a4SJacob Faibussowitsch {
245c4762a1bSJed Brown   PetscScalar        v[3];
246c4762a1bSJed Brown   const PetscScalar *tmp;
247c4762a1bSJed Brown   PetscInt           idx[3], i;
248c4762a1bSJed Brown 
2493ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
2509371c9d4SSatish Balay   idx[0] = 0;
2519371c9d4SSatish Balay   idx[1] = 1;
2529371c9d4SSatish Balay   idx[2] = 2;
2539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &tmp));
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   i    = 0;
2569371c9d4SSatish Balay   v[0] = 2.0;
2579371c9d4SSatish Balay   v[1] = 1.0;
2589371c9d4SSatish Balay   v[2] = 0.0;
2599566063dSJacob Faibussowitsch   PetscCall(MatSetValues(BB, 1, &i, 3, idx, v, INSERT_VALUES));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   i    = 1;
2629371c9d4SSatish Balay   v[0] = 1.0;
2639371c9d4SSatish Balay   v[1] = 2.0;
2649371c9d4SSatish Balay   v[2] = 1.0;
2659566063dSJacob Faibussowitsch   PetscCall(MatSetValues(BB, 1, &i, 3, idx, v, INSERT_VALUES));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   i    = 2;
2689371c9d4SSatish Balay   v[0] = 0.0;
2699371c9d4SSatish Balay   v[1] = 1.0;
2709371c9d4SSatish Balay   v[2] = 2.0;
2719566063dSJacob Faibussowitsch   PetscCall(MatSetValues(BB, 1, &i, 3, idx, v, INSERT_VALUES));
272c4762a1bSJed Brown 
2739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(BB, MAT_FINAL_ASSEMBLY));
2749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(BB, MAT_FINAL_ASSEMBLY));
275c4762a1bSJed Brown 
276c4762a1bSJed Brown   if (A != BB) {
2779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2789566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
279c4762a1bSJed Brown   }
2809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &tmp));
281c4762a1bSJed Brown 
2823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
283c4762a1bSJed Brown }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown /*
286c4762a1bSJed Brown       The exact solutions
287c4762a1bSJed Brown */
288d71ae5a4SJacob Faibussowitsch PetscReal solx(PetscReal t)
289d71ae5a4SJacob Faibussowitsch {
2909371c9d4SSatish Balay   return PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / 2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / (2.0 * PetscSqrtReal(2.0)) + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / 2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / (2.0 * PetscSqrtReal(2.0));
291c4762a1bSJed Brown }
292c4762a1bSJed Brown 
293d71ae5a4SJacob Faibussowitsch PetscReal soly(PetscReal t)
294d71ae5a4SJacob Faibussowitsch {
2959371c9d4SSatish Balay   return PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / 2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / PetscSqrtReal(2.0) + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / 2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / PetscSqrtReal(2.0);
296c4762a1bSJed Brown }
297c4762a1bSJed Brown 
298d71ae5a4SJacob Faibussowitsch PetscReal solz(PetscReal t)
299d71ae5a4SJacob Faibussowitsch {
3009371c9d4SSatish Balay   return PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / 2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / (2.0 * PetscSqrtReal(2.0)) + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / 2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / (2.0 * PetscSqrtReal(2.0));
301c4762a1bSJed Brown }
302c4762a1bSJed Brown 
303c4762a1bSJed Brown /*TEST
304c4762a1bSJed Brown 
305c4762a1bSJed Brown     test:
306c4762a1bSJed Brown       suffix: euler
307*eaf6e66dSStefano Zampini       args: -ts_type euler -nest {{0 1}}
308c4762a1bSJed Brown       requires: !single
309c4762a1bSJed Brown 
310c4762a1bSJed Brown     test:
311c4762a1bSJed Brown       suffix: beuler
312*eaf6e66dSStefano Zampini       args: -ts_type beuler -nest {{0 1}}
313*eaf6e66dSStefano Zampini       requires: !single
314*eaf6e66dSStefano Zampini 
315*eaf6e66dSStefano Zampini     test:
316*eaf6e66dSStefano Zampini       suffix: rk
317*eaf6e66dSStefano Zampini       args: -ts_type rk -nest {{0 1}} -ts_adapt_monitor
318c4762a1bSJed Brown       requires: !single
319c4762a1bSJed Brown 
320c4762a1bSJed Brown TEST*/
321