Actual source code: ex49.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2: static char help[] = "Solves the van der Pol equation.\n\
  3: Input parameters include:\n";

  5: /*
  6:    Concepts: TS^time-dependent nonlinear problems
  7:    Concepts: TS^van der Pol equation DAE equivalent
  8:    Processors: 1
  9: */
 10: /* ------------------------------------------------------------------------

 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.

 21:    This is a copy and modification of ex20.c to exactly match a test 
 22:    problem that comes with the Radau5 integrator package.

 24:   ------------------------------------------------------------------------- */

 26:  #include <petscts.h>

 28: typedef struct _n_User *User;
 29: struct _n_User {
 30:   PetscReal mu;
 31:   PetscReal next_output;
 32: };


 35: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 36: {
 37:   PetscErrorCode    ierr;
 38:   User              user = (User)ctx;
 39:   const PetscScalar *x,*xdot;
 40:   PetscScalar       *f;

 43:   VecGetArrayRead(X,&x);
 44:   VecGetArrayRead(Xdot,&xdot);
 45:   VecGetArray(F,&f);
 46:   f[0] = xdot[0] - x[1];
 47:   f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]);
 48:   VecRestoreArrayRead(X,&x);
 49:   VecRestoreArrayRead(Xdot,&xdot);
 50:   VecRestoreArray(F,&f);
 51:   return(0);
 52: }

 54: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 55: {
 56:   PetscErrorCode    ierr;
 57:   User              user     = (User)ctx;
 58:   PetscInt          rowcol[] = {0,1};
 59:   const PetscScalar *x;
 60:   PetscScalar       J[2][2];

 63:   VecGetArrayRead(X,&x);
 64:   J[0][0] = a;     J[0][1] = -1.0;
 65:   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:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 67:   VecRestoreArrayRead(X,&x);

 69:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 70:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 71:   if (A != B) {
 72:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 73:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 74:   }
 75:   return(0);
 76: }


 79: int main(int argc,char **argv)
 80: {
 81:   TS             ts;            /* nonlinear solver */
 82:   Vec            x;             /* solution, residual vectors */
 83:   Mat            A;             /* Jacobian matrix */
 84:   PetscInt       steps;
 85:   PetscReal      ftime   = 2;
 86:   PetscScalar    *x_ptr;
 87:   PetscMPIInt    size;
 88:   struct _n_User user;

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 92:      Initialize program
 93:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 94:   PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
 95:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 96:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

 98:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 99:     Set runtime options
100:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101:   user.next_output = 0.0;
102:   user.mu          = 1.0e6;
103:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);
104:   PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);
105:   PetscOptionsEnd();

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:     Create necessary matrix and vectors, solve same ODE on every process
109:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110:   MatCreate(PETSC_COMM_WORLD,&A);
111:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
112:   MatSetFromOptions(A);
113:   MatSetUp(A);

115:   MatCreateVecs(A,&x,NULL);

117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118:      Create timestepping solver context
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120:   TSCreate(PETSC_COMM_WORLD,&ts);
121:   TSSetType(ts,TSBEULER);
122:   TSSetIFunction(ts,NULL,IFunction,&user);
123:   TSSetIJacobian(ts,A,A,IJacobian,&user);

125:   TSSetMaxTime(ts,ftime);
126:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
127:   TSSetTolerances(ts,1.e-4,NULL,1.e-4,NULL);
128:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129:      Set initial conditions
130:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131:   VecGetArray(x,&x_ptr);
132:   x_ptr[0] = 2.0;   x_ptr[1] = -6.6e-01;
133:   VecRestoreArray(x,&x_ptr);
134:   TSSetTimeStep(ts,.000001);

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:      Set runtime options
138:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139:   TSSetFromOptions(ts);

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142:      Solve nonlinear system
143:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144:   TSSolve(ts,x);
145:   TSGetSolveTime(ts,&ftime);
146:   TSGetStepNumber(ts,&steps);
147:   PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
148:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151:      Free work space.  All PETSc objects should be destroyed when they
152:      are no longer needed.
153:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154:   MatDestroy(&A);
155:   VecDestroy(&x);
156:   TSDestroy(&ts);

158:   PetscFinalize();
159:   return(ierr);
160: }

162: /*TEST

164:     build:
165:       requires: double !complex !define(PETSC_USE_64BIT_INDICES) radau5

167:     test:
168:       args: -ts_monitor_solution -ts_type radau5

170: TEST*/