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