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