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; 47*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 48*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot,&xdot)); 49*9566063dSJacob Faibussowitsch PetscCall(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]) ; 52*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 53*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot,&xdot)); 54*9566063dSJacob Faibussowitsch PetscCall(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; 66*9566063dSJacob Faibussowitsch PetscCall(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]); 69*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 70*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 71c4762a1bSJed Brown 72*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 73*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 74c4762a1bSJed Brown if (A != B) { 75*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 76*9566063dSJacob Faibussowitsch PetscCall(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; 90*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 91c4762a1bSJed Brown J[0][0] = 0; 92c4762a1bSJed Brown J[1][0] = (1.-x[0]*x[0])*x[1]-x[0]; 93*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES)); 94*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 95c4762a1bSJed Brown 96*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 97*9566063dSJacob Faibussowitsch PetscCall(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; 110*9566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 111*9566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts,&tfinal)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 114*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&interpolatedX)); 115*9566063dSJacob Faibussowitsch PetscCall(TSInterpolate(ts,user->next_output,interpolatedX)); 116*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolatedX,&x)); 117*9566063dSJacob Faibussowitsch PetscCall(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]), 119b122ec5aSJacob Faibussowitsch (double)PetscRealPart(x[1]))); 120*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolatedX,&x)); 121*9566063dSJacob Faibussowitsch PetscCall(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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 140c4762a1bSJed Brown 141*9566063dSJacob Faibussowitsch PetscCallMPI(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; 152*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 153*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL)); 154*9566063dSJacob Faibussowitsch PetscCall(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; 161*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Jac)); 162*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.Jac,PETSC_DECIDE,PETSC_DECIDE,2,2)); 163*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.Jac)); 164*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.Jac)); 165*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jac,&user.x,NULL)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 168c4762a1bSJed Brown Create timestepping solver context 169c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 170*9566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 171*9566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSBEULER)); 172*9566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,IFunction,&user)); 173*9566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,user.Jac,user.Jac,IJacobian,&user)); 174*9566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 175*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,user.ftime)); 176c4762a1bSJed Brown if (monitor) { 177*9566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,Monitor,&user,NULL)); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Set initial conditions 182c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 183*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.x,&x_ptr)); 184c4762a1bSJed Brown x_ptr[0] = 2.0; x_ptr[1] = -0.66666654321; 185*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.x,&x_ptr)); 186*9566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,1.0/1024.0)); 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* Set up forward sensitivity */ 189*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Jacp)); 190*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,rows,cols)); 191*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.Jacp)); 192*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.Jacp)); 193*9566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,rows,cols,NULL,&user.sp)); 194c4762a1bSJed Brown if (user.combined) { 195*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(user.sp)); 196*9566063dSJacob Faibussowitsch PetscCall(MatShift(user.sp,1.0)); 197c4762a1bSJed Brown } else { 198*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(user.sp)); 199c4762a1bSJed Brown } 200*9566063dSJacob Faibussowitsch PetscCall(TSForwardSetSensitivities(ts,cols,user.sp)); 201*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 204c4762a1bSJed Brown Set runtime options 205c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 206*9566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 207c4762a1bSJed Brown 208*9566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,user.x)); 209*9566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&user.ftime)); 210*9566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&user.steps)); 211*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime)); 212*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n ode solution \n")); 213*9566063dSJacob Faibussowitsch PetscCall(VecView(user.x,PETSC_VIEWER_STDOUT_WORLD)); 214c4762a1bSJed Brown 215c4762a1bSJed Brown if (user.combined) { 216*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n")); 217c4762a1bSJed Brown } else { 218*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[mu]\n")); 219c4762a1bSJed Brown } 220*9566063dSJacob Faibussowitsch PetscCall(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 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 226*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.Jac)); 227*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.sp)); 228*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.Jacp)); 229*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.x)); 230*9566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 231c4762a1bSJed Brown 232*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 233b122ec5aSJacob 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