1 static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n\ 2 Input parameters include:\n\ 3 -mu : stiffness parameter\n\n"; 4 5 /* 6 Concepts: TS^time-dependent nonlinear problems 7 Concepts: TS^van der Pol equation 8 Concepts: TS^adjoint sensitivity analysis 9 Processors: 1 10 */ 11 /* ------------------------------------------------------------------------ 12 13 This program solves the van der Pol equation 14 y'' - \mu (1-y^2)*y' + y = 0 (1) 15 on the domain 0 <= x <= 1, with the boundary conditions 16 y(0) = 2, y'(0) = 0, 17 and computes the sensitivities of the final solution w.r.t. initial conditions and parameter \mu with an explicit Runge-Kutta method and its discrete tangent linear model. 18 19 Notes: 20 This code demonstrates the TSForward interface to a system of ordinary differential equations (ODEs) in the form of u_t = f(u,t). 21 22 (1) can be turned into a system of first order ODEs 23 [ y' ] = [ z ] 24 [ z' ] [ \mu (1 - y^2) z - y ] 25 26 which then we can write as a vector equation 27 28 [ u_1' ] = [ u_2 ] (2) 29 [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ] 30 31 which is now in the form of u_t = F(u,t). 32 33 The user provides the right-hand-side function 34 35 [ f(u,t) ] = [ u_2 ] 36 [ \mu (1 - u_1^2) u_2 - u_1 ] 37 38 the Jacobian function 39 40 df [ 0 ; 1 ] 41 -- = [ ] 42 du [ -2 \mu u_1*u_2 - 1; \mu (1 - u_1^2) ] 43 44 and the JacobainP (the Jacobian w.r.t. parameter) function 45 46 df [ 0; 0; 0 ] 47 --- = [ ] 48 d\mu [ 0; 0; (1 - u_1^2) u_2 ] 49 50 ------------------------------------------------------------------------- */ 51 52 #include <petscts.h> 53 #include <petscmat.h> 54 typedef struct _n_User *User; 55 struct _n_User { 56 PetscReal mu; 57 PetscReal next_output; 58 PetscReal tprev; 59 }; 60 61 /* 62 User-defined routines 63 */ 64 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 65 { 66 User user = (User)ctx; 67 PetscScalar *f; 68 const PetscScalar *x; 69 70 PetscFunctionBeginUser; 71 PetscCall(VecGetArrayRead(X,&x)); 72 PetscCall(VecGetArray(F,&f)); 73 f[0] = x[1]; 74 f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0]; 75 PetscCall(VecRestoreArrayRead(X,&x)); 76 PetscCall(VecRestoreArray(F,&f)); 77 PetscFunctionReturn(0); 78 } 79 80 static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx) 81 { 82 User user = (User)ctx; 83 PetscReal mu = user->mu; 84 PetscInt rowcol[] = {0,1}; 85 PetscScalar J[2][2]; 86 const PetscScalar *x; 87 88 PetscFunctionBeginUser; 89 PetscCall(VecGetArrayRead(X,&x)); 90 J[0][0] = 0; 91 J[1][0] = -2.*mu*x[1]*x[0]-1.; 92 J[0][1] = 1.0; 93 J[1][1] = mu*(1.0-x[0]*x[0]); 94 PetscCall(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 95 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 96 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 97 if (A != B) { 98 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 99 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 100 } 101 PetscCall(VecRestoreArrayRead(X,&x)); 102 PetscFunctionReturn(0); 103 } 104 105 static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx) 106 { 107 PetscInt row[] = {0,1},col[]={2}; 108 PetscScalar J[2][1]; 109 const PetscScalar *x; 110 111 PetscFunctionBeginUser; 112 PetscCall(VecGetArrayRead(X,&x)); 113 J[0][0] = 0; 114 J[1][0] = (1.-x[0]*x[0])*x[1]; 115 PetscCall(VecRestoreArrayRead(X,&x)); 116 PetscCall(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES)); 117 118 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 119 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 120 PetscFunctionReturn(0); 121 } 122 123 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 124 static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 125 { 126 const PetscScalar *x; 127 PetscReal tfinal, dt, tprev; 128 User user = (User)ctx; 129 130 PetscFunctionBeginUser; 131 PetscCall(TSGetTimeStep(ts,&dt)); 132 PetscCall(TSGetMaxTime(ts,&tfinal)); 133 PetscCall(TSGetPrevTime(ts,&tprev)); 134 PetscCall(VecGetArrayRead(X,&x)); 135 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D 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]))); 136 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev)); 137 PetscCall(VecRestoreArrayRead(X,&x)); 138 PetscFunctionReturn(0); 139 } 140 141 int main(int argc,char **argv) 142 { 143 TS ts; /* nonlinear solver */ 144 Vec x; /* solution, residual vectors */ 145 Mat A; /* Jacobian matrix */ 146 Mat Jacp; /* JacobianP matrix */ 147 PetscInt steps; 148 PetscReal ftime =0.5; 149 PetscBool monitor = PETSC_FALSE; 150 PetscScalar *x_ptr; 151 PetscMPIInt size; 152 struct _n_User user; 153 Mat sp; 154 155 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156 Initialize program 157 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 158 PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 159 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 160 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 161 162 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 163 Set runtime options 164 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 165 user.mu = 1; 166 user.next_output = 0.0; 167 168 PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL)); 169 PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 170 171 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 172 Create necessary matrix and vectors, solve same ODE on every process 173 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 174 PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 175 PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2)); 176 PetscCall(MatSetFromOptions(A)); 177 PetscCall(MatSetUp(A)); 178 PetscCall(MatCreateVecs(A,&x,NULL)); 179 180 PetscCall(MatCreate(PETSC_COMM_WORLD,&Jacp)); 181 PetscCall(MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,3)); 182 PetscCall(MatSetFromOptions(Jacp)); 183 PetscCall(MatSetUp(Jacp)); 184 185 PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,3,NULL,&sp)); 186 PetscCall(MatZeroEntries(sp)); 187 PetscCall(MatShift(sp,1.0)); 188 189 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 190 Create timestepping solver context 191 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 192 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 193 PetscCall(TSSetType(ts,TSRK)); 194 PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,&user)); 195 /* Set RHS Jacobian for the adjoint integration */ 196 PetscCall(TSSetRHSJacobian(ts,A,A,RHSJacobian,&user)); 197 PetscCall(TSSetMaxTime(ts,ftime)); 198 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 199 if (monitor) { 200 PetscCall(TSMonitorSet(ts,Monitor,&user,NULL)); 201 } 202 PetscCall(TSForwardSetSensitivities(ts,3,sp)); 203 PetscCall(TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user)); 204 205 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 206 Set initial conditions 207 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 208 PetscCall(VecGetArray(x,&x_ptr)); 209 210 x_ptr[0] = 2; x_ptr[1] = 0.66666654321; 211 PetscCall(VecRestoreArray(x,&x_ptr)); 212 PetscCall(TSSetTimeStep(ts,.001)); 213 214 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 215 Set runtime options 216 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 217 PetscCall(TSSetFromOptions(ts)); 218 219 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 220 Solve nonlinear system 221 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 222 PetscCall(TSSolve(ts,x)); 223 PetscCall(TSGetSolveTime(ts,&ftime)); 224 PetscCall(TSGetStepNumber(ts,&steps)); 225 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime)); 226 PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 227 228 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n")); 229 PetscCall(MatView(sp,PETSC_VIEWER_STDOUT_WORLD)); 230 231 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 232 Free work space. All PETSc objects should be destroyed when they 233 are no longer needed. 234 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 235 PetscCall(MatDestroy(&A)); 236 PetscCall(MatDestroy(&Jacp)); 237 PetscCall(VecDestroy(&x)); 238 PetscCall(MatDestroy(&sp)); 239 PetscCall(TSDestroy(&ts)); 240 PetscCall(PetscFinalize()); 241 return 0; 242 } 243 244 /*TEST 245 246 test: 247 args: -monitor 0 -ts_adapt_type none 248 249 TEST*/ 250