xref: /petsc/src/ts/tests/ex2.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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;
37c4762a1bSJed Brown 
38327415f7SBarry Smith   PetscFunctionBeginUser;
399566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-time", &time_steps, NULL));
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* set initial conditions */
439566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &global));
449566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(global, PETSC_DECIDE, 3));
459566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(global));
469566063dSJacob Faibussowitsch   PetscCall(Initial(global, NULL));
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   /* make timestep context */
499566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
509566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
519566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, Monitor, NULL, NULL));
52c4762a1bSJed Brown   dt = 0.001;
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /*
55c4762a1bSJed Brown     The user provides the RHS and Jacobian
56c4762a1bSJed Brown   */
579566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL));
589566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
599566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 3, 3));
609566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
619566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
629566063dSJacob Faibussowitsch   PetscCall(RHSJacobian(ts, 0.0, global, A, A, NULL));
639566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL));
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, 3, 3, 3, 3, NULL, &S));
669566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S, MATOP_MULT, (void (*)(void))MyMatMult));
679566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, S, A, RHSJacobian, NULL));
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
709566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
739566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps));
749566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 1));
759566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, global));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, global));
789566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
799566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* free the memories */
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&S));
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
89b122ec5aSJacob Faibussowitsch   return 0;
90c4762a1bSJed Brown }
91c4762a1bSJed Brown 
92d71ae5a4SJacob Faibussowitsch PetscErrorCode MyMatMult(Mat S, Vec x, Vec y)
93d71ae5a4SJacob Faibussowitsch {
94c4762a1bSJed Brown   const PetscScalar *inptr;
95c4762a1bSJed Brown   PetscScalar       *outptr;
96c4762a1bSJed Brown 
977510d9b0SBarry Smith   PetscFunctionBeginUser;
989566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &inptr));
999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y, &outptr));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   outptr[0] = 2.0 * inptr[0] + inptr[1];
102c4762a1bSJed Brown   outptr[1] = inptr[0] + 2.0 * inptr[1] + inptr[2];
103c4762a1bSJed Brown   outptr[2] = inptr[1] + 2.0 * inptr[2];
104c4762a1bSJed Brown 
1059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &inptr));
1069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y, &outptr));
107*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108c4762a1bSJed Brown }
109c4762a1bSJed Brown 
110c4762a1bSJed Brown /* this test problem has initial values (1,1,1).                      */
111d71ae5a4SJacob Faibussowitsch PetscErrorCode Initial(Vec global, void *ctx)
112d71ae5a4SJacob Faibussowitsch {
113c4762a1bSJed Brown   PetscScalar *localptr;
114c4762a1bSJed Brown   PetscInt     i, mybase, myend, locsize;
115c4762a1bSJed Brown 
116*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
117c4762a1bSJed Brown   /* determine starting point of each processor */
1189566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(global, &mybase, &myend));
1199566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global, &locsize));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* Initialize the array */
1229566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(global, &localptr));
123c4762a1bSJed Brown   for (i = 0; i < locsize; i++) localptr[i] = 1.0;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   if (mybase == 0) localptr[0] = 1.0;
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(global, &localptr));
128*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec global, void *ctx)
132d71ae5a4SJacob Faibussowitsch {
133c4762a1bSJed Brown   VecScatter         scatter;
134c4762a1bSJed Brown   IS                 from, to;
135c4762a1bSJed Brown   PetscInt           i, n, *idx;
136c4762a1bSJed Brown   Vec                tmp_vec;
137c4762a1bSJed Brown   const PetscScalar *tmp;
138c4762a1bSJed Brown 
139*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
140c4762a1bSJed Brown   /* Get the size of the vector */
1419566063dSJacob Faibussowitsch   PetscCall(VecGetSize(global, &n));
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* Set the index sets */
1449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idx));
145c4762a1bSJed Brown   for (i = 0; i < n; i++) idx[i] = i;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* Create local sequential vectors */
1489566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &tmp_vec));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* Create scatter context */
1519566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &from));
1529566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &to));
1539566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global, from, tmp_vec, to, &scatter));
1549566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD));
1559566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD));
156c4762a1bSJed Brown 
1579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_vec, &tmp));
1589371c9d4SSatish 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])));
1599371c9d4SSatish 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))));
1609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_vec, &tmp));
1619566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
1629566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
1639566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
1649566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
1659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_vec));
166*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx)
170d71ae5a4SJacob Faibussowitsch {
171c4762a1bSJed Brown   PetscScalar       *outptr;
172c4762a1bSJed Brown   const PetscScalar *inptr;
173c4762a1bSJed Brown   PetscInt           i, n, *idx;
174c4762a1bSJed Brown   IS                 from, to;
175c4762a1bSJed Brown   VecScatter         scatter;
176c4762a1bSJed Brown   Vec                tmp_in, tmp_out;
177c4762a1bSJed Brown 
178*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
179c4762a1bSJed Brown   /* Get the length of parallel vector */
1809566063dSJacob Faibussowitsch   PetscCall(VecGetSize(globalin, &n));
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /* Set the index sets */
1839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idx));
184c4762a1bSJed Brown   for (i = 0; i < n; i++) idx[i] = i;
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* Create local sequential vectors */
1879566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &tmp_in));
1889566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tmp_in, &tmp_out));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /* Create scatter context */
1919566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &from));
1929566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &to));
1939566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(globalin, from, tmp_in, to, &scatter));
1949566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD));
1959566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD));
1969566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /*Extract income array */
1999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_in, &inptr));
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* Extract outcome array*/
2029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(tmp_out, &outptr));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   outptr[0] = 2.0 * inptr[0] + inptr[1];
205c4762a1bSJed Brown   outptr[1] = inptr[0] + 2.0 * inptr[1] + inptr[2];
206c4762a1bSJed Brown   outptr[2] = inptr[1] + 2.0 * inptr[2];
207c4762a1bSJed Brown 
2089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_in, &inptr));
2099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(tmp_out, &outptr));
210c4762a1bSJed Brown 
2119566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(tmp_out, from, globalout, to, &scatter));
2129566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD));
2139566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD));
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* Destroy idx aand scatter */
2169566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
2179566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
2189566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
2199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_in));
2209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_out));
2219566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
222*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
223c4762a1bSJed Brown }
224c4762a1bSJed Brown 
225d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat BB, void *ctx)
226d71ae5a4SJacob Faibussowitsch {
227c4762a1bSJed Brown   PetscScalar        v[3];
228c4762a1bSJed Brown   const PetscScalar *tmp;
229c4762a1bSJed Brown   PetscInt           idx[3], i;
230c4762a1bSJed Brown 
231*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
2329371c9d4SSatish Balay   idx[0] = 0;
2339371c9d4SSatish Balay   idx[1] = 1;
2349371c9d4SSatish Balay   idx[2] = 2;
2359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &tmp));
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   i    = 0;
2389371c9d4SSatish Balay   v[0] = 2.0;
2399371c9d4SSatish Balay   v[1] = 1.0;
2409371c9d4SSatish Balay   v[2] = 0.0;
2419566063dSJacob Faibussowitsch   PetscCall(MatSetValues(BB, 1, &i, 3, idx, v, INSERT_VALUES));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   i    = 1;
2449371c9d4SSatish Balay   v[0] = 1.0;
2459371c9d4SSatish Balay   v[1] = 2.0;
2469371c9d4SSatish Balay   v[2] = 1.0;
2479566063dSJacob Faibussowitsch   PetscCall(MatSetValues(BB, 1, &i, 3, idx, v, INSERT_VALUES));
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   i    = 2;
2509371c9d4SSatish Balay   v[0] = 0.0;
2519371c9d4SSatish Balay   v[1] = 1.0;
2529371c9d4SSatish Balay   v[2] = 2.0;
2539566063dSJacob Faibussowitsch   PetscCall(MatSetValues(BB, 1, &i, 3, idx, v, INSERT_VALUES));
254c4762a1bSJed Brown 
2559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(BB, MAT_FINAL_ASSEMBLY));
2569566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(BB, MAT_FINAL_ASSEMBLY));
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   if (A != BB) {
2599566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2609566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
261c4762a1bSJed Brown   }
2629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &tmp));
263c4762a1bSJed Brown 
264*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
265c4762a1bSJed Brown }
266c4762a1bSJed Brown 
267c4762a1bSJed Brown /*
268c4762a1bSJed Brown       The exact solutions
269c4762a1bSJed Brown */
270d71ae5a4SJacob Faibussowitsch PetscReal solx(PetscReal t)
271d71ae5a4SJacob Faibussowitsch {
2729371c9d4SSatish 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));
273c4762a1bSJed Brown }
274c4762a1bSJed Brown 
275d71ae5a4SJacob Faibussowitsch PetscReal soly(PetscReal t)
276d71ae5a4SJacob Faibussowitsch {
2779371c9d4SSatish 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);
278c4762a1bSJed Brown }
279c4762a1bSJed Brown 
280d71ae5a4SJacob Faibussowitsch PetscReal solz(PetscReal t)
281d71ae5a4SJacob Faibussowitsch {
2829371c9d4SSatish 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));
283c4762a1bSJed Brown }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown /*TEST
286c4762a1bSJed Brown 
287c4762a1bSJed Brown     test:
288c4762a1bSJed Brown       suffix: euler
289c4762a1bSJed Brown       args: -ts_type euler
290c4762a1bSJed Brown       requires: !single
291c4762a1bSJed Brown 
292c4762a1bSJed Brown     test:
293c4762a1bSJed Brown       suffix: beuler
294c4762a1bSJed Brown       args:   -ts_type beuler
295c4762a1bSJed Brown       requires: !single
296c4762a1bSJed Brown 
297c4762a1bSJed Brown TEST*/
298