1c4762a1bSJed Brown /* 2c4762a1bSJed Brown * ex_vdp.c 3c4762a1bSJed Brown * 4c4762a1bSJed Brown * Created on: Jun 1, 2012 5c4762a1bSJed Brown * Author: Hong Zhang 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation. \n Input parameters include:\n"; 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* 10c4762a1bSJed Brown * Processors:1 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown /* 14c4762a1bSJed Brown * This program solves the van der Pol equation 15c4762a1bSJed Brown * y' = z (1) 16c4762a1bSJed Brown * z' = (((1-y^2)*z-y)/eps (2) 17c4762a1bSJed Brown * on the domain 0<=x<=0.5, with the initial conditions 18c4762a1bSJed Brown * y(0) = 2, 19c4762a1bSJed Brown * z(0) = -2/3 + 10/81*eps - 292/2187*eps^2-1814/19683*eps^3 20c4762a1bSJed Brown * IMEX schemes are applied to the splitted equation 21c4762a1bSJed Brown * [y'] = [z] + [0 ] 22c4762a1bSJed Brown * [z'] [0] [(((1-y^2)*z-y)/eps] 23c4762a1bSJed Brown * 24c4762a1bSJed Brown * F(x)= [z] 25c4762a1bSJed Brown * [0] 26c4762a1bSJed Brown * 27c4762a1bSJed Brown * G(x)= [y'] - [0 ] 28c4762a1bSJed Brown * [z'] [(((1-y^2)*z-y)/eps] 29c4762a1bSJed Brown * 30c4762a1bSJed Brown * JG(x) = G_x + a G_xdot 31c4762a1bSJed Brown */ 32c4762a1bSJed Brown 33c4762a1bSJed Brown #include <petscdmda.h> 34c4762a1bSJed Brown #include <petscts.h> 35c4762a1bSJed Brown 36c4762a1bSJed Brown typedef struct _User *User; 37c4762a1bSJed Brown struct _User { 38c4762a1bSJed Brown PetscReal mu; /*stiffness control coefficient: epsilon*/ 39c4762a1bSJed Brown }; 40c4762a1bSJed Brown 41c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 42c4762a1bSJed Brown static PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*); 43c4762a1bSJed Brown static PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 44c4762a1bSJed Brown 45c4762a1bSJed Brown int main(int argc, char **argv) 46c4762a1bSJed Brown { 47c4762a1bSJed Brown TS ts; 48c4762a1bSJed Brown Vec x; /*solution vector*/ 49c4762a1bSJed Brown Mat A; /*Jacobian*/ 50c4762a1bSJed Brown PetscInt steps,mx,eimex_rowcol[2],two; 51c4762a1bSJed Brown PetscErrorCode ierr; 52c4762a1bSJed Brown PetscScalar *x_ptr; 53c4762a1bSJed Brown PetscReal ftime,dt,norm; 54c4762a1bSJed Brown Vec ref; 55c4762a1bSJed Brown struct _User user; /* user-defined work context */ 56c4762a1bSJed Brown PetscViewer viewer; 57c4762a1bSJed Brown 58c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 59c4762a1bSJed Brown /* Initialize user application context */ 60*1e1ea65dSPierre Jolivet ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"van der Pol options","");CHKERRQ(ierr); 61c4762a1bSJed Brown user.mu = 1e0; 62c4762a1bSJed Brown ierr = PetscOptionsReal("-eps","Stiffness controller","",user.mu,&user.mu,NULL);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 66c4762a1bSJed Brown Set runtime options 67c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 68c4762a1bSJed Brown /* 69c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); 70c4762a1bSJed Brown */ 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 73c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 74c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 75c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 78c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr); 80c4762a1bSJed Brown 81c4762a1bSJed Brown ierr = MatCreateVecs(A,&ref,NULL);CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = VecGetArray(ref,&x_ptr);CHKERRQ(ierr); 83c4762a1bSJed Brown /* 84c4762a1bSJed Brown * [0,1], mu=10^-3 85c4762a1bSJed Brown */ 86c4762a1bSJed Brown /* 87c4762a1bSJed Brown x_ptr[0] = -1.8881254106283; 88c4762a1bSJed Brown x_ptr[1] = 0.7359074233370;*/ 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* 91c4762a1bSJed Brown * [0,0.5],mu=10^-3 92c4762a1bSJed Brown */ 93c4762a1bSJed Brown /* 94c4762a1bSJed Brown x_ptr[0] = 1.596980778659137; 95c4762a1bSJed Brown x_ptr[1] = -1.029103015879544; 96c4762a1bSJed Brown */ 97c4762a1bSJed Brown /* 98c4762a1bSJed Brown * [0,0.5],mu=1 99c4762a1bSJed Brown */ 100c4762a1bSJed Brown x_ptr[0] = 1.619084329683235; 101c4762a1bSJed Brown x_ptr[1] = -0.803530465176385; 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 104c4762a1bSJed Brown Create timestepping solver context 105c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 106c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = TSSetType(ts,TSEIMEX);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,A,IJacobian,&user);CHKERRQ(ierr); 111c4762a1bSJed Brown 112c4762a1bSJed Brown dt = 0.00001; 113c4762a1bSJed Brown ftime = 1.1; 114c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 117c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 118c4762a1bSJed Brown Set initial conditions 119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 120c4762a1bSJed Brown ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr); 121c4762a1bSJed Brown x_ptr[0] = 2.; 122c4762a1bSJed Brown x_ptr[1] = -2./3. + 10./81.*(user.mu) - 292./2187.* (user.mu) * (user.mu) 123c4762a1bSJed Brown -1814./19683.*(user.mu)*(user.mu)*(user.mu); 124c4762a1bSJed Brown ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = VecGetSize(x,&mx);CHKERRQ(ierr); 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 128c4762a1bSJed Brown Set runtime options 129c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 130c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 133c4762a1bSJed Brown Solve nonlinear system 134c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 135c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 136c4762a1bSJed Brown ierr = TSGetTime(ts,&ftime);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 138c4762a1bSJed Brown 139c4762a1bSJed Brown ierr = VecAXPY(x,-1.0,ref);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 142c4762a1bSJed Brown 143c4762a1bSJed Brown eimex_rowcol[0] = 0; eimex_rowcol[1] = 0; two = 2; 144c4762a1bSJed Brown ierr = PetscOptionsGetIntArray(NULL,NULL,"-ts_eimex_row_col",eimex_rowcol,&two,NULL);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"order %11s %18s %37s\n","dt","norm","final solution components 0 and 1");CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"(%D,%D) %10.8f %18.15f %18.15f %18.15f\n",eimex_rowcol[0],eimex_rowcol[1],(double)dt,(double)norm,(double)PetscRealPart(x_ptr[0]),(double)PetscRealPart(x_ptr[1]));CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* Write line in convergence log */ 151c4762a1bSJed Brown ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr); 152c4762a1bSJed Brown ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = PetscViewerFileSetMode(viewer,FILE_MODE_APPEND);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = PetscViewerFileSetName(viewer,"eimex_nonstiff_vdp.txt");CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"%D %D %10.8f %18.15f\n",eimex_rowcol[0],eimex_rowcol[1],(double)dt,(double)norm);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 159c4762a1bSJed Brown Free work space. 160c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 161c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = VecDestroy(&ref);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = PetscFinalize(); 166c4762a1bSJed Brown return ierr; 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 170c4762a1bSJed Brown { 171c4762a1bSJed Brown PetscErrorCode ierr; 172c4762a1bSJed Brown PetscScalar *f; 173c4762a1bSJed Brown const PetscScalar *x; 174c4762a1bSJed Brown 175c4762a1bSJed Brown PetscFunctionBegin; 176c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 178c4762a1bSJed Brown f[0] = x[1]; 179c4762a1bSJed Brown f[1] = 0.0; 180c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 182c4762a1bSJed Brown PetscFunctionReturn(0); 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr) 186c4762a1bSJed Brown { 187c4762a1bSJed Brown User user = (User)ptr; 188c4762a1bSJed Brown PetscScalar *f; 189c4762a1bSJed Brown const PetscScalar *x,*xdot; 190c4762a1bSJed Brown PetscErrorCode ierr; 191c4762a1bSJed Brown 192c4762a1bSJed Brown PetscFunctionBegin; 193c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 196c4762a1bSJed Brown f[0] = xdot[0]; 197c4762a1bSJed Brown f[1] = xdot[1]-((1.-x[0]*x[0])*x[1]-x[0])/user->mu; 198c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 201c4762a1bSJed Brown PetscFunctionReturn(0); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown 204c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ptr) 205c4762a1bSJed Brown { 206c4762a1bSJed Brown PetscErrorCode ierr; 207c4762a1bSJed Brown User user = (User)ptr; 208c4762a1bSJed Brown PetscReal mu = user->mu; 209c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 210c4762a1bSJed Brown PetscScalar J[2][2]; 211c4762a1bSJed Brown const PetscScalar *x; 212c4762a1bSJed Brown 213c4762a1bSJed Brown PetscFunctionBegin; 214c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 215c4762a1bSJed Brown J[0][0] = a; 216c4762a1bSJed Brown J[0][1] = 0; 217c4762a1bSJed Brown J[1][0] = (2.*x[0]*x[1]+1.)/mu; 218c4762a1bSJed Brown J[1][1] = a - (1. - x[0]*x[0])/mu; 219c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 221c4762a1bSJed Brown 222c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 224c4762a1bSJed Brown if (A != B) { 225c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown PetscFunctionReturn(0); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown 231c4762a1bSJed Brown /*TEST 232c4762a1bSJed Brown 233c4762a1bSJed Brown test: 234c4762a1bSJed Brown args: -ts_type eimex -ts_adapt_type none -pc_type lu -ts_dt 0.01 -ts_max_time 10 -ts_eimex_row_col 3,3 -ts_monitor_lg_solution 235c4762a1bSJed Brown requires: x 236c4762a1bSJed Brown 237c4762a1bSJed Brown test: 238c4762a1bSJed Brown suffix: adapt 239c4762a1bSJed Brown args: -ts_type eimex -ts_adapt_type none -pc_type lu -ts_dt 0.01 -ts_max_time 10 -ts_eimex_order_adapt -ts_eimex_max_rows 7 -ts_monitor_lg_solution 240c4762a1bSJed Brown requires: x 241c4762a1bSJed Brown 242c4762a1bSJed Brown test: 243c4762a1bSJed Brown suffix: loop 244c4762a1bSJed Brown args: -ts_type eimex -ts_adapt_type none -pc_type lu -ts_dt {{0.005 0.001 0.0005}separate output} -ts_max_steps {{100 500 1000}separate output} -ts_eimex_row_col {{1,1 2,1 3,1 2,2 3,2 3,3}separate output} 245c4762a1bSJed Brown requires: x 246c4762a1bSJed Brown 247c4762a1bSJed Brown TEST*/ 248