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