xref: /petsc/src/ts/tutorials/ex20fwd.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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.
162d4ee042Sprj-    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.
172d4ee042Sprj-    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 /*
380e3d61c9SBarry Smith    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   User              user = (User)ctx;
43c4762a1bSJed Brown   const PetscScalar *x,*xdot;
44c4762a1bSJed Brown   PetscScalar       *f;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   PetscFunctionBeginUser;
475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Xdot,&xdot));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
50c4762a1bSJed Brown   f[0] = xdot[0] - x[1];
51c4762a1bSJed Brown   f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ;
525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Xdot,&xdot));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
55c4762a1bSJed Brown   PetscFunctionReturn(0);
56c4762a1bSJed Brown }
57c4762a1bSJed Brown 
58c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
59c4762a1bSJed Brown {
60c4762a1bSJed Brown   User              user     = (User)ctx;
61c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
62c4762a1bSJed Brown   PetscScalar       J[2][2];
63c4762a1bSJed Brown   const PetscScalar *x;
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   PetscFunctionBeginUser;
665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
67c4762a1bSJed Brown   J[0][0] = a;     J[0][1] =  -1.0;
68c4762a1bSJed 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]);
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
71c4762a1bSJed Brown 
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
74c4762a1bSJed Brown   if (A != B) {
755f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
765f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
77c4762a1bSJed Brown   }
78c4762a1bSJed Brown   PetscFunctionReturn(0);
79c4762a1bSJed Brown }
80c4762a1bSJed Brown 
81c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
82c4762a1bSJed Brown {
83c4762a1bSJed Brown   User              user = (User)ctx;
84c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[]={0};
85c4762a1bSJed Brown   PetscScalar       J[2][1];
86c4762a1bSJed Brown   const PetscScalar *x;
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   PetscFunctionBeginUser;
89c4762a1bSJed Brown   if (user->combined) col[0] = 2;
905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
91c4762a1bSJed Brown   J[0][0] = 0;
92c4762a1bSJed Brown   J[1][0] = (1.-x[0]*x[0])*x[1]-x[0];
935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
95c4762a1bSJed Brown 
965f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
98c4762a1bSJed Brown   PetscFunctionReturn(0);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
102c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
103c4762a1bSJed Brown {
104c4762a1bSJed Brown   const PetscScalar *x;
105c4762a1bSJed Brown   PetscReal         tfinal, dt;
106c4762a1bSJed Brown   User              user = (User)ctx;
107c4762a1bSJed Brown   Vec               interpolatedX;
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   PetscFunctionBeginUser;
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTimeStep(ts,&dt));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetMaxTime(ts,&tfinal));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(X,&interpolatedX));
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(TSInterpolate(ts,user->next_output,interpolatedX));
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(interpolatedX,&x));
117*b122ec5aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
118c4762a1bSJed Brown                         user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]),
119*b122ec5aSJacob Faibussowitsch                         (double)PetscRealPart(x[1])));
1205f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(interpolatedX,&x));
1215f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&interpolatedX));
122c4762a1bSJed Brown     user->next_output += 0.1;
123c4762a1bSJed Brown   }
124c4762a1bSJed Brown   PetscFunctionReturn(0);
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown int main(int argc,char **argv)
128c4762a1bSJed Brown {
129c4762a1bSJed Brown   TS             ts;
130c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
131c4762a1bSJed Brown   PetscScalar    *x_ptr;
132c4762a1bSJed Brown   PetscMPIInt    size;
133c4762a1bSJed Brown   struct _n_User user;
134c4762a1bSJed Brown   PetscInt       rows,cols;
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137c4762a1bSJed Brown      Initialize program
138c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
140c4762a1bSJed Brown 
1415f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1423c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown     Set runtime options
146c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147c4762a1bSJed Brown   user.next_output = 0.0;
148c4762a1bSJed Brown   user.mu          = 1.0e6;
149c4762a1bSJed Brown   user.steps       = 0;
150c4762a1bSJed Brown   user.ftime       = 0.5;
151c4762a1bSJed Brown   user.combined    = PETSC_FALSE;
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-combined",&user.combined,NULL));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
158c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159c4762a1bSJed Brown   rows = 2;
160c4762a1bSJed Brown   cols = user.combined ? 3 : 1;
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.Jac));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.Jac,PETSC_DECIDE,PETSC_DECIDE,2,2));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user.Jac));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.Jac));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jac,&user.x,NULL));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168c4762a1bSJed Brown      Create timestepping solver context
169c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSBEULER));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,IFunction,&user));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts,user.Jac,user.Jac,IJacobian,&user));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,user.ftime));
176c4762a1bSJed Brown   if (monitor) {
1775f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(ts,Monitor,&user,NULL));
178c4762a1bSJed Brown   }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181c4762a1bSJed Brown      Set initial conditions
182c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user.x,&x_ptr));
184c4762a1bSJed Brown   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user.x,&x_ptr));
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,1.0/1024.0));
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* Set up forward sensitivity */
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.Jacp));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,rows,cols));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user.Jacp));
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.Jacp));
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,rows,cols,NULL,&user.sp));
194c4762a1bSJed Brown   if (user.combined) {
1955f80ce2aSJacob Faibussowitsch     CHKERRQ(MatZeroEntries(user.sp));
1965f80ce2aSJacob Faibussowitsch     CHKERRQ(MatShift(user.sp,1.0));
197c4762a1bSJed Brown   } else {
1985f80ce2aSJacob Faibussowitsch     CHKERRQ(MatZeroEntries(user.sp));
199c4762a1bSJed Brown   }
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(TSForwardSetSensitivities(ts,cols,user.sp));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204c4762a1bSJed Brown      Set runtime options
205c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
207c4762a1bSJed Brown 
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,user.x));
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&user.ftime));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&user.steps));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n ode solution \n"));
2135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(user.x,PETSC_VIEWER_STDOUT_WORLD));
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   if (user.combined) {
2165f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n"));
217c4762a1bSJed Brown   } else {
2185f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[mu]\n"));
219c4762a1bSJed Brown   }
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(user.sp,PETSC_VIEWER_STDOUT_WORLD));
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
224c4762a1bSJed Brown      are no longer needed.
225c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.Jac));
2275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.sp));
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.Jacp));
2295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.x));
2305f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
231c4762a1bSJed Brown 
232*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
233*b122ec5aSJacob Faibussowitsch   return 0;
234c4762a1bSJed Brown }
235c4762a1bSJed Brown 
236c4762a1bSJed Brown /*TEST
237c4762a1bSJed Brown 
238c4762a1bSJed Brown     test:
239c4762a1bSJed Brown       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -combined
240c4762a1bSJed Brown       requires:  !complex !single
241c4762a1bSJed Brown 
242c4762a1bSJed Brown     test:
243c4762a1bSJed Brown       suffix: 2
244c4762a1bSJed Brown       requires: !complex !single
245c4762a1bSJed Brown       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5
246c4762a1bSJed Brown 
247c4762a1bSJed Brown TEST*/
248