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