xref: /petsc/src/ts/tests/ex2.c (revision 9566063d113dddea24716c546802770db7481bc0)
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 
30c4762a1bSJed Brown int main(int argc,char **argv)
31c4762a1bSJed Brown {
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 
38*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
39*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* set initial conditions */
42*9566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&global));
43*9566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(global,PETSC_DECIDE,3));
44*9566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(global));
45*9566063dSJacob Faibussowitsch   PetscCall(Initial(global,NULL));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* make timestep context */
48*9566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
49*9566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
50*9566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts,Monitor,NULL,NULL));
51c4762a1bSJed Brown   dt = 0.001;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /*
54c4762a1bSJed Brown     The user provides the RHS and Jacobian
55c4762a1bSJed Brown   */
56*9566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,NULL));
57*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
58*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3));
59*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
60*9566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
61*9566063dSJacob Faibussowitsch   PetscCall(RHSJacobian(ts,0.0,global,A,A,NULL));
62*9566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,A,A,RHSJacobian,NULL));
63c4762a1bSJed Brown 
64*9566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,3,3,3,3,NULL,&S));
65*9566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MyMatMult));
66*9566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,S,A,RHSJacobian,NULL));
67c4762a1bSJed Brown 
68*9566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
69*9566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
70c4762a1bSJed Brown 
71*9566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
72*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts,time_steps));
73*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,1));
74*9566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,global));
75c4762a1bSJed Brown 
76*9566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,global));
77*9566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
78*9566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* free the memories */
81c4762a1bSJed Brown 
82*9566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
83*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
84*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
85*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&S));
86c4762a1bSJed Brown 
87*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
88b122ec5aSJacob Faibussowitsch   return 0;
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat S,Vec x,Vec y)
92c4762a1bSJed Brown {
93c4762a1bSJed Brown   const PetscScalar  *inptr;
94c4762a1bSJed Brown   PetscScalar        *outptr;
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   PetscFunctionBegin;
97*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&inptr));
98*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y,&outptr));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   outptr[0] = 2.0*inptr[0]+inptr[1];
101c4762a1bSJed Brown   outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2];
102c4762a1bSJed Brown   outptr[2] = inptr[1]+2.0*inptr[2];
103c4762a1bSJed Brown 
104*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&inptr));
105*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y,&outptr));
106c4762a1bSJed Brown   PetscFunctionReturn(0);
107c4762a1bSJed Brown }
108c4762a1bSJed Brown 
109c4762a1bSJed Brown /* this test problem has initial values (1,1,1).                      */
110c4762a1bSJed Brown PetscErrorCode Initial(Vec global,void *ctx)
111c4762a1bSJed Brown {
112c4762a1bSJed Brown   PetscScalar    *localptr;
113c4762a1bSJed Brown   PetscInt       i,mybase,myend,locsize;
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* determine starting point of each processor */
116*9566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(global,&mybase,&myend));
117*9566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global,&locsize));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* Initialize the array */
120*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(global,&localptr));
121c4762a1bSJed Brown   for (i=0; i<locsize; i++) localptr[i] = 1.0;
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   if (mybase == 0) localptr[0]=1.0;
124c4762a1bSJed Brown 
125*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(global,&localptr));
126c4762a1bSJed Brown   return 0;
127c4762a1bSJed Brown }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx)
130c4762a1bSJed Brown {
131c4762a1bSJed Brown   VecScatter        scatter;
132c4762a1bSJed Brown   IS                from,to;
133c4762a1bSJed Brown   PetscInt          i,n,*idx;
134c4762a1bSJed Brown   Vec               tmp_vec;
135c4762a1bSJed Brown   PetscErrorCode    ierr;
136c4762a1bSJed Brown   const PetscScalar *tmp;
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /* Get the size of the vector */
139*9566063dSJacob Faibussowitsch   PetscCall(VecGetSize(global,&n));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Set the index sets */
142*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&idx));
143c4762a1bSJed Brown   for (i=0; i<n; i++) idx[i]=i;
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* Create local sequential vectors */
146*9566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Create scatter context */
149*9566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from));
150*9566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to));
151*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global,from,tmp_vec,to,&scatter));
152*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD));
153*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD));
154c4762a1bSJed Brown 
155*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_vec,&tmp));
156c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e u = %14.6e  %14.6e  %14.6e \n",
157*9566063dSJacob Faibussowitsch                      (double)time,(double)PetscRealPart(tmp[0]),(double)PetscRealPart(tmp[1]),(double)PetscRealPart(tmp[2]));PetscCall(ierr);
158c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e errors = %14.6e  %14.6e  %14.6e \n",
159*9566063dSJacob Faibussowitsch                      (double)time,(double)PetscRealPart(tmp[0]-solx(time)),(double)PetscRealPart(tmp[1]-soly(time)),(double)PetscRealPart(tmp[2]-solz(time)));PetscCall(ierr);
160*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_vec,&tmp));
161*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
162*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
163*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
164*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
165*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_vec));
166c4762a1bSJed Brown   return 0;
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
170c4762a1bSJed Brown {
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 
178c4762a1bSJed Brown   /* Get the length of parallel vector */
179*9566063dSJacob Faibussowitsch   PetscCall(VecGetSize(globalin,&n));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* Set the index sets */
182*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&idx));
183c4762a1bSJed Brown   for (i=0; i<n; i++) idx[i]=i;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   /* Create local sequential vectors */
186*9566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&tmp_in));
187*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tmp_in,&tmp_out));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* Create scatter context */
190*9566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from));
191*9566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to));
192*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(globalin,from,tmp_in,to,&scatter));
193*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD));
194*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD));
195*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /*Extract income array */
198*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_in,&inptr));
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* Extract outcome array*/
201*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(tmp_out,&outptr));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   outptr[0] = 2.0*inptr[0]+inptr[1];
204c4762a1bSJed Brown   outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2];
205c4762a1bSJed Brown   outptr[2] = inptr[1]+2.0*inptr[2];
206c4762a1bSJed Brown 
207*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_in,&inptr));
208*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(tmp_out,&outptr));
209c4762a1bSJed Brown 
210*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(tmp_out,from,globalout,to,&scatter));
211*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD));
212*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /* Destroy idx aand scatter */
215*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
216*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
217*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
218*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_in));
219*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_out));
220*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
221c4762a1bSJed Brown   return 0;
222c4762a1bSJed Brown }
223c4762a1bSJed Brown 
224c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ctx)
225c4762a1bSJed Brown {
226c4762a1bSJed Brown   PetscScalar       v[3];
227c4762a1bSJed Brown   const PetscScalar *tmp;
228c4762a1bSJed Brown   PetscInt          idx[3],i;
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   idx[0]=0; idx[1]=1; idx[2]=2;
231*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&tmp));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   i    = 0;
234c4762a1bSJed Brown   v[0] = 2.0; v[1] = 1.0; v[2] = 0.0;
235*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES));
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   i    = 1;
238c4762a1bSJed Brown   v[0] = 1.0; v[1] = 2.0; v[2] = 1.0;
239*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES));
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   i    = 2;
242c4762a1bSJed Brown   v[0] = 0.0; v[1] = 1.0; v[2] = 2.0;
243*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES));
244c4762a1bSJed Brown 
245*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY));
246*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY));
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   if (A != BB) {
249*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
250*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
251c4762a1bSJed Brown   }
252*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&tmp));
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   return 0;
255c4762a1bSJed Brown }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown /*
258c4762a1bSJed Brown       The exact solutions
259c4762a1bSJed Brown */
260c4762a1bSJed Brown PetscReal solx(PetscReal t)
261c4762a1bSJed Brown {
262c4762a1bSJed Brown   return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) +
263c4762a1bSJed Brown          PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0));
264c4762a1bSJed Brown }
265c4762a1bSJed Brown 
266c4762a1bSJed Brown PetscReal soly(PetscReal t)
267c4762a1bSJed Brown {
268c4762a1bSJed Brown   return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0) +
269c4762a1bSJed Brown          PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0);
270c4762a1bSJed Brown }
271c4762a1bSJed Brown 
272c4762a1bSJed Brown PetscReal solz(PetscReal t)
273c4762a1bSJed Brown {
274c4762a1bSJed Brown   return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) +
275c4762a1bSJed Brown          PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0));
276c4762a1bSJed Brown }
277c4762a1bSJed Brown 
278c4762a1bSJed Brown /*TEST
279c4762a1bSJed Brown 
280c4762a1bSJed Brown     test:
281c4762a1bSJed Brown       suffix: euler
282c4762a1bSJed Brown       args: -ts_type euler
283c4762a1bSJed Brown       requires: !single
284c4762a1bSJed Brown 
285c4762a1bSJed Brown     test:
286c4762a1bSJed Brown       suffix: beuler
287c4762a1bSJed Brown       args:   -ts_type beuler
288c4762a1bSJed Brown       requires: !single
289c4762a1bSJed Brown 
290c4762a1bSJed Brown TEST*/
291