1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\ 3*c4762a1bSJed Brown Input parameters include:\n"; 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown /* 6*c4762a1bSJed Brown Concepts: TS^time-dependent nonlinear problems 7*c4762a1bSJed Brown Concepts: TS^van der Pol equation DAE equivalent 8*c4762a1bSJed Brown Processors: 1 9*c4762a1bSJed Brown */ 10*c4762a1bSJed Brown /* ------------------------------------------------------------------------ 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown This program solves the van der Pol DAE ODE equivalent 13*c4762a1bSJed Brown y' = z (1) 14*c4762a1bSJed Brown z' = mu[(1-y^2)z-y] 15*c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions 16*c4762a1bSJed Brown y(0) = 2, y'(0) = -6.6e-01, 17*c4762a1bSJed Brown and 18*c4762a1bSJed Brown mu = 10^6. 19*c4762a1bSJed Brown This is a nonlinear equation. 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown This is a copy and modification of ex20.c to exactly match a test 22*c4762a1bSJed Brown problem that comes with the Radau5 integrator package. 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown ------------------------------------------------------------------------- */ 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown #include <petscts.h> 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown typedef struct _n_User *User; 29*c4762a1bSJed Brown struct _n_User { 30*c4762a1bSJed Brown PetscReal mu; 31*c4762a1bSJed Brown PetscReal next_output; 32*c4762a1bSJed Brown }; 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 36*c4762a1bSJed Brown { 37*c4762a1bSJed Brown PetscErrorCode ierr; 38*c4762a1bSJed Brown User user = (User)ctx; 39*c4762a1bSJed Brown const PetscScalar *x,*xdot; 40*c4762a1bSJed Brown PetscScalar *f; 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown PetscFunctionBeginUser; 43*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 44*c4762a1bSJed Brown ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 46*c4762a1bSJed Brown f[0] = xdot[0] - x[1]; 47*c4762a1bSJed Brown f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]); 48*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 51*c4762a1bSJed Brown PetscFunctionReturn(0); 52*c4762a1bSJed Brown } 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 55*c4762a1bSJed Brown { 56*c4762a1bSJed Brown PetscErrorCode ierr; 57*c4762a1bSJed Brown User user = (User)ctx; 58*c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 59*c4762a1bSJed Brown const PetscScalar *x; 60*c4762a1bSJed Brown PetscScalar J[2][2]; 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown PetscFunctionBeginUser; 63*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 64*c4762a1bSJed Brown J[0][0] = a; J[0][1] = -1.0; 65*c4762a1bSJed Brown J[1][0] = user->mu*(1.0 + 2.0*x[0]*x[1]); J[1][1] = a - user->mu*(1.0-x[0]*x[0]); 66*c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71*c4762a1bSJed Brown if (A != B) { 72*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 74*c4762a1bSJed Brown } 75*c4762a1bSJed Brown PetscFunctionReturn(0); 76*c4762a1bSJed Brown } 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown int main(int argc,char **argv) 80*c4762a1bSJed Brown { 81*c4762a1bSJed Brown TS ts; /* nonlinear solver */ 82*c4762a1bSJed Brown Vec x; /* solution, residual vectors */ 83*c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 84*c4762a1bSJed Brown PetscInt steps; 85*c4762a1bSJed Brown PetscReal ftime = 2; 86*c4762a1bSJed Brown PetscScalar *x_ptr; 87*c4762a1bSJed Brown PetscMPIInt size; 88*c4762a1bSJed Brown struct _n_User user; 89*c4762a1bSJed Brown PetscErrorCode ierr; 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92*c4762a1bSJed Brown Initialize program 93*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 94*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 95*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 96*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 99*c4762a1bSJed Brown Set runtime options 100*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 101*c4762a1bSJed Brown user.next_output = 0.0; 102*c4762a1bSJed Brown user.mu = 1.0e6; 103*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL); 104*c4762a1bSJed Brown ierr = PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 108*c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 109*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 110*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 114*c4762a1bSJed Brown 115*c4762a1bSJed Brown ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 118*c4762a1bSJed Brown Create timestepping solver context 119*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 120*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,A,IJacobian,&user);CHKERRQ(ierr); 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = TSSetTolerances(ts,1.e-4,NULL,1.e-4,NULL);CHKERRQ(ierr); 128*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129*c4762a1bSJed Brown Set initial conditions 130*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 131*c4762a1bSJed Brown ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr); 132*c4762a1bSJed Brown x_ptr[0] = 2.0; x_ptr[1] = -6.6e-01; 133*c4762a1bSJed Brown ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,.000001);CHKERRQ(ierr); 135*c4762a1bSJed Brown 136*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 137*c4762a1bSJed Brown Set runtime options 138*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 139*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142*c4762a1bSJed Brown Solve nonlinear system 143*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144*c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 152*c4762a1bSJed Brown are no longer needed. 153*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 154*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown ierr = PetscFinalize(); 159*c4762a1bSJed Brown return(ierr); 160*c4762a1bSJed Brown } 161*c4762a1bSJed Brown 162*c4762a1bSJed Brown /*TEST 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown build: 165*c4762a1bSJed Brown requires: double !complex !define(PETSC_USE_64BIT_INDICES) radau5 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown test: 168*c4762a1bSJed Brown args: -ts_monitor_solution -ts_type radau5 169*c4762a1bSJed Brown 170*c4762a1bSJed Brown TEST*/ 171