1 2 static char help[] = "Solves the van der Pol DAE.\n\ 3 Input parameters include:\n"; 4 5 /* ------------------------------------------------------------------------ 6 7 This program solves the van der Pol DAE 8 y' = -z = f(y,z) (1) 9 0 = y-(z^3/3 - z) = g(y,z) 10 on the domain 0 <= x <= 1, with the boundary conditions 11 y(0) = -2, y'(0) = -2.355301397608119909925287735864250951918 12 This is a nonlinear equation. 13 14 Notes: 15 This code demonstrates the TS solver interface with the Van der Pol DAE, 16 namely it is the case when there is no RHS (meaning the RHS == 0), and the 17 equations are converted to two variants of linear problems, u_t = f(u,t), 18 namely turning (1) into a vector equation in terms of u, 19 20 [ y' + z ] = [ 0 ] 21 [ (z^3/3 - z) - y ] [ 0 ] 22 23 which then we can write as a vector equation 24 25 [ u_1' + u_2 ] = [ 0 ] (2) 26 [ (u_2^3/3 - u_2) - u_1 ] [ 0 ] 27 28 which is now in the desired form of u_t = f(u,t). As this is a DAE, and 29 there is no u_2', there is no need for a split, 30 31 so 32 33 [ F(u',u,t) ] = [ u_1' ] + [ u_2 ] 34 [ 0 ] [ (u_2^3/3 - u_2) - u_1 ] 35 36 Using the definition of the Jacobian of F (from the PETSc user manual), 37 in the equation F(u',u,t) = G(u,t), 38 39 dF dF 40 J(F) = a * -- - -- 41 du' du 42 43 where d is the partial derivative. In this example, 44 45 dF [ 1 ; 0 ] 46 -- = [ ] 47 du' [ 0 ; 0 ] 48 49 dF [ 0 ; 1 ] 50 -- = [ ] 51 du [ -1 ; 1 - u_2^2 ] 52 53 Hence, 54 55 [ a ; -1 ] 56 J(F) = [ ] 57 [ 1 ; u_2^2 - 1 ] 58 59 ------------------------------------------------------------------------- */ 60 61 #include <petscts.h> 62 63 typedef struct _n_User *User; 64 struct _n_User { 65 PetscReal next_output; 66 }; 67 68 /* 69 User-defined routines 70 */ 71 72 static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 73 { 74 PetscScalar *f; 75 const PetscScalar *x,*xdot; 76 77 PetscFunctionBeginUser; 78 PetscCall(VecGetArrayRead(X,&x)); 79 PetscCall(VecGetArrayRead(Xdot,&xdot)); 80 PetscCall(VecGetArray(F,&f)); 81 f[0] = xdot[0] + x[1]; 82 f[1] = (x[1]*x[1]*x[1]/3.0 - x[1])-x[0]; 83 PetscCall(VecRestoreArrayRead(X,&x)); 84 PetscCall(VecRestoreArrayRead(Xdot,&xdot)); 85 PetscCall(VecRestoreArray(F,&f)); 86 PetscFunctionReturn(0); 87 } 88 89 static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 90 { 91 PetscInt rowcol[] = {0,1}; 92 PetscScalar J[2][2]; 93 const PetscScalar *x; 94 95 PetscFunctionBeginUser; 96 PetscCall(VecGetArrayRead(X,&x)); 97 J[0][0] = a; J[0][1] = -1.; 98 J[1][0] = 1.; J[1][1] = -1. + x[1]*x[1]; 99 PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 100 PetscCall(VecRestoreArrayRead(X,&x)); 101 102 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 103 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 104 if (A != B) { 105 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 106 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 107 } 108 PetscFunctionReturn(0); 109 } 110 111 static PetscErrorCode RegisterMyARK2(void) 112 { 113 PetscFunctionBeginUser; 114 { 115 const PetscReal 116 A[3][3] = {{0,0,0}, 117 {0.41421356237309504880,0,0}, 118 {0.75,0.25,0}}, 119 At[3][3] = {{0,0,0}, 120 {0.12132034355964257320,0.29289321881345247560,0}, 121 {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}}, 122 *bembedt = NULL,*bembed = NULL; 123 PetscCall(TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembed,0,NULL,NULL)); 124 } 125 PetscFunctionReturn(0); 126 } 127 128 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 129 static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 130 { 131 const PetscScalar *x; 132 PetscReal tfinal, dt; 133 User user = (User)ctx; 134 Vec interpolatedX; 135 136 PetscFunctionBeginUser; 137 PetscCall(TSGetTimeStep(ts,&dt)); 138 PetscCall(TSGetMaxTime(ts,&tfinal)); 139 140 while (user->next_output <= t && user->next_output <= tfinal) { 141 PetscCall(VecDuplicate(X,&interpolatedX)); 142 PetscCall(TSInterpolate(ts,user->next_output,interpolatedX)); 143 PetscCall(VecGetArrayRead(interpolatedX,&x)); 144 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %3" PetscInt_FMT " 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]))); 145 PetscCall(VecRestoreArrayRead(interpolatedX,&x)); 146 PetscCall(VecDestroy(&interpolatedX)); 147 user->next_output += PetscRealConstant(0.1); 148 } 149 PetscFunctionReturn(0); 150 } 151 152 int main(int argc,char **argv) 153 { 154 TS ts; /* nonlinear solver */ 155 Vec x; /* solution, residual vectors */ 156 Mat A; /* Jacobian matrix */ 157 PetscInt steps; 158 PetscReal ftime = 0.5; 159 PetscBool monitor = PETSC_FALSE; 160 PetscScalar *x_ptr; 161 PetscMPIInt size; 162 struct _n_User user; 163 164 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 165 Initialize program 166 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 167 PetscFunctionBeginUser; 168 PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 169 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 170 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 171 172 PetscCall(RegisterMyARK2()); 173 174 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 175 Set runtime options 176 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 177 178 user.next_output = 0.0; 179 PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 180 181 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 182 Create necessary matrix and vectors, solve same ODE on every process 183 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 184 PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 185 PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2)); 186 PetscCall(MatSetFromOptions(A)); 187 PetscCall(MatSetUp(A)); 188 PetscCall(MatCreateVecs(A,&x,NULL)); 189 190 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 191 Create timestepping solver context 192 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 193 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 194 PetscCall(TSSetType(ts,TSBEULER)); 195 PetscCall(TSSetIFunction(ts,NULL,IFunction,&user)); 196 PetscCall(TSSetIJacobian(ts,A,A,IJacobian,&user)); 197 PetscCall(TSSetMaxTime(ts,ftime)); 198 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 199 if (monitor) { 200 PetscCall(TSMonitorSet(ts,Monitor,&user,NULL)); 201 } 202 203 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 204 Set initial conditions 205 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 206 PetscCall(VecGetArray(x,&x_ptr)); 207 x_ptr[0] = -2; x_ptr[1] = -2.355301397608119909925287735864250951918; 208 PetscCall(VecRestoreArray(x,&x_ptr)); 209 PetscCall(TSSetTimeStep(ts,.001)); 210 211 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 212 Set runtime options 213 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 214 PetscCall(TSSetFromOptions(ts)); 215 216 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 217 Solve nonlinear system 218 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 219 PetscCall(TSSolve(ts,x)); 220 PetscCall(TSGetSolveTime(ts,&ftime)); 221 PetscCall(TSGetStepNumber(ts,&steps)); 222 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"steps %3" PetscInt_FMT ", ftime %g\n",steps,(double)ftime)); 223 PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 224 225 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 226 Free work space. All PETSc objects should be destroyed when they 227 are no longer needed. 228 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 229 PetscCall(MatDestroy(&A)); 230 PetscCall(VecDestroy(&x)); 231 PetscCall(TSDestroy(&ts)); 232 233 PetscCall(PetscFinalize()); 234 return 0; 235 } 236 237 /*TEST 238 239 test: 240 requires: !single 241 suffix: a 242 args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp 243 output_file: output/ex19_pi42.out 244 245 test: 246 requires: !single 247 suffix: b 248 args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp -ts_adapt_dsp_filter PI42 249 output_file: output/ex19_pi42.out 250 251 test: 252 requires: !single 253 suffix: c 254 args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp -ts_adapt_dsp_pid 0.4,0.2 255 output_file: output/ex19_pi42.out 256 257 test: 258 requires: !single 259 suffix: bdf_reject 260 args: -ts_type bdf -ts_dt 0.5 -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false -ts_adapt_monitor 261 262 TEST*/ 263