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