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