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