xref: /petsc/src/ts/tutorials/ex20fwd.c (revision 2d4ee0429272f36cabb4acc35b5dc00ca17b1c83)
1c4762a1bSJed Brown #define c11 1.0
2c4762a1bSJed Brown #define c12 0
3c4762a1bSJed Brown #define c21 2.0
4c4762a1bSJed Brown #define c22 1.0
5c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\
6c4762a1bSJed Brown Input parameters include:\n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown /*
9c4762a1bSJed Brown    Concepts: TS^forward sensitivity analysis for time-dependent nonlinear problems
10c4762a1bSJed Brown    Concepts: TS^van der Pol equation DAE equivalent
11c4762a1bSJed Brown    Processors: 1
12c4762a1bSJed Brown */
13c4762a1bSJed Brown /* ------------------------------------------------------------------------
14c4762a1bSJed Brown 
15c4762a1bSJed Brown    This code demonstrates how to compute trajectory sensitivties w.r.t. the stiffness parameter mu.
16*2d4ee042Sprj-    1) Use two vectors s and sp for sensitivities w.r.t. initial values and paraeters respectively. This case requires the original Jacobian matrix and a JacobianP matrix for the only parameter mu.
17*2d4ee042Sprj-    2) Consider the initial values to be parameters as well. Then there are three parameters in total. The JacobianP matrix will be combined matrix of the Jacobian matrix and JacobianP matrix in the previous case. This choice can be selected by using command line option '-combined'
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ------------------------------------------------------------------------- */
20c4762a1bSJed Brown #include <petscts.h>
21c4762a1bSJed Brown #include <petsctao.h>
22c4762a1bSJed Brown 
23c4762a1bSJed Brown typedef struct _n_User *User;
24c4762a1bSJed Brown struct _n_User {
25c4762a1bSJed Brown   PetscReal mu;
26c4762a1bSJed Brown   PetscReal next_output;
27c4762a1bSJed Brown   PetscBool combined;
28c4762a1bSJed Brown   /* Sensitivity analysis support */
29c4762a1bSJed Brown   PetscInt  steps;
30c4762a1bSJed Brown   PetscReal ftime;
31c4762a1bSJed Brown   Mat       Jac;                    /* Jacobian matrix */
32c4762a1bSJed Brown   Mat       Jacp;                   /* JacobianP matrix */
33c4762a1bSJed Brown   Vec       x;
34c4762a1bSJed Brown   Mat       sp;                     /* forward sensitivity variables */
35c4762a1bSJed Brown };
36c4762a1bSJed Brown 
37c4762a1bSJed Brown /*
38c4762a1bSJed Brown *  User-defined routines
39c4762a1bSJed Brown */
40c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
41c4762a1bSJed Brown {
42c4762a1bSJed Brown   PetscErrorCode    ierr;
43c4762a1bSJed Brown   User              user = (User)ctx;
44c4762a1bSJed Brown   const PetscScalar *x,*xdot;
45c4762a1bSJed Brown   PetscScalar       *f;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   PetscFunctionBeginUser;
48c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
51c4762a1bSJed Brown   f[0] = xdot[0] - x[1];
52c4762a1bSJed Brown   f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ;
53c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
54c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
55c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
56c4762a1bSJed Brown   PetscFunctionReturn(0);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   PetscErrorCode    ierr;
62c4762a1bSJed Brown   User              user     = (User)ctx;
63c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
64c4762a1bSJed Brown   PetscScalar       J[2][2];
65c4762a1bSJed Brown   const PetscScalar *x;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   PetscFunctionBeginUser;
68c4762a1bSJed Brown   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
69c4762a1bSJed Brown   J[0][0] = a;     J[0][1] =  -1.0;
70c4762a1bSJed Brown   J[1][0] = c21*a + user->mu*(1.0 + 2.0*x[0]*x[1]);   J[1][1] = -c21 + a - user->mu*(1.0-x[0]*x[0]);
71c4762a1bSJed Brown   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
72c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76c4762a1bSJed Brown   if (A != B) {
77c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79c4762a1bSJed Brown   }
80c4762a1bSJed Brown   PetscFunctionReturn(0);
81c4762a1bSJed Brown }
82c4762a1bSJed Brown 
83c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
84c4762a1bSJed Brown {
85c4762a1bSJed Brown   User              user = (User)ctx;
86c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[]={0};
87c4762a1bSJed Brown   PetscScalar       J[2][1];
88c4762a1bSJed Brown   const PetscScalar *x;
89c4762a1bSJed Brown   PetscErrorCode    ierr;
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   PetscFunctionBeginUser;
92c4762a1bSJed Brown   if (user->combined) col[0] = 2;
93c4762a1bSJed Brown   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
94c4762a1bSJed Brown   J[0][0] = 0;
95c4762a1bSJed Brown   J[1][0] = (1.-x[0]*x[0])*x[1]-x[0];
96c4762a1bSJed Brown   ierr    = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
97c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
100c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
105c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
106c4762a1bSJed Brown {
107c4762a1bSJed Brown   PetscErrorCode    ierr;
108c4762a1bSJed Brown   const PetscScalar *x;
109c4762a1bSJed Brown   PetscReal         tfinal, dt;
110c4762a1bSJed Brown   User              user = (User)ctx;
111c4762a1bSJed Brown   Vec               interpolatedX;
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   PetscFunctionBeginUser;
114c4762a1bSJed Brown   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr);
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
118c4762a1bSJed Brown     ierr = VecDuplicate(X,&interpolatedX);CHKERRQ(ierr);
119c4762a1bSJed Brown     ierr = TSInterpolate(ts,user->next_output,interpolatedX);CHKERRQ(ierr);
120c4762a1bSJed Brown     ierr = VecGetArrayRead(interpolatedX,&x);CHKERRQ(ierr);
121c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
122c4762a1bSJed Brown                        user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]),
123c4762a1bSJed Brown                        (double)PetscRealPart(x[1]));CHKERRQ(ierr);
124c4762a1bSJed Brown     ierr = VecRestoreArrayRead(interpolatedX,&x);CHKERRQ(ierr);
125c4762a1bSJed Brown     ierr = VecDestroy(&interpolatedX);CHKERRQ(ierr);
126c4762a1bSJed Brown     user->next_output += 0.1;
127c4762a1bSJed Brown   }
128c4762a1bSJed Brown   PetscFunctionReturn(0);
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown int main(int argc,char **argv)
132c4762a1bSJed Brown {
133c4762a1bSJed Brown   TS             ts;
134c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
135c4762a1bSJed Brown   PetscScalar    *x_ptr;
136c4762a1bSJed Brown   PetscMPIInt    size;
137c4762a1bSJed Brown   struct _n_User user;
138c4762a1bSJed Brown   PetscInt       rows,cols;
139c4762a1bSJed Brown   PetscErrorCode ierr;
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142c4762a1bSJed Brown      Initialize program
143c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
147c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150c4762a1bSJed Brown     Set runtime options
151c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152c4762a1bSJed Brown   user.next_output = 0.0;
153c4762a1bSJed Brown   user.mu          = 1.0e6;
154c4762a1bSJed Brown   user.steps       = 0;
155c4762a1bSJed Brown   user.ftime       = 0.5;
156c4762a1bSJed Brown   user.combined    = PETSC_FALSE;
157c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr);
158c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr);
159c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-combined",&user.combined,NULL);CHKERRQ(ierr);
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
163c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164c4762a1bSJed Brown   rows = 2;
165c4762a1bSJed Brown   cols = user.combined ? 3 : 1;
166c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&user.Jac);CHKERRQ(ierr);
167c4762a1bSJed Brown   ierr = MatSetSizes(user.Jac,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
168c4762a1bSJed Brown   ierr = MatSetFromOptions(user.Jac);CHKERRQ(ierr);
169c4762a1bSJed Brown   ierr = MatSetUp(user.Jac);CHKERRQ(ierr);
170c4762a1bSJed Brown   ierr = MatCreateVecs(user.Jac,&user.x,NULL);CHKERRQ(ierr);
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173c4762a1bSJed Brown      Create timestepping solver context
174c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
176c4762a1bSJed Brown   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
177c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr);
178c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,user.Jac,user.Jac,IJacobian,&user);CHKERRQ(ierr);
179c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
180c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,user.ftime);CHKERRQ(ierr);
181c4762a1bSJed Brown   if (monitor) {
182c4762a1bSJed Brown     ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr);
183c4762a1bSJed Brown   }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186c4762a1bSJed Brown      Set initial conditions
187c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
188c4762a1bSJed Brown   ierr = VecGetArray(user.x,&x_ptr);CHKERRQ(ierr);
189c4762a1bSJed Brown   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
190c4762a1bSJed Brown   ierr = VecRestoreArray(user.x,&x_ptr);CHKERRQ(ierr);
191c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,1.0/1024.0);CHKERRQ(ierr);
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   /* Set up forward sensitivity */
194c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&user.Jacp);CHKERRQ(ierr);
195c4762a1bSJed Brown   ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr);
196c4762a1bSJed Brown   ierr = MatSetFromOptions(user.Jacp);CHKERRQ(ierr);
197c4762a1bSJed Brown   ierr = MatSetUp(user.Jacp);CHKERRQ(ierr);
198c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,rows,cols,NULL,&user.sp);CHKERRQ(ierr);
199c4762a1bSJed Brown   if (user.combined) {
200c4762a1bSJed Brown     ierr = MatZeroEntries(user.sp);CHKERRQ(ierr);
201c4762a1bSJed Brown     ierr = MatShift(user.sp,1.0);CHKERRQ(ierr);
202c4762a1bSJed Brown   } else {
203c4762a1bSJed Brown     ierr = MatZeroEntries(user.sp);CHKERRQ(ierr);
204c4762a1bSJed Brown   }
205c4762a1bSJed Brown   ierr = TSForwardSetSensitivities(ts,cols,user.sp);CHKERRQ(ierr);
206c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user);CHKERRQ(ierr);
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209c4762a1bSJed Brown      Set runtime options
210c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   ierr = TSSolve(ts,user.x);CHKERRQ(ierr);
214c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&user.ftime);CHKERRQ(ierr);
215c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&user.steps);CHKERRQ(ierr);
216c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime);CHKERRQ(ierr);
217c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n ode solution \n");CHKERRQ(ierr);
218c4762a1bSJed Brown   ierr = VecView(user.x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   if (user.combined) {
221c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n");CHKERRQ(ierr);
222c4762a1bSJed Brown   } else {
223c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[mu]\n");CHKERRQ(ierr);
224c4762a1bSJed Brown   }
225c4762a1bSJed Brown   ierr = MatView(user.sp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
229c4762a1bSJed Brown      are no longer needed.
230c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
231c4762a1bSJed Brown   ierr = MatDestroy(&user.Jac);CHKERRQ(ierr);
232c4762a1bSJed Brown   ierr = MatDestroy(&user.sp);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = MatDestroy(&user.Jacp);CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = VecDestroy(&user.x);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   ierr = PetscFinalize();
238c4762a1bSJed Brown   return(ierr);
239c4762a1bSJed Brown }
240c4762a1bSJed Brown 
241c4762a1bSJed Brown /*TEST
242c4762a1bSJed Brown 
243c4762a1bSJed Brown     test:
244c4762a1bSJed Brown       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -combined
245c4762a1bSJed Brown       requires:  !complex !single
246c4762a1bSJed Brown 
247c4762a1bSJed Brown     test:
248c4762a1bSJed Brown       suffix: 2
249c4762a1bSJed Brown       requires: !complex !single
250c4762a1bSJed Brown       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5
251c4762a1bSJed Brown 
252c4762a1bSJed Brown TEST*/
253