xref: /petsc/src/ts/tutorials/ex49.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
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;
42c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
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]);
47c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
48c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
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;
62c4762a1bSJed Brown   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
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]);
65c4762a1bSJed Brown   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
66c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70c4762a1bSJed Brown   if (A != B) {
71c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
72c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
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      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
93ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
94c4762a1bSJed Brown   if (size != 1) SETERRQ(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;
1011e1ea65dSPierre Jolivet   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);CHKERRQ(ierr);
103c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
107c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
109c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
110c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
111c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr);
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116c4762a1bSJed Brown      Create timestepping solver context
117c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr);
121c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,A,A,IJacobian,&user);CHKERRQ(ierr);
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = TSSetTolerances(ts,1.e-4,NULL,1.e-4,NULL);CHKERRQ(ierr);
126c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown      Set initial conditions
128c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129c4762a1bSJed Brown   ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr);
130c4762a1bSJed Brown   x_ptr[0] = 2.0;   x_ptr[1] = -6.6e-01;
131c4762a1bSJed Brown   ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr);
132c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,.000001);CHKERRQ(ierr);
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135c4762a1bSJed Brown      Set runtime options
136c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140c4762a1bSJed Brown      Solve nonlinear system
141c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142c4762a1bSJed Brown   ierr = TSSolve(ts,x);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);CHKERRQ(ierr);
146c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
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    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
153c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
154c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   ierr = PetscFinalize();
157c4762a1bSJed Brown   return(ierr);
158c4762a1bSJed Brown }
159c4762a1bSJed Brown 
160c4762a1bSJed Brown /*TEST
161c4762a1bSJed Brown 
162c4762a1bSJed Brown     build:
163*dfd57a17SPierre 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