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