1 2 static char help[] = "Solves the van der Pol equation.\n\ 3 Input parameters include:\n"; 4 5 /* ------------------------------------------------------------------------ 6 7 This program solves the van der Pol DAE ODE equivalent 8 y' = z (1) 9 z' = \mu ((1-y^2)z-y) 10 on the domain 0 <= x <= 1, with the boundary conditions 11 y(0) = 2, y'(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2), 12 and 13 \mu = 10^6 ( y'(0) ~ -0.6666665432100101). 14 This is a nonlinear equation. The well prepared initial condition gives errors that are not dominated by the first few steps of the method when \mu is large. 15 16 Notes: 17 This code demonstrates the TS solver interface to an ODE -- RHSFunction for explicit form and IFunction for implicit form. 18 19 ------------------------------------------------------------------------- */ 20 21 #include <petscts.h> 22 23 typedef struct _n_User *User; 24 struct _n_User { 25 PetscReal mu; 26 PetscReal next_output; 27 }; 28 29 /* 30 User-defined routines 31 */ 32 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 33 { 34 User user = (User)ctx; 35 PetscScalar *f; 36 const PetscScalar *x; 37 38 PetscFunctionBeginUser; 39 PetscCall(VecGetArrayRead(X,&x)); 40 PetscCall(VecGetArray(F,&f)); 41 f[0] = x[1]; 42 f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0]; 43 PetscCall(VecRestoreArrayRead(X,&x)); 44 PetscCall(VecRestoreArray(F,&f)); 45 PetscFunctionReturn(0); 46 } 47 48 static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 49 { 50 User user = (User)ctx; 51 const PetscScalar *x,*xdot; 52 PetscScalar *f; 53 54 PetscFunctionBeginUser; 55 PetscCall(VecGetArrayRead(X,&x)); 56 PetscCall(VecGetArrayRead(Xdot,&xdot)); 57 PetscCall(VecGetArray(F,&f)); 58 f[0] = xdot[0] - x[1]; 59 f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]); 60 PetscCall(VecRestoreArrayRead(X,&x)); 61 PetscCall(VecRestoreArrayRead(Xdot,&xdot)); 62 PetscCall(VecRestoreArray(F,&f)); 63 PetscFunctionReturn(0); 64 } 65 66 static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 67 { 68 User user = (User)ctx; 69 PetscInt rowcol[] = {0,1}; 70 const PetscScalar *x; 71 PetscScalar J[2][2]; 72 73 PetscFunctionBeginUser; 74 PetscCall(VecGetArrayRead(X,&x)); 75 J[0][0] = a; J[0][1] = -1.0; 76 J[1][0] = user->mu*(2.0*x[0]*x[1] + 1.0); J[1][1] = a - user->mu*(1.0-x[0]*x[0]); 77 PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 78 PetscCall(VecRestoreArrayRead(X,&x)); 79 80 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 81 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 82 if (A != B) { 83 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 84 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 85 } 86 PetscFunctionReturn(0); 87 } 88 89 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 90 static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 91 { 92 const PetscScalar *x; 93 PetscReal tfinal, dt; 94 User user = (User)ctx; 95 Vec interpolatedX; 96 97 PetscFunctionBeginUser; 98 PetscCall(TSGetTimeStep(ts,&dt)); 99 PetscCall(TSGetMaxTime(ts,&tfinal)); 100 101 while (user->next_output <= t && user->next_output <= tfinal) { 102 PetscCall(VecDuplicate(X,&interpolatedX)); 103 PetscCall(TSInterpolate(ts,user->next_output,interpolatedX)); 104 PetscCall(VecGetArrayRead(interpolatedX,&x)); 105 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %" PetscInt_FMT " TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n", 106 (double)user->next_output,step,(double)t,(double)dt, 107 (double)PetscRealPart(x[0]),(double)PetscRealPart(x[1]))); 108 PetscCall(VecRestoreArrayRead(interpolatedX,&x)); 109 PetscCall(VecDestroy(&interpolatedX)); 110 user->next_output += 0.1; 111 } 112 PetscFunctionReturn(0); 113 } 114 115 int main(int argc,char **argv) 116 { 117 TS ts; /* nonlinear solver */ 118 Vec x; /* solution, residual vectors */ 119 Mat A; /* Jacobian matrix */ 120 PetscInt steps; 121 PetscReal ftime = 0.5; 122 PetscBool monitor = PETSC_FALSE,implicitform = PETSC_TRUE; 123 PetscScalar *x_ptr; 124 PetscMPIInt size; 125 struct _n_User user; 126 127 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 128 Initialize program 129 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 130 PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 131 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 132 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 133 134 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135 Set runtime options 136 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137 user.next_output = 0.0; 138 user.mu = 1.0e3; 139 PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 140 PetscCall(PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL)); 141 PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL); 142 PetscCall(PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL)); 143 PetscOptionsEnd(); 144 145 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146 Create necessary matrix and vectors, solve same ODE on every process 147 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 148 PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 149 PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2)); 150 PetscCall(MatSetFromOptions(A)); 151 PetscCall(MatSetUp(A)); 152 153 PetscCall(MatCreateVecs(A,&x,NULL)); 154 155 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156 Create timestepping solver context 157 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 158 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 159 if (implicitform) { 160 PetscCall(TSSetIFunction(ts,NULL,IFunction,&user)); 161 PetscCall(TSSetIJacobian(ts,A,A,IJacobian,&user)); 162 PetscCall(TSSetType(ts,TSBEULER)); 163 } else { 164 PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,&user)); 165 PetscCall(TSSetType(ts,TSRK)); 166 } 167 PetscCall(TSSetMaxTime(ts,ftime)); 168 PetscCall(TSSetTimeStep(ts,0.001)); 169 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 170 if (monitor) { 171 PetscCall(TSMonitorSet(ts,Monitor,&user,NULL)); 172 } 173 174 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 175 Set initial conditions 176 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 177 PetscCall(VecGetArray(x,&x_ptr)); 178 x_ptr[0] = 2.0; 179 x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu); 180 PetscCall(VecRestoreArray(x,&x_ptr)); 181 182 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 183 Set runtime options 184 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 185 PetscCall(TSSetFromOptions(ts)); 186 187 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 188 Solve nonlinear system 189 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 190 PetscCall(TSSolve(ts,x)); 191 PetscCall(TSGetSolveTime(ts,&ftime)); 192 PetscCall(TSGetStepNumber(ts,&steps)); 193 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"steps %" PetscInt_FMT ", ftime %g\n",steps,(double)ftime)); 194 PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 195 196 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 197 Free work space. All PETSc objects should be destroyed when they 198 are no longer needed. 199 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 200 PetscCall(MatDestroy(&A)); 201 PetscCall(VecDestroy(&x)); 202 PetscCall(TSDestroy(&ts)); 203 204 PetscCall(PetscFinalize()); 205 return(0); 206 } 207 208 /*TEST 209 210 test: 211 requires: !single 212 args: -mu 1e6 213 214 test: 215 requires: !single 216 suffix: 2 217 args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp 218 219 test: 220 requires: !single 221 suffix: 3 222 args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp -ts_adapt_dsp_filter H0312 223 224 TEST*/ 225