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