1*c4762a1bSJed Brown #define c11 1.0 2*c4762a1bSJed Brown #define c12 0 3*c4762a1bSJed Brown #define c21 2.0 4*c4762a1bSJed Brown #define c22 1.0 5*c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\ 6*c4762a1bSJed Brown Input parameters include:\n"; 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown /* 9*c4762a1bSJed Brown Concepts: TS^forward sensitivity analysis for time-dependent nonlinear problems 10*c4762a1bSJed Brown Concepts: TS^van der Pol equation DAE equivalent 11*c4762a1bSJed Brown Processors: 1 12*c4762a1bSJed Brown */ 13*c4762a1bSJed Brown /* ------------------------------------------------------------------------ 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown This code demonstrates how to compute trajectory sensitivties w.r.t. the stiffness parameter mu. 16*c4762a1bSJed Brown 1) Use two vectors s and sp for sensitivities w.r.t. intial values and paraeters respectively. This case requires the original Jacobian matrix and a JacobianP matrix for the only parameter mu. 17*c4762a1bSJed Brown 2) Consider the intial 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' 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown ------------------------------------------------------------------------- */ 20*c4762a1bSJed Brown #include <petscts.h> 21*c4762a1bSJed Brown #include <petsctao.h> 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown typedef struct _n_User *User; 24*c4762a1bSJed Brown struct _n_User { 25*c4762a1bSJed Brown PetscReal mu; 26*c4762a1bSJed Brown PetscReal next_output; 27*c4762a1bSJed Brown PetscBool combined; 28*c4762a1bSJed Brown /* Sensitivity analysis support */ 29*c4762a1bSJed Brown PetscInt steps; 30*c4762a1bSJed Brown PetscReal ftime; 31*c4762a1bSJed Brown Mat Jac; /* Jacobian matrix */ 32*c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix */ 33*c4762a1bSJed Brown Vec x; 34*c4762a1bSJed Brown Mat sp; /* forward sensitivity variables */ 35*c4762a1bSJed Brown }; 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown /* 38*c4762a1bSJed Brown * User-defined routines 39*c4762a1bSJed Brown */ 40*c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 41*c4762a1bSJed Brown { 42*c4762a1bSJed Brown PetscErrorCode ierr; 43*c4762a1bSJed Brown User user = (User)ctx; 44*c4762a1bSJed Brown const PetscScalar *x,*xdot; 45*c4762a1bSJed Brown PetscScalar *f; 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown PetscFunctionBeginUser; 48*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 51*c4762a1bSJed Brown f[0] = xdot[0] - x[1]; 52*c4762a1bSJed Brown f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ; 53*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 56*c4762a1bSJed Brown PetscFunctionReturn(0); 57*c4762a1bSJed Brown } 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 60*c4762a1bSJed Brown { 61*c4762a1bSJed Brown PetscErrorCode ierr; 62*c4762a1bSJed Brown User user = (User)ctx; 63*c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 64*c4762a1bSJed Brown PetscScalar J[2][2]; 65*c4762a1bSJed Brown const PetscScalar *x; 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown PetscFunctionBeginUser; 68*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 69*c4762a1bSJed Brown J[0][0] = a; J[0][1] = -1.0; 70*c4762a1bSJed 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]); 71*c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 76*c4762a1bSJed Brown if (A != B) { 77*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 79*c4762a1bSJed Brown } 80*c4762a1bSJed Brown PetscFunctionReturn(0); 81*c4762a1bSJed Brown } 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx) 84*c4762a1bSJed Brown { 85*c4762a1bSJed Brown User user = (User)ctx; 86*c4762a1bSJed Brown PetscInt row[] = {0,1},col[]={0}; 87*c4762a1bSJed Brown PetscScalar J[2][1]; 88*c4762a1bSJed Brown const PetscScalar *x; 89*c4762a1bSJed Brown PetscErrorCode ierr; 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown PetscFunctionBeginUser; 92*c4762a1bSJed Brown if (user->combined) col[0] = 2; 93*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 94*c4762a1bSJed Brown J[0][0] = 0; 95*c4762a1bSJed Brown J[1][0] = (1.-x[0]*x[0])*x[1]-x[0]; 96*c4762a1bSJed Brown ierr = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 98*c4762a1bSJed Brown 99*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 101*c4762a1bSJed Brown PetscFunctionReturn(0); 102*c4762a1bSJed Brown } 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 105*c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 106*c4762a1bSJed Brown { 107*c4762a1bSJed Brown PetscErrorCode ierr; 108*c4762a1bSJed Brown const PetscScalar *x; 109*c4762a1bSJed Brown PetscReal tfinal, dt; 110*c4762a1bSJed Brown User user = (User)ctx; 111*c4762a1bSJed Brown Vec interpolatedX; 112*c4762a1bSJed Brown 113*c4762a1bSJed Brown PetscFunctionBeginUser; 114*c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 118*c4762a1bSJed Brown ierr = VecDuplicate(X,&interpolatedX);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = TSInterpolate(ts,user->next_output,interpolatedX);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = VecGetArrayRead(interpolatedX,&x);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n", 122*c4762a1bSJed Brown user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]), 123*c4762a1bSJed Brown (double)PetscRealPart(x[1]));CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = VecRestoreArrayRead(interpolatedX,&x);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = VecDestroy(&interpolatedX);CHKERRQ(ierr); 126*c4762a1bSJed Brown user->next_output += 0.1; 127*c4762a1bSJed Brown } 128*c4762a1bSJed Brown PetscFunctionReturn(0); 129*c4762a1bSJed Brown } 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown int main(int argc,char **argv) 132*c4762a1bSJed Brown { 133*c4762a1bSJed Brown TS ts; 134*c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 135*c4762a1bSJed Brown PetscScalar *x_ptr; 136*c4762a1bSJed Brown PetscMPIInt size; 137*c4762a1bSJed Brown struct _n_User user; 138*c4762a1bSJed Brown PetscInt rows,cols; 139*c4762a1bSJed Brown PetscErrorCode ierr; 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142*c4762a1bSJed Brown Initialize program 143*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 145*c4762a1bSJed Brown 146*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 147*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 150*c4762a1bSJed Brown Set runtime options 151*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 152*c4762a1bSJed Brown user.next_output = 0.0; 153*c4762a1bSJed Brown user.mu = 1.0e6; 154*c4762a1bSJed Brown user.steps = 0; 155*c4762a1bSJed Brown user.ftime = 0.5; 156*c4762a1bSJed Brown user.combined = PETSC_FALSE; 157*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-combined",&user.combined,NULL);CHKERRQ(ierr); 160*c4762a1bSJed Brown 161*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 162*c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 163*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 164*c4762a1bSJed Brown rows = 2; 165*c4762a1bSJed Brown cols = user.combined ? 3 : 1; 166*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&user.Jac);CHKERRQ(ierr); 167*c4762a1bSJed Brown ierr = MatSetSizes(user.Jac,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 168*c4762a1bSJed Brown ierr = MatSetFromOptions(user.Jac);CHKERRQ(ierr); 169*c4762a1bSJed Brown ierr = MatSetUp(user.Jac);CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = MatCreateVecs(user.Jac,&user.x,NULL);CHKERRQ(ierr); 171*c4762a1bSJed Brown 172*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 173*c4762a1bSJed Brown Create timestepping solver context 174*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 175*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr); 178*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,user.Jac,user.Jac,IJacobian,&user);CHKERRQ(ierr); 179*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 180*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,user.ftime);CHKERRQ(ierr); 181*c4762a1bSJed Brown if (monitor) { 182*c4762a1bSJed Brown ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr); 183*c4762a1bSJed Brown } 184*c4762a1bSJed Brown 185*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 186*c4762a1bSJed Brown Set initial conditions 187*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 188*c4762a1bSJed Brown ierr = VecGetArray(user.x,&x_ptr);CHKERRQ(ierr); 189*c4762a1bSJed Brown x_ptr[0] = 2.0; x_ptr[1] = -0.66666654321; 190*c4762a1bSJed Brown ierr = VecRestoreArray(user.x,&x_ptr);CHKERRQ(ierr); 191*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,1.0/1024.0);CHKERRQ(ierr); 192*c4762a1bSJed Brown 193*c4762a1bSJed Brown /* Set up forward sensitivity */ 194*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&user.Jacp);CHKERRQ(ierr); 195*c4762a1bSJed Brown ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr); 196*c4762a1bSJed Brown ierr = MatSetFromOptions(user.Jacp);CHKERRQ(ierr); 197*c4762a1bSJed Brown ierr = MatSetUp(user.Jacp);CHKERRQ(ierr); 198*c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,rows,cols,NULL,&user.sp);CHKERRQ(ierr); 199*c4762a1bSJed Brown if (user.combined) { 200*c4762a1bSJed Brown ierr = MatZeroEntries(user.sp);CHKERRQ(ierr); 201*c4762a1bSJed Brown ierr = MatShift(user.sp,1.0);CHKERRQ(ierr); 202*c4762a1bSJed Brown } else { 203*c4762a1bSJed Brown ierr = MatZeroEntries(user.sp);CHKERRQ(ierr); 204*c4762a1bSJed Brown } 205*c4762a1bSJed Brown ierr = TSForwardSetSensitivities(ts,cols,user.sp);CHKERRQ(ierr); 206*c4762a1bSJed Brown ierr = TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user);CHKERRQ(ierr); 207*c4762a1bSJed Brown 208*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 209*c4762a1bSJed Brown Set runtime options 210*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 211*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 212*c4762a1bSJed Brown 213*c4762a1bSJed Brown ierr = TSSolve(ts,user.x);CHKERRQ(ierr); 214*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&user.ftime);CHKERRQ(ierr); 215*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&user.steps);CHKERRQ(ierr); 216*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime);CHKERRQ(ierr); 217*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n ode solution \n");CHKERRQ(ierr); 218*c4762a1bSJed Brown ierr = VecView(user.x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 219*c4762a1bSJed Brown 220*c4762a1bSJed Brown if (user.combined) { 221*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n");CHKERRQ(ierr); 222*c4762a1bSJed Brown } else { 223*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[mu]\n");CHKERRQ(ierr); 224*c4762a1bSJed Brown } 225*c4762a1bSJed Brown ierr = MatView(user.sp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 226*c4762a1bSJed Brown 227*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 228*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 229*c4762a1bSJed Brown are no longer needed. 230*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 231*c4762a1bSJed Brown ierr = MatDestroy(&user.Jac);CHKERRQ(ierr); 232*c4762a1bSJed Brown ierr = MatDestroy(&user.sp);CHKERRQ(ierr); 233*c4762a1bSJed Brown ierr = MatDestroy(&user.Jacp);CHKERRQ(ierr); 234*c4762a1bSJed Brown ierr = VecDestroy(&user.x);CHKERRQ(ierr); 235*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 236*c4762a1bSJed Brown 237*c4762a1bSJed Brown ierr = PetscFinalize(); 238*c4762a1bSJed Brown return(ierr); 239*c4762a1bSJed Brown } 240*c4762a1bSJed Brown 241*c4762a1bSJed Brown /*TEST 242*c4762a1bSJed Brown 243*c4762a1bSJed Brown test: 244*c4762a1bSJed Brown args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -combined 245*c4762a1bSJed Brown requires: !complex !single 246*c4762a1bSJed Brown 247*c4762a1bSJed Brown test: 248*c4762a1bSJed Brown suffix: 2 249*c4762a1bSJed Brown requires: !complex !single 250*c4762a1bSJed Brown args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 251*c4762a1bSJed Brown 252*c4762a1bSJed Brown TEST*/ 253