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 405ab1ac2bSHong Zhang [ G(u,t) ] = [ u_2 ] 41c4762a1bSJed Brown [ 0 ] 42c4762a1bSJed Brown 43c4762a1bSJed Brown and 44c4762a1bSJed Brown 455ab1ac2bSHong Zhang [ F(u',u,t) ] = [ u_1' ] - [ 0 ] 46c4762a1bSJed Brown [ u_2' ] [ \mu ((1 - u_1^2) u_2 - u_1) ] 47c4762a1bSJed Brown 485ab1ac2bSHong Zhang Using the definition of the Jacobian of F (from the PETSc user manual), 495ab1ac2bSHong Zhang in the equation F(u',u,t) = G(u,t), 50c4762a1bSJed Brown 515ab1ac2bSHong Zhang dF dF 525ab1ac2bSHong Zhang J(F) = a * -- - -- 53c4762a1bSJed Brown du' du 54c4762a1bSJed Brown 55c4762a1bSJed Brown where d is the partial derivative. In this example, 56c4762a1bSJed Brown 575ab1ac2bSHong Zhang dF [ 1 ; 0 ] 58c4762a1bSJed Brown -- = [ ] 59c4762a1bSJed Brown du' [ 0 ; 1 ] 60c4762a1bSJed Brown 615ab1ac2bSHong 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 ] 685ab1ac2bSHong 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 /* 830e3d61c9SBarry Smith User-defined routines 84c4762a1bSJed Brown */ 85c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 86c4762a1bSJed Brown { 87c4762a1bSJed Brown User user = (User)ctx; 88c4762a1bSJed Brown PetscScalar *f; 89c4762a1bSJed Brown const PetscScalar *x; 90c4762a1bSJed Brown 91c4762a1bSJed Brown PetscFunctionBeginUser; 92*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 93*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 94c4762a1bSJed Brown f[0] = (user->imex ? x[1] : 0); 95c4762a1bSJed Brown f[1] = 0.0; 96*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 97*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 98c4762a1bSJed Brown PetscFunctionReturn(0); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 102c4762a1bSJed Brown { 103c4762a1bSJed Brown User user = (User)ctx; 104c4762a1bSJed Brown const PetscScalar *x,*xdot; 105c4762a1bSJed Brown PetscScalar *f; 106c4762a1bSJed Brown 107c4762a1bSJed Brown PetscFunctionBeginUser; 108*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 109*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot,&xdot)); 110*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 111c4762a1bSJed Brown f[0] = xdot[0] + (user->imex ? 0 : x[1]); 112c4762a1bSJed Brown f[1] = xdot[1] - user->mu*((1. - x[0]*x[0])*x[1] - x[0]); 113*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 114*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot,&xdot)); 115*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 116c4762a1bSJed Brown PetscFunctionReturn(0); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 120c4762a1bSJed Brown { 121c4762a1bSJed Brown User user = (User)ctx; 122c4762a1bSJed Brown PetscReal mu = user->mu; 123c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 124c4762a1bSJed Brown const PetscScalar *x; 125c4762a1bSJed Brown PetscScalar J[2][2]; 126c4762a1bSJed Brown 127c4762a1bSJed Brown PetscFunctionBeginUser; 128*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 129c4762a1bSJed Brown J[0][0] = a; J[0][1] = (user->imex ? 0 : 1.); 130c4762a1bSJed Brown J[1][0] = mu*(2.*x[0]*x[1]+1.); J[1][1] = a - mu*(1. - x[0]*x[0]); 131*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 132*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 133c4762a1bSJed Brown 134*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 135*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 136c4762a1bSJed Brown if (A != B) { 137*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 138*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown PetscFunctionReturn(0); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 143c4762a1bSJed Brown static PetscErrorCode RegisterMyARK2(void) 144c4762a1bSJed Brown { 145c4762a1bSJed Brown PetscFunctionBeginUser; 146c4762a1bSJed Brown { 147c4762a1bSJed Brown const PetscReal 148c4762a1bSJed Brown A[3][3] = {{0,0,0}, 149c4762a1bSJed Brown {0.41421356237309504880,0,0}, 150c4762a1bSJed Brown {0.75,0.25,0}}, 151c4762a1bSJed Brown At[3][3] = {{0,0,0}, 152c4762a1bSJed Brown {0.12132034355964257320,0.29289321881345247560,0}, 153c4762a1bSJed Brown {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}}, 154c4762a1bSJed Brown *bembedt = NULL,*bembed = NULL; 155*9566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembed,0,NULL,NULL)); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown PetscFunctionReturn(0); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 161c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 162c4762a1bSJed Brown { 163c4762a1bSJed Brown const PetscScalar *x; 164c4762a1bSJed Brown PetscReal tfinal, dt; 165c4762a1bSJed Brown User user = (User)ctx; 166c4762a1bSJed Brown Vec interpolatedX; 167c4762a1bSJed Brown 168c4762a1bSJed Brown PetscFunctionBeginUser; 169*9566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 170*9566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts,&tfinal)); 171c4762a1bSJed Brown 172c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 173*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&interpolatedX)); 174*9566063dSJacob Faibussowitsch PetscCall(TSInterpolate(ts,user->next_output,interpolatedX)); 175*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolatedX,&x)); 176*9566063dSJacob Faibussowitsch PetscCall(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]))); 177*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolatedX,&x)); 178*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&interpolatedX)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown user->next_output += 0.1; 181c4762a1bSJed Brown } 182c4762a1bSJed Brown PetscFunctionReturn(0); 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185c4762a1bSJed Brown int main(int argc,char **argv) 186c4762a1bSJed Brown { 187c4762a1bSJed Brown TS ts; /* nonlinear solver */ 188c4762a1bSJed Brown Vec x; /* solution, residual vectors */ 189c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 190c4762a1bSJed Brown PetscInt steps; 191c4762a1bSJed Brown PetscReal ftime = 0.5; 192c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 193c4762a1bSJed Brown PetscScalar *x_ptr; 194c4762a1bSJed Brown PetscMPIInt size; 195c4762a1bSJed Brown struct _n_User user; 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 198c4762a1bSJed Brown Initialize program 199c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 200*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 201*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 2023c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 203c4762a1bSJed Brown 204*9566063dSJacob Faibussowitsch PetscCall(RegisterMyARK2()); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 207c4762a1bSJed Brown Set runtime options 208c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 209c4762a1bSJed Brown user.mu = 1000.0; 210c4762a1bSJed Brown user.imex = PETSC_TRUE; 211c4762a1bSJed Brown user.next_output = 0.0; 212c4762a1bSJed Brown 213*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL)); 214*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-imex",&user.imex,NULL)); 215*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 218c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 219c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 220*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 221*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2)); 222*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 223*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 224*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&x,NULL)); 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 227c4762a1bSJed Brown Create timestepping solver context 228c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 229*9566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 230*9566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSBEULER)); 231*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,&user)); 232*9566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,IFunction,&user)); 233*9566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,A,A,IJacobian,&user)); 234*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,ftime)); 235*9566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 236c4762a1bSJed Brown if (monitor) { 237*9566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,Monitor,&user,NULL)); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 241c4762a1bSJed Brown Set initial conditions 242c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 243*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&x_ptr)); 244c4762a1bSJed Brown x_ptr[0] = 2.0; 245c4762a1bSJed Brown x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu); 246*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&x_ptr)); 247*9566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,0.01)); 248c4762a1bSJed Brown 249c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 250c4762a1bSJed Brown Set runtime options 251c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 252*9566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 255c4762a1bSJed Brown Solve nonlinear system 256c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 257*9566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,x)); 258*9566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ftime)); 259*9566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 260*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime)); 261*9566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 264c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 265c4762a1bSJed Brown are no longer needed. 266c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 267*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 268*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 269*9566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 270c4762a1bSJed Brown 271*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 272b122ec5aSJacob Faibussowitsch return 0; 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275c4762a1bSJed Brown /*TEST 276c4762a1bSJed Brown 277c4762a1bSJed Brown test: 278c4762a1bSJed Brown args: -ts_type arkimex -ts_arkimex_type myark2 -ts_adapt_type none 279c4762a1bSJed Brown requires: !single 280c4762a1bSJed Brown 281c4762a1bSJed Brown TEST*/ 282