Actual source code: ex19.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: static char help[] = "Solves the van der Pol DAE.\n\
  3: Input parameters include:\n";

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

 12:    This program solves the van der Pol DAE
 13:        y' = -z = f(y,z)        (1)
 14:        0  = y-(z^3/3 - z) = g(y,z)
 15:    on the domain 0 <= x <= 1, with the boundary conditions
 16:        y(0) = -2, y'(0) = -2.355301397608119909925287735864250951918
 17:    This is a nonlinear equation.

 19:    Notes:
 20:    This code demonstrates the TS solver interface with the Van der Pol DAE,
 21:    namely it is the case when there is no RHS (meaning the RHS == 0), and the
 22:    equations are converted to two variants of linear problems, u_t = f(u,t),
 23:    namely turning (1) into a vector equation in terms of u,

 25:    [     y' + z      ] = [ 0 ]
 26:    [ (z^3/3 - z) - y ]   [ 0 ]

 28:    which then we can write as a vector equation

 30:    [      u_1' + u_2       ] = [ 0 ]  (2)
 31:    [ (u_2^3/3 - u_2) - u_1 ]   [ 0 ]

 33:    which is now in the desired form of u_t = f(u,t). As this is a DAE, and
 34:    there is no u_2', there is no need for a split,

 36:    so

 38:    [ G(u',u,t) ] = [ u_1' ] + [         u_2           ]
 39:                    [  0   ]   [ (u_2^3/3 - u_2) - u_1 ]

 41:    Using the definition of the Jacobian of G (from the PETSc user manual),
 42:    in the equation G(u',u,t) = F(u,t),

 44:               dG   dG
 45:    J(G) = a * -- - --
 46:               du'  du

 48:    where d is the partial derivative. In this example,

 50:    dG   [ 1 ; 0 ]
 51:    -- = [       ]
 52:    du'  [ 0 ; 0 ]

 54:    dG   [  0 ;      1     ]
 55:    -- = [                 ]
 56:    du   [ -1 ; 1 - u_2^2  ]

 58:    Hence,

 60:           [ a ;    -1     ]
 61:    J(G) = [               ]
 62:           [ 1 ; u_2^2 - 1 ]

 64:   ------------------------------------------------------------------------- */

 66: #include <petscts.h>

 68: typedef struct _n_User *User;
 69: struct _n_User {
 70:   PetscReal next_output;
 71: };

 73: /*
 74: *  User-defined routines
 75: */

 79: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 80: {
 82:   PetscScalar    *x,*xdot,*f;

 85:   VecGetArray(X,&x);
 86:   VecGetArray(Xdot,&xdot);
 87:   VecGetArray(F,&f);
 88:   f[0] = xdot[0] + x[1];
 89:   f[1] = (x[1]*x[1]*x[1]/3 - x[1])-x[0];
 90:   VecRestoreArray(X,&x);
 91:   VecRestoreArray(Xdot,&xdot);
 92:   VecRestoreArray(F,&f);
 93:   return(0);
 94: }

 98: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 99: {
101:   PetscInt       rowcol[] = {0,1};
102:   PetscScalar    *x,J[2][2];

105:   VecGetArray(X,&x);
106:   J[0][0] = a;    J[0][1] = -1.;
107:   J[1][0] = 1.;   J[1][1] = -1. + x[1]*x[1];
108:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
109:   VecRestoreArray(X,&x);

111:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
112:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
113:   if (A != B) {
114:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
115:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
116:   }
117:   return(0);
118: }

122: static PetscErrorCode RegisterMyARK2(void)
123: {

127:   {
128:     const PetscReal
129:       A[3][3] = {{0,0,0},
130:                  {0.41421356237309504880,0,0},
131:                  {0.75,0.25,0}},
132:       At[3][3] = {{0,0,0},
133:                   {0.12132034355964257320,0.29289321881345247560,0},
134:                   {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}},
135:     *bembedt = NULL,*bembed = NULL;
136:     TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembed,0,NULL,NULL);
137:   }
138:   return(0);
139: }

143: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
144: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
145: {
146:   PetscErrorCode    ierr;
147:   const PetscScalar *x;
148:   PetscReal         tfinal, dt;
149:   User              user = (User)ctx;
150:   Vec               interpolatedX;

153:   TSGetTimeStep(ts,&dt);
154:   TSGetDuration(ts,NULL,&tfinal);

156:   while (user->next_output <= t && user->next_output <= tfinal) {
157:     VecDuplicate(X,&interpolatedX);
158:     TSInterpolate(ts,user->next_output,interpolatedX);
159:     VecGetArrayRead(interpolatedX,&x);
160:     PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",user->next_output,step,t,dt,(double)PetscRealPart(x[0]),(double)PetscRealPart(x[1]));
161:     VecRestoreArrayRead(interpolatedX,&x);
162:     VecDestroy(&interpolatedX);
163:     user->next_output += 0.1;
164:   }
165:   return(0);
166: }

170: int main(int argc,char **argv)
171: {
172:   TS             ts;            /* nonlinear solver */
173:   Vec            x;             /* solution, residual vectors */
174:   Mat            A;             /* Jacobian matrix */
175:   PetscInt       steps;
176:   PetscReal      ftime   = 0.5;
177:   PetscBool      monitor = PETSC_FALSE;
178:   PetscScalar    *x_ptr;
179:   PetscMPIInt    size;
180:   struct _n_User user;

183:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184:      Initialize program
185:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186:   PetscInitialize(&argc,&argv,NULL,help);

188:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
189:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

191:   RegisterMyARK2();

193:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194:     Set runtime options
195:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

197:   user.next_output = 0.0;
198:   PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);

200:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201:     Create necessary matrix and vectors, solve same ODE on every process
202:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203:   MatCreate(PETSC_COMM_WORLD,&A);
204:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
205:   MatSetFromOptions(A);
206:   MatSetUp(A);
207:   MatGetVecs(A,&x,NULL);

209:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210:      Create timestepping solver context
211:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212:   TSCreate(PETSC_COMM_WORLD,&ts);
213:   TSSetType(ts,TSBEULER);
214:   TSSetIFunction(ts,NULL,IFunction,&user);
215:   TSSetIJacobian(ts,A,A,IJacobian,&user);
216:   TSSetDuration(ts,PETSC_DEFAULT,ftime);
217:   if (monitor) {
218:     TSMonitorSet(ts,Monitor,&user,NULL);
219:   }

221:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222:      Set initial conditions
223:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224:   VecGetArray(x,&x_ptr);
225:   x_ptr[0] = -2;   x_ptr[1] = -2.355301397608119909925287735864250951918;
226:   VecRestoreArray(x,&x_ptr);
227:   TSSetInitialTimeStep(ts,0.0,.001);

229:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230:      Set runtime options
231:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
232:   TSSetFromOptions(ts);

234:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
235:      Solve nonlinear system
236:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
237:   TSSolve(ts,x);
238:   TSGetSolveTime(ts,&ftime);
239:   TSGetTimeStepNumber(ts,&steps);
240:   PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
241:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

243:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244:      Free work space.  All PETSc objects should be destroyed when they
245:      are no longer needed.
246:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
247:   MatDestroy(&A);
248:   VecDestroy(&x);
249:   TSDestroy(&ts);

251:   PetscFinalize();
252:   return(0);
253: }