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 PetscErrorCode ierr; 67 User user = (User)ctx; 68 PetscScalar *f; 69 const PetscScalar *x; 70 71 PetscFunctionBeginUser; 72 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 73 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 74 f[0] = x[1]; 75 f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0]; 76 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 77 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 78 PetscFunctionReturn(0); 79 } 80 81 static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx) 82 { 83 PetscErrorCode ierr; 84 User user = (User)ctx; 85 PetscReal mu = user->mu; 86 PetscInt rowcol[] = {0,1}; 87 PetscScalar J[2][2]; 88 const PetscScalar *x; 89 90 PetscFunctionBeginUser; 91 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 92 J[0][0] = 0; 93 J[1][0] = -2.*mu*x[1]*x[0]-1.; 94 J[0][1] = 1.0; 95 J[1][1] = mu*(1.0-x[0]*x[0]); 96 ierr = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 97 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 98 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 99 if (A != B) { 100 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 101 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 102 } 103 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 104 PetscFunctionReturn(0); 105 } 106 107 static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx) 108 { 109 PetscErrorCode ierr; 110 PetscInt row[] = {0,1},col[]={2}; 111 PetscScalar J[2][1]; 112 const PetscScalar *x; 113 114 PetscFunctionBeginUser; 115 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 116 J[0][0] = 0; 117 J[1][0] = (1.-x[0]*x[0])*x[1]; 118 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 119 ierr = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 120 121 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 122 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 123 PetscFunctionReturn(0); 124 } 125 126 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 127 static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 128 { 129 PetscErrorCode ierr; 130 const PetscScalar *x; 131 PetscReal tfinal, dt, tprev; 132 User user = (User)ctx; 133 134 PetscFunctionBeginUser; 135 ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 136 ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr); 137 ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr); 138 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 139 ierr = 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]));CHKERRQ(ierr); 140 ierr = PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);CHKERRQ(ierr); 141 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 142 PetscFunctionReturn(0); 143 } 144 145 int main(int argc,char **argv) 146 { 147 TS ts; /* nonlinear solver */ 148 Vec x; /* solution, residual vectors */ 149 Mat A; /* Jacobian matrix */ 150 Mat Jacp; /* JacobianP matrix */ 151 PetscInt steps; 152 PetscReal ftime =0.5; 153 PetscBool monitor = PETSC_FALSE; 154 PetscScalar *x_ptr; 155 PetscMPIInt size; 156 struct _n_User user; 157 PetscErrorCode ierr; 158 Mat sp; 159 160 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 161 Initialize program 162 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 163 ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 164 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 165 if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 166 167 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 168 Set runtime options 169 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 170 user.mu = 1; 171 user.next_output = 0.0; 172 173 ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr); 174 ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); 175 176 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 177 Create necessary matrix and vectors, solve same ODE on every process 178 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 179 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 180 ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 181 ierr = MatSetFromOptions(A);CHKERRQ(ierr); 182 ierr = MatSetUp(A);CHKERRQ(ierr); 183 ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr); 184 185 ierr = MatCreate(PETSC_COMM_WORLD,&Jacp);CHKERRQ(ierr); 186 ierr = MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,3);CHKERRQ(ierr); 187 ierr = MatSetFromOptions(Jacp);CHKERRQ(ierr); 188 ierr = MatSetUp(Jacp);CHKERRQ(ierr); 189 190 ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,3,NULL,&sp);CHKERRQ(ierr); 191 ierr = MatZeroEntries(sp);CHKERRQ(ierr); 192 ierr = MatShift(sp,1.0);CHKERRQ(ierr); 193 194 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 195 Create timestepping solver context 196 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 197 ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 198 ierr = TSSetType(ts,TSRK);CHKERRQ(ierr); 199 ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr); 200 /* Set RHS Jacobian for the adjoint integration */ 201 ierr = TSSetRHSJacobian(ts,A,A,RHSJacobian,&user);CHKERRQ(ierr); 202 ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); 203 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 204 if (monitor) { 205 ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr); 206 } 207 ierr = TSForwardSetSensitivities(ts,3,sp);CHKERRQ(ierr); 208 ierr = TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user);CHKERRQ(ierr); 209 210 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 211 Set initial conditions 212 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 213 ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr); 214 215 x_ptr[0] = 2; x_ptr[1] = 0.66666654321; 216 ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr); 217 ierr = TSSetTimeStep(ts,.001);CHKERRQ(ierr); 218 219 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 220 Set runtime options 221 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 222 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 223 224 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 225 Solve nonlinear system 226 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 227 ierr = TSSolve(ts,x);CHKERRQ(ierr); 228 ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 229 ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 230 ierr = PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);CHKERRQ(ierr); 231 ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 232 233 ierr = PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n");CHKERRQ(ierr); 234 ierr = MatView(sp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 235 236 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 237 Free work space. All PETSc objects should be destroyed when they 238 are no longer needed. 239 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 240 ierr = MatDestroy(&A);CHKERRQ(ierr); 241 ierr = MatDestroy(&Jacp);CHKERRQ(ierr); 242 ierr = VecDestroy(&x);CHKERRQ(ierr); 243 ierr = MatDestroy(&sp);CHKERRQ(ierr); 244 ierr = TSDestroy(&ts);CHKERRQ(ierr); 245 ierr = PetscFinalize(); 246 return ierr; 247 } 248 249 /*TEST 250 251 test: 252 args: -monitor 0 -ts_adapt_type none 253 254 TEST*/ 255