xref: /petsc/src/ts/tutorials/ex51.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[] = "Small ODE to test TS accuracy.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   The ODE
5c4762a1bSJed Brown                   u1_t = cos(t),
6c4762a1bSJed Brown                   u2_t = sin(u2)
7c4762a1bSJed Brown   with analytical solution
8c4762a1bSJed Brown                   u1(t) = sin(t),
9c4762a1bSJed Brown                   u2(t) = 2 * atan(exp(t) * tan(0.5))
10c4762a1bSJed Brown   is used to test the accuracy of TS schemes.
11c4762a1bSJed Brown */
12c4762a1bSJed Brown 
13c4762a1bSJed Brown #include <petscts.h>
14c4762a1bSJed Brown 
15c4762a1bSJed Brown /*
16c4762a1bSJed Brown      Defines the ODE passed to the ODE solver in explicit form: U_t = F(U)
17c4762a1bSJed Brown */
18c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *s)
19c4762a1bSJed Brown {
20c4762a1bSJed Brown   PetscScalar       *f;
21c4762a1bSJed Brown   const PetscScalar *u;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   PetscFunctionBegin;
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   f[0] = PetscCosReal(t);
28c4762a1bSJed Brown   f[1] = PetscSinReal(u[1]);
29c4762a1bSJed Brown 
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
32c4762a1bSJed Brown   PetscFunctionReturn(0);
33c4762a1bSJed Brown }
34c4762a1bSJed Brown 
35c4762a1bSJed Brown /*
36c4762a1bSJed Brown      Defines the exact solution.
37c4762a1bSJed Brown */
38c4762a1bSJed Brown static PetscErrorCode ExactSolution(PetscReal t, Vec U)
39c4762a1bSJed Brown {
40c4762a1bSJed Brown   PetscScalar       *u;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   PetscFunctionBegin;
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(U,&u));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   u[0] = PetscSinReal(t);
46c4762a1bSJed Brown   u[1] = 2 * PetscAtanReal(PetscExpReal(t) * PetscTanReal(0.5));
47c4762a1bSJed Brown 
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(U,&u));
49c4762a1bSJed Brown   PetscFunctionReturn(0);
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
52c4762a1bSJed Brown int main(int argc,char **argv)
53c4762a1bSJed Brown {
54c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
55c4762a1bSJed Brown   Vec            U;             /* numerical solution will be stored here */
56c4762a1bSJed Brown   Vec            Uex;           /* analytical (exact) solution will be stored here */
57c4762a1bSJed Brown   PetscErrorCode ierr;
58c4762a1bSJed Brown   PetscMPIInt    size;
59c4762a1bSJed Brown   PetscInt       n = 2;
60c4762a1bSJed Brown   PetscScalar    *u;
61c4762a1bSJed Brown   PetscReal      t, final_time = 1.0, dt = 0.25;
62c4762a1bSJed Brown   PetscReal      error;
63c4762a1bSJed Brown   TSAdapt        adapt;
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66c4762a1bSJed Brown      Initialize program
67c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
69*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
703c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73c4762a1bSJed Brown      Create timestepping solver context
74c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSROSW));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79c4762a1bSJed Brown      Set ODE routines
80c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,NULL));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85c4762a1bSJed Brown      Set initial conditions
86c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&U));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(U,n,PETSC_DETERMINE));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetUp(U));
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(U,&u));
91c4762a1bSJed Brown   u[0] = 0.0;
92c4762a1bSJed Brown   u[1] = 1.0;
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(U,&u));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,U));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97c4762a1bSJed Brown      Set solver options
98c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(ts));
100*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,final_time));
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
103a5b23f4aSJose E. Roman   /* The adaptive time step controller is forced to take constant time steps. */
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetAdapt(ts,&adapt));
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdaptSetType(adapt,TSADAPTNONE));
106c4762a1bSJed Brown 
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110c4762a1bSJed Brown      Run timestepping solver and compute error
111c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,U));
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts,&t));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   if (PetscAbsReal(t-final_time)>100*PETSC_MACHINE_EPSILON) {
116*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Note: There is a difference of %g between the prescribed final time %g and the actual final time.\n",(double)(final_time-t),(double)final_time));
117c4762a1bSJed Brown   }
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(U,&Uex));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ExactSolution(t,Uex));
120c4762a1bSJed Brown 
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAYPX(Uex,-1.0,U));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(Uex,NORM_2,&error));
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error at final time: %.2E\n",(double)error));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
127c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Uex));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   ierr = PetscFinalize();
133c4762a1bSJed Brown   return ierr;
134c4762a1bSJed Brown }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown /*TEST
137c4762a1bSJed Brown 
138c4762a1bSJed Brown     test:
139c4762a1bSJed Brown       suffix: 3bs
140c4762a1bSJed Brown       args: -ts_type rk -ts_rk_type 3bs
141c4762a1bSJed Brown       requires: !single
142c4762a1bSJed Brown 
143c4762a1bSJed Brown     test:
144c4762a1bSJed Brown       suffix: 5bs
145c4762a1bSJed Brown       args: -ts_type rk -ts_rk_type 5bs
146c4762a1bSJed Brown       requires: !single
147c4762a1bSJed Brown 
148c4762a1bSJed Brown     test:
149c4762a1bSJed Brown       suffix: 5dp
150c4762a1bSJed Brown       args: -ts_type rk -ts_rk_type 5dp
151c4762a1bSJed Brown       requires: !single
152c4762a1bSJed Brown 
153c4762a1bSJed Brown     test:
154c4762a1bSJed Brown       suffix: 6vr
155c4762a1bSJed Brown       args: -ts_type rk -ts_rk_type 6vr
156c4762a1bSJed Brown       requires: !single
157c4762a1bSJed Brown 
158c4762a1bSJed Brown     test:
159c4762a1bSJed Brown       suffix: 7vr
160c4762a1bSJed Brown       args: -ts_type rk -ts_rk_type 7vr
161c4762a1bSJed Brown       requires: !single
162c4762a1bSJed Brown 
163c4762a1bSJed Brown     test:
164c4762a1bSJed Brown       suffix: 8vr
165c4762a1bSJed Brown       args: -ts_type rk -ts_rk_type 8vr
166c4762a1bSJed Brown       requires: !single
167c4762a1bSJed Brown 
168c4762a1bSJed Brown TEST*/
169