xref: /petsc/src/ts/tutorials/ex19.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves the van der Pol DAE.\n\
3c4762a1bSJed Brown Input parameters include:\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
7c4762a1bSJed Brown    Concepts: TS^van der Pol DAE
8c4762a1bSJed Brown    Processors: 1
9c4762a1bSJed Brown */
10c4762a1bSJed Brown /* ------------------------------------------------------------------------
11c4762a1bSJed Brown 
12c4762a1bSJed Brown    This program solves the van der Pol DAE
13c4762a1bSJed Brown        y' = -z = f(y,z)        (1)
14c4762a1bSJed Brown        0  = y-(z^3/3 - z) = g(y,z)
15c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
16c4762a1bSJed Brown        y(0) = -2, y'(0) = -2.355301397608119909925287735864250951918
17c4762a1bSJed Brown    This is a nonlinear equation.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown    Notes:
20c4762a1bSJed Brown    This code demonstrates the TS solver interface with the Van der Pol DAE,
21c4762a1bSJed Brown    namely it is the case when there is no RHS (meaning the RHS == 0), and the
22c4762a1bSJed Brown    equations are converted to two variants of linear problems, u_t = f(u,t),
23c4762a1bSJed Brown    namely turning (1) into a vector equation in terms of u,
24c4762a1bSJed Brown 
25c4762a1bSJed Brown    [     y' + z      ] = [ 0 ]
26c4762a1bSJed Brown    [ (z^3/3 - z) - y ]   [ 0 ]
27c4762a1bSJed Brown 
28c4762a1bSJed Brown    which then we can write as a vector equation
29c4762a1bSJed Brown 
30c4762a1bSJed Brown    [      u_1' + u_2       ] = [ 0 ]  (2)
31c4762a1bSJed Brown    [ (u_2^3/3 - u_2) - u_1 ]   [ 0 ]
32c4762a1bSJed Brown 
33c4762a1bSJed Brown    which is now in the desired form of u_t = f(u,t). As this is a DAE, and
34c4762a1bSJed Brown    there is no u_2', there is no need for a split,
35c4762a1bSJed Brown 
36c4762a1bSJed Brown    so
37c4762a1bSJed Brown 
385ab1ac2bSHong Zhang    [ F(u',u,t) ] = [ u_1' ] + [         u_2           ]
39c4762a1bSJed Brown                    [  0   ]   [ (u_2^3/3 - u_2) - u_1 ]
40c4762a1bSJed Brown 
415ab1ac2bSHong Zhang    Using the definition of the Jacobian of F (from the PETSc user manual),
425ab1ac2bSHong Zhang    in the equation F(u',u,t) = G(u,t),
43c4762a1bSJed Brown 
445ab1ac2bSHong Zhang               dF   dF
455ab1ac2bSHong Zhang    J(F) = a * -- - --
46c4762a1bSJed Brown               du'  du
47c4762a1bSJed Brown 
48c4762a1bSJed Brown    where d is the partial derivative. In this example,
49c4762a1bSJed Brown 
505ab1ac2bSHong Zhang    dF   [ 1 ; 0 ]
51c4762a1bSJed Brown    -- = [       ]
52c4762a1bSJed Brown    du'  [ 0 ; 0 ]
53c4762a1bSJed Brown 
545ab1ac2bSHong Zhang    dF   [  0 ;      1     ]
55c4762a1bSJed Brown    -- = [                 ]
56c4762a1bSJed Brown    du   [ -1 ; 1 - u_2^2  ]
57c4762a1bSJed Brown 
58c4762a1bSJed Brown    Hence,
59c4762a1bSJed Brown 
60c4762a1bSJed Brown           [ a ;    -1     ]
615ab1ac2bSHong Zhang    J(F) = [               ]
62c4762a1bSJed Brown           [ 1 ; u_2^2 - 1 ]
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   ------------------------------------------------------------------------- */
65c4762a1bSJed Brown 
66c4762a1bSJed Brown #include <petscts.h>
67c4762a1bSJed Brown 
68c4762a1bSJed Brown typedef struct _n_User *User;
69c4762a1bSJed Brown struct _n_User {
70c4762a1bSJed Brown   PetscReal next_output;
71c4762a1bSJed Brown };
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /*
740e3d61c9SBarry Smith    User-defined routines
75c4762a1bSJed Brown */
76c4762a1bSJed Brown 
77c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
78c4762a1bSJed Brown {
79c4762a1bSJed Brown   PetscScalar       *f;
80c4762a1bSJed Brown   const PetscScalar *x,*xdot;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   PetscFunctionBeginUser;
83*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
84*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot,&xdot));
85*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
86c4762a1bSJed Brown   f[0] = xdot[0] + x[1];
87c4762a1bSJed Brown   f[1] = (x[1]*x[1]*x[1]/3.0 - x[1])-x[0];
88*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
89*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot,&xdot));
90*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
91c4762a1bSJed Brown   PetscFunctionReturn(0);
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
95c4762a1bSJed Brown {
96c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
97c4762a1bSJed Brown   PetscScalar       J[2][2];
98c4762a1bSJed Brown   const PetscScalar *x;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   PetscFunctionBeginUser;
101*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
102c4762a1bSJed Brown   J[0][0] = a;    J[0][1] = -1.;
103c4762a1bSJed Brown   J[1][0] = 1.;   J[1][1] = -1. + x[1]*x[1];
104*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
105*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
106c4762a1bSJed Brown 
107*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
108*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
109c4762a1bSJed Brown   if (A != B) {
110*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
111*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
112c4762a1bSJed Brown   }
113c4762a1bSJed Brown   PetscFunctionReturn(0);
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown static PetscErrorCode RegisterMyARK2(void)
117c4762a1bSJed Brown {
118c4762a1bSJed Brown   PetscFunctionBeginUser;
119c4762a1bSJed Brown   {
120c4762a1bSJed Brown     const PetscReal
121c4762a1bSJed Brown       A[3][3] = {{0,0,0},
122c4762a1bSJed Brown                  {0.41421356237309504880,0,0},
123c4762a1bSJed Brown                  {0.75,0.25,0}},
124c4762a1bSJed Brown       At[3][3] = {{0,0,0},
125c4762a1bSJed Brown                   {0.12132034355964257320,0.29289321881345247560,0},
126c4762a1bSJed Brown                   {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}},
127c4762a1bSJed Brown     *bembedt = NULL,*bembed = NULL;
128*9566063dSJacob Faibussowitsch     PetscCall(TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembed,0,NULL,NULL));
129c4762a1bSJed Brown   }
130c4762a1bSJed Brown   PetscFunctionReturn(0);
131c4762a1bSJed Brown }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
134c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
135c4762a1bSJed Brown {
136c4762a1bSJed Brown   const PetscScalar *x;
137c4762a1bSJed Brown   PetscReal         tfinal, dt;
138c4762a1bSJed Brown   User              user = (User)ctx;
139c4762a1bSJed Brown   Vec               interpolatedX;
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   PetscFunctionBeginUser;
142*9566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts,&dt));
143*9566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts,&tfinal));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
146*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X,&interpolatedX));
147*9566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts,user->next_output,interpolatedX));
148*9566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedX,&x));
149*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %3D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",(double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]),(double)PetscRealPart(x[1])));
150*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedX,&x));
151*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedX));
152c4762a1bSJed Brown     user->next_output += PetscRealConstant(0.1);
153c4762a1bSJed Brown   }
154c4762a1bSJed Brown   PetscFunctionReturn(0);
155c4762a1bSJed Brown }
156c4762a1bSJed Brown 
157c4762a1bSJed Brown int main(int argc,char **argv)
158c4762a1bSJed Brown {
159c4762a1bSJed Brown   TS             ts;            /* nonlinear solver */
160c4762a1bSJed Brown   Vec            x;             /* solution, residual vectors */
161c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
162c4762a1bSJed Brown   PetscInt       steps;
163c4762a1bSJed Brown   PetscReal      ftime   = 0.5;
164c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
165c4762a1bSJed Brown   PetscScalar    *x_ptr;
166c4762a1bSJed Brown   PetscMPIInt    size;
167c4762a1bSJed Brown   struct _n_User user;
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170c4762a1bSJed Brown      Initialize program
171c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
172*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
173*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1743c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
175c4762a1bSJed Brown 
176*9566063dSJacob Faibussowitsch   PetscCall(RegisterMyARK2());
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179c4762a1bSJed Brown     Set runtime options
180c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   user.next_output = 0.0;
183*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
187c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
188*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
189*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2));
190*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
191*9566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
192*9566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&x,NULL));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195c4762a1bSJed Brown      Create timestepping solver context
196c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197*9566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
198*9566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSBEULER));
199*9566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,IFunction,&user));
200*9566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,A,A,IJacobian,&user));
201*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime));
202*9566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
203c4762a1bSJed Brown   if (monitor) {
204*9566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,Monitor,&user,NULL));
205c4762a1bSJed Brown   }
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208c4762a1bSJed Brown      Set initial conditions
209c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
210*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&x_ptr));
211c4762a1bSJed Brown   x_ptr[0] = -2;   x_ptr[1] = -2.355301397608119909925287735864250951918;
212*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&x_ptr));
213*9566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,.001));
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216c4762a1bSJed Brown      Set runtime options
217c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
218*9566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221c4762a1bSJed Brown      Solve nonlinear system
222c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
223*9566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,x));
224*9566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
225*9566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
226*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"steps %3D, ftime %g\n",steps,(double)ftime));
227*9566063dSJacob Faibussowitsch   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
231c4762a1bSJed Brown      are no longer needed.
232c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
233*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
234*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
235*9566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
236c4762a1bSJed Brown 
237*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
238b122ec5aSJacob Faibussowitsch   return 0;
239c4762a1bSJed Brown }
240c4762a1bSJed Brown 
241c4762a1bSJed Brown /*TEST
242c4762a1bSJed Brown 
243c4762a1bSJed Brown    test:
244c4762a1bSJed Brown       requires: !single
245c4762a1bSJed Brown       suffix: a
246c4762a1bSJed Brown       args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp
247c4762a1bSJed Brown       output_file: output/ex19_pi42.out
248c4762a1bSJed Brown 
249c4762a1bSJed Brown    test:
250c4762a1bSJed Brown       requires: !single
251c4762a1bSJed Brown       suffix: b
252c4762a1bSJed Brown       args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp -ts_adapt_dsp_filter PI42
253c4762a1bSJed Brown       output_file: output/ex19_pi42.out
254c4762a1bSJed Brown 
255c4762a1bSJed Brown    test:
256c4762a1bSJed Brown       requires: !single
257c4762a1bSJed Brown       suffix: c
258c4762a1bSJed Brown       args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp -ts_adapt_dsp_pid 0.4,0.2
259c4762a1bSJed Brown       output_file: output/ex19_pi42.out
260c4762a1bSJed Brown 
261e5b8ffdfSLisandro Dalcin    test:
262e5b8ffdfSLisandro Dalcin       requires: !single
263e5b8ffdfSLisandro Dalcin       suffix: bdf_reject
264e5b8ffdfSLisandro Dalcin       args: -ts_type bdf -ts_dt 0.5 -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false -ts_adapt_monitor
265e5b8ffdfSLisandro Dalcin 
266c4762a1bSJed Brown TEST*/
267