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*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(Xdot,&xdot)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(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*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(Xdot,&xdot)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(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*5f80ce2aSJacob Faibussowitsch CHKERRQ(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*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 71c4762a1bSJed Brown 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 74c4762a1bSJed Brown if (A != B) { 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(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*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 91c4762a1bSJed Brown J[0][0] = 0; 92c4762a1bSJed Brown J[1][0] = (1.-x[0]*x[0])*x[1]-x[0]; 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 95c4762a1bSJed Brown 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 PetscErrorCode ierr; 105c4762a1bSJed Brown const PetscScalar *x; 106c4762a1bSJed Brown PetscReal tfinal, dt; 107c4762a1bSJed Brown User user = (User)ctx; 108c4762a1bSJed Brown Vec interpolatedX; 109c4762a1bSJed Brown 110c4762a1bSJed Brown PetscFunctionBeginUser; 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTimeStep(ts,&dt)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetMaxTime(ts,&tfinal)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&interpolatedX)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSInterpolate(ts,user->next_output,interpolatedX)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(interpolatedX,&x)); 118c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n", 119c4762a1bSJed Brown user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]), 120c4762a1bSJed Brown (double)PetscRealPart(x[1]));CHKERRQ(ierr); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(interpolatedX,&x)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&interpolatedX)); 123c4762a1bSJed Brown user->next_output += 0.1; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown PetscFunctionReturn(0); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown int main(int argc,char **argv) 129c4762a1bSJed Brown { 130c4762a1bSJed Brown TS ts; 131c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 132c4762a1bSJed Brown PetscScalar *x_ptr; 133c4762a1bSJed Brown PetscMPIInt size; 134c4762a1bSJed Brown struct _n_User user; 135c4762a1bSJed Brown PetscInt rows,cols; 136c4762a1bSJed Brown PetscErrorCode ierr; 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 139c4762a1bSJed Brown Initialize program 140c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 141c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 142c4762a1bSJed Brown 143*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1443c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 147c4762a1bSJed Brown Set runtime options 148c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 149c4762a1bSJed Brown user.next_output = 0.0; 150c4762a1bSJed Brown user.mu = 1.0e6; 151c4762a1bSJed Brown user.steps = 0; 152c4762a1bSJed Brown user.ftime = 0.5; 153c4762a1bSJed Brown user.combined = PETSC_FALSE; 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL)); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-combined",&user.combined,NULL)); 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 159c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 160c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 161c4762a1bSJed Brown rows = 2; 162c4762a1bSJed Brown cols = user.combined ? 3 : 1; 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.Jac)); 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(user.Jac,PETSC_DECIDE,PETSC_DECIDE,2,2)); 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(user.Jac)); 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(user.Jac)); 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(user.Jac,&user.x,NULL)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 170c4762a1bSJed Brown Create timestepping solver context 171c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSBEULER)); 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,IFunction,&user)); 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,user.Jac,user.Jac,IJacobian,&user)); 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,user.ftime)); 178c4762a1bSJed Brown if (monitor) { 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts,Monitor,&user,NULL)); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 183c4762a1bSJed Brown Set initial conditions 184c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 185*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user.x,&x_ptr)); 186c4762a1bSJed Brown x_ptr[0] = 2.0; x_ptr[1] = -0.66666654321; 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user.x,&x_ptr)); 188*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,1.0/1024.0)); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* Set up forward sensitivity */ 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.Jacp)); 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,rows,cols)); 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(user.Jacp)); 194*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(user.Jacp)); 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,rows,cols,NULL,&user.sp)); 196c4762a1bSJed Brown if (user.combined) { 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(user.sp)); 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(user.sp,1.0)); 199c4762a1bSJed Brown } else { 200*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(user.sp)); 201c4762a1bSJed Brown } 202*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSForwardSetSensitivities(ts,cols,user.sp)); 203*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user)); 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 206c4762a1bSJed Brown Set runtime options 207c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 208*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 209c4762a1bSJed Brown 210*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,user.x)); 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&user.ftime)); 212*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&user.steps)); 213*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime)); 214*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n ode solution \n")); 215*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(user.x,PETSC_VIEWER_STDOUT_WORLD)); 216c4762a1bSJed Brown 217c4762a1bSJed Brown if (user.combined) { 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n")); 219c4762a1bSJed Brown } else { 220*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[mu]\n")); 221c4762a1bSJed Brown } 222*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(user.sp,PETSC_VIEWER_STDOUT_WORLD)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 225c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 226c4762a1bSJed Brown are no longer needed. 227c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user.Jac)); 229*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user.sp)); 230*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user.Jacp)); 231*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user.x)); 232*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 233c4762a1bSJed Brown 234c4762a1bSJed Brown ierr = PetscFinalize(); 235c4762a1bSJed Brown return(ierr); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown /*TEST 239c4762a1bSJed Brown 240c4762a1bSJed Brown test: 241c4762a1bSJed Brown args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -combined 242c4762a1bSJed Brown requires: !complex !single 243c4762a1bSJed Brown 244c4762a1bSJed Brown test: 245c4762a1bSJed Brown suffix: 2 246c4762a1bSJed Brown requires: !complex !single 247c4762a1bSJed Brown args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 248c4762a1bSJed Brown 249c4762a1bSJed Brown TEST*/ 250