static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n\ Input parameters include:\n\ -mu : stiffness parameter\n\n"; /* Concepts: TS^time-dependent nonlinear problems Concepts: TS^van der Pol equation Concepts: TS^adjoint sensitivity analysis Processors: 1 */ /* ------------------------------------------------------------------------ This program solves the van der Pol equation y'' - \mu (1-y^2)*y' + y = 0 (1) on the domain 0 <= x <= 1, with the boundary conditions y(0) = 2, y'(0) = 0, 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. Notes: This code demonstrates the TSForward interface to a system of ordinary differential equations (ODEs) in the form of u_t = f(u,t). (1) can be turned into a system of first order ODEs [ y' ] = [ z ] [ z' ] [ \mu (1 - y^2) z - y ] which then we can write as a vector equation [ u_1' ] = [ u_2 ] (2) [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ] which is now in the form of u_t = F(u,t). The user provides the right-hand-side function [ f(u,t) ] = [ u_2 ] [ \mu (1 - u_1^2) u_2 - u_1 ] the Jacobian function df [ 0 ; 1 ] -- = [ ] du [ -2 \mu u_1*u_2 - 1; \mu (1 - u_1^2) ] and the JacobainP (the Jacobian w.r.t. parameter) function df [ 0; 0; 0 ] --- = [ ] d\mu [ 0; 0; (1 - u_1^2) u_2 ] ------------------------------------------------------------------------- */ #include #include typedef struct _n_User *User; struct _n_User { PetscReal mu; PetscReal next_output; PetscReal tprev; }; /* User-defined routines */ static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx) { User user = (User)ctx; PetscScalar *f; const PetscScalar *x; PetscFunctionBeginUser; CHKERRQ(VecGetArrayRead(X,&x)); CHKERRQ(VecGetArray(F,&f)); f[0] = x[1]; f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0]; CHKERRQ(VecRestoreArrayRead(X,&x)); CHKERRQ(VecRestoreArray(F,&f)); PetscFunctionReturn(0); } static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx) { User user = (User)ctx; PetscReal mu = user->mu; PetscInt rowcol[] = {0,1}; PetscScalar J[2][2]; const PetscScalar *x; PetscFunctionBeginUser; CHKERRQ(VecGetArrayRead(X,&x)); J[0][0] = 0; J[1][0] = -2.*mu*x[1]*x[0]-1.; J[0][1] = 1.0; J[1][1] = mu*(1.0-x[0]*x[0]); CHKERRQ(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); if (A != B) { CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); } CHKERRQ(VecRestoreArrayRead(X,&x)); PetscFunctionReturn(0); } static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx) { PetscInt row[] = {0,1},col[]={2}; PetscScalar J[2][1]; const PetscScalar *x; PetscFunctionBeginUser; CHKERRQ(VecGetArrayRead(X,&x)); J[0][0] = 0; J[1][0] = (1.-x[0]*x[0])*x[1]; CHKERRQ(VecRestoreArrayRead(X,&x)); CHKERRQ(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES)); CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); PetscFunctionReturn(0); } /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) { const PetscScalar *x; PetscReal tfinal, dt, tprev; User user = (User)ctx; PetscFunctionBeginUser; CHKERRQ(TSGetTimeStep(ts,&dt)); CHKERRQ(TSGetMaxTime(ts,&tfinal)); CHKERRQ(TSGetPrevTime(ts,&tprev)); CHKERRQ(VecGetArrayRead(X,&x)); CHKERRQ(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(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev)); CHKERRQ(VecRestoreArrayRead(X,&x)); PetscFunctionReturn(0); } int main(int argc,char **argv) { TS ts; /* nonlinear solver */ Vec x; /* solution, residual vectors */ Mat A; /* Jacobian matrix */ Mat Jacp; /* JacobianP matrix */ PetscInt steps; PetscReal ftime =0.5; PetscBool monitor = PETSC_FALSE; PetscScalar *x_ptr; PetscMPIInt size; struct _n_User user; PetscErrorCode ierr; Mat sp; /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Initialize program - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set runtime options - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ user.mu = 1; user.next_output = 0.0; CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL)); CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create necessary matrix and vectors, solve same ODE on every process - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2)); CHKERRQ(MatSetFromOptions(A)); CHKERRQ(MatSetUp(A)); CHKERRQ(MatCreateVecs(A,&x,NULL)); CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Jacp)); CHKERRQ(MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,3)); CHKERRQ(MatSetFromOptions(Jacp)); CHKERRQ(MatSetUp(Jacp)); CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,3,NULL,&sp)); CHKERRQ(MatZeroEntries(sp)); CHKERRQ(MatShift(sp,1.0)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create timestepping solver context - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); CHKERRQ(TSSetType(ts,TSRK)); CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,&user)); /* Set RHS Jacobian for the adjoint integration */ CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSJacobian,&user)); CHKERRQ(TSSetMaxTime(ts,ftime)); CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); if (monitor) { CHKERRQ(TSMonitorSet(ts,Monitor,&user,NULL)); } CHKERRQ(TSForwardSetSensitivities(ts,3,sp)); CHKERRQ(TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set initial conditions - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(VecGetArray(x,&x_ptr)); x_ptr[0] = 2; x_ptr[1] = 0.66666654321; CHKERRQ(VecRestoreArray(x,&x_ptr)); CHKERRQ(TSSetTimeStep(ts,.001)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set runtime options - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(TSSetFromOptions(ts)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Solve nonlinear system - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(TSSolve(ts,x)); CHKERRQ(TSGetSolveTime(ts,&ftime)); CHKERRQ(TSGetStepNumber(ts,&steps)); CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime)); CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n")); CHKERRQ(MatView(sp,PETSC_VIEWER_STDOUT_WORLD)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Free work space. All PETSc objects should be destroyed when they are no longer needed. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(MatDestroy(&A)); CHKERRQ(MatDestroy(&Jacp)); CHKERRQ(VecDestroy(&x)); CHKERRQ(MatDestroy(&sp)); CHKERRQ(TSDestroy(&ts)); ierr = PetscFinalize(); return ierr; } /*TEST test: args: -monitor 0 -ts_adapt_type none TEST*/