static char help[] = "Solves the van der Pol equation and demonstrate IMEX.\n\ Input parameters include:\n\ -mu : stiffness parameter\n\n"; /* Concepts: TS^time-dependent nonlinear problems Concepts: TS^van der Pol equation Processors: 1 */ /* ------------------------------------------------------------------------ This program solves the van der Pol equation y'' - \mu ((1-y^2)*y' - y) = 0 (1) on the domain 0 <= x <= 1, with the boundary conditions y(0) = 2, y'(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2), 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. Notes: This code demonstrates the TS solver interface to two variants of linear problems, u_t = f(u,t), namely turning (1) into a system of first order differential equations, [ y' ] = [ z ] [ z' ] [ \mu ((1 - y^2) z - y) ] which then we can write as a vector equation [ u_1' ] = [ u_2 ] (2) [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ] which is now in the desired form of u_t = f(u,t). One way that we can split f(u,t) in (2) is to split by component, [ u_1' ] = [ u_2 ] + [ 0 ] [ u_2' ] [ 0 ] [ \mu ((1 - u_1^2) u_2 - u_1) ] where [ F(u,t) ] = [ u_2 ] [ 0 ] and [ G(u',u,t) ] = [ u_1' ] - [ 0 ] [ u_2' ] [ \mu ((1 - u_1^2) u_2 - u_1) ] Using the definition of the Jacobian of G (from the PETSc user manual), in the equation G(u',u,t) = F(u,t), dG dG J(G) = a * -- - -- du' du where d is the partial derivative. In this example, dG [ 1 ; 0 ] -- = [ ] du' [ 0 ; 1 ] dG [ 0 ; 0 ] -- = [ ] du [ -\mu (2*u_1*u_2 + 1); \mu (1 - u_1^2) ] Hence, [ a ; 0 ] J(G) = [ ] [ \mu (2*u_1*u_2 + 1); a - \mu (1 - u_1^2) ] ------------------------------------------------------------------------- */ #include typedef struct _n_User *User; struct _n_User { PetscReal mu; PetscBool imex; PetscReal next_output; }; /* * User-defined routines */ static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx) { PetscErrorCode ierr; User user = (User)ctx; PetscScalar *f; const PetscScalar *x; PetscFunctionBeginUser; ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); ierr = VecGetArray(F,&f);CHKERRQ(ierr); f[0] = (user->imex ? x[1] : 0); f[1] = 0.0; ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); PetscFunctionReturn(0); } static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) { PetscErrorCode ierr; User user = (User)ctx; const PetscScalar *x,*xdot; PetscScalar *f; PetscFunctionBeginUser; ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); ierr = VecGetArray(F,&f);CHKERRQ(ierr); f[0] = xdot[0] + (user->imex ? 0 : x[1]); f[1] = xdot[1] - user->mu*((1. - x[0]*x[0])*x[1] - x[0]); ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); PetscFunctionReturn(0); } static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) { PetscErrorCode ierr; User user = (User)ctx; PetscReal mu = user->mu; PetscInt rowcol[] = {0,1}; const PetscScalar *x; PetscScalar J[2][2]; PetscFunctionBeginUser; ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); J[0][0] = a; J[0][1] = (user->imex ? 0 : 1.); J[1][0] = mu*(2.*x[0]*x[1]+1.); J[1][1] = a - mu*(1. - x[0]*x[0]); ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); if (A != B) { ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); } PetscFunctionReturn(0); } static PetscErrorCode RegisterMyARK2(void) { PetscErrorCode ierr; PetscFunctionBeginUser; { const PetscReal A[3][3] = {{0,0,0}, {0.41421356237309504880,0,0}, {0.75,0.25,0}}, At[3][3] = {{0,0,0}, {0.12132034355964257320,0.29289321881345247560,0}, {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}}, *bembedt = NULL,*bembed = NULL; ierr = TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembed,0,NULL,NULL);CHKERRQ(ierr); } PetscFunctionReturn(0); } /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) { PetscErrorCode ierr; const PetscScalar *x; PetscReal tfinal, dt; User user = (User)ctx; Vec interpolatedX; PetscFunctionBeginUser; ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr); while (user->next_output <= t && user->next_output <= tfinal) { ierr = VecDuplicate(X,&interpolatedX);CHKERRQ(ierr); ierr = TSInterpolate(ts,user->next_output,interpolatedX);CHKERRQ(ierr); ierr = VecGetArrayRead(interpolatedX,&x);CHKERRQ(ierr); 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); ierr = VecRestoreArrayRead(interpolatedX,&x);CHKERRQ(ierr); ierr = VecDestroy(&interpolatedX);CHKERRQ(ierr); user->next_output += 0.1; } PetscFunctionReturn(0); } int main(int argc,char **argv) { TS ts; /* nonlinear solver */ Vec x; /* solution, residual vectors */ Mat A; /* Jacobian matrix */ PetscInt steps; PetscReal ftime = 0.5; PetscBool monitor = PETSC_FALSE; PetscScalar *x_ptr; PetscMPIInt size; struct _n_User user; PetscErrorCode ierr; /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Initialize program - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); ierr = RegisterMyARK2();CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set runtime options - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ user.mu = 1000.0; user.imex = PETSC_TRUE; user.next_output = 0.0; ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,NULL,"-imex",&user.imex,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create necessary matrix and vectors, solve same ODE on every process - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); ierr = MatSetFromOptions(A);CHKERRQ(ierr); ierr = MatSetUp(A);CHKERRQ(ierr); ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create timestepping solver context - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr); ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr); ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr); ierr = TSSetIJacobian(ts,A,A,IJacobian,&user);CHKERRQ(ierr); ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); if (monitor) { ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr); } /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set initial conditions - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr); x_ptr[0] = 2.0; x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu); ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr); ierr = TSSetTimeStep(ts,0.01);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set runtime options - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Solve nonlinear system - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSSolve(ts,x);CHKERRQ(ierr); ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);CHKERRQ(ierr); ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Free work space. All PETSc objects should be destroyed when they are no longer needed. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = MatDestroy(&A);CHKERRQ(ierr); ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = TSDestroy(&ts);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; } /*TEST test: args: -ts_type arkimex -ts_arkimex_type myark2 -ts_adapt_type none requires: !single TEST*/