xref: /petsc/src/ts/tutorials/multirate/ex3.c (revision f56ea12d3d7e20aed2e24c9530b477864b138a71)
1c4762a1bSJed Brown static char help[] = "Basic problem for multi-rate method.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*F
4c4762a1bSJed Brown 
5c4762a1bSJed Brown \begin{eqnarray}
6c4762a1bSJed Brown                  ys' = -2.0*\frac{-1.0+ys^2.0-\cos(t)}{2.0*ys}+0.05*\frac{-2.0+yf^2-\cos(5.0*t)}{2.0*yf}-\frac{\sin(t)}{2.0*ys}\\
7c4762a1bSJed Brown                  yf' = 0.05*\frac{-1.0+ys^2-\cos(t)}{2.0*ys}-\frac{-2.0+yf^2-\cos(5.0*t)}{2.0*yf}-5.0*\frac{\sin(5.0*t)}{2.0*yf}\\
8c4762a1bSJed Brown \end{eqnarray}
9c4762a1bSJed Brown 
10c4762a1bSJed Brown F*/
11c4762a1bSJed Brown 
12c4762a1bSJed Brown #include <petscts.h>
13c4762a1bSJed Brown 
14c4762a1bSJed Brown typedef struct {
15c4762a1bSJed Brown   PetscReal Tf,dt;
16c4762a1bSJed Brown } AppCtx;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
19c4762a1bSJed Brown {
20c4762a1bSJed Brown   PetscErrorCode    ierr;
21c4762a1bSJed Brown   const PetscScalar *u;
22c4762a1bSJed Brown   PetscScalar       *f;
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   PetscFunctionBegin;
25c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
26c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
27c4762a1bSJed Brown   f[0] = -2.0*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])+0.05*(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-PetscSinScalar(t)/(2.0*u[0]);
28c4762a1bSJed Brown   f[1] = 0.05*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])-(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-5.0*PetscSinScalar(5.0*t)/(2.0*u[1]);
29c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
31c4762a1bSJed Brown   PetscFunctionReturn(0);
32c4762a1bSJed Brown }
33c4762a1bSJed Brown 
34c4762a1bSJed Brown static PetscErrorCode RHSFunctionslow(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
35c4762a1bSJed Brown {
36c4762a1bSJed Brown   PetscErrorCode    ierr;
37c4762a1bSJed Brown   const PetscScalar *u;
38c4762a1bSJed Brown   PetscScalar       *f;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   PetscFunctionBegin;
41c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
43c4762a1bSJed Brown   f[0] = -2.0*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])+0.05*(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-PetscSinScalar(t)/(2.0*u[0]);
44c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
46c4762a1bSJed Brown   PetscFunctionReturn(0);
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown static PetscErrorCode RHSFunctionfast(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
50c4762a1bSJed Brown {
51c4762a1bSJed Brown   PetscErrorCode    ierr;
52c4762a1bSJed Brown   const PetscScalar *u;
53c4762a1bSJed Brown   PetscScalar       *f;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   PetscFunctionBegin;
56c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
57c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
58c4762a1bSJed Brown   f[0] = 0.05*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])-(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-5.0*PetscSinScalar(5.0*t)/(2.0*u[1]);
59c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
61c4762a1bSJed Brown   PetscFunctionReturn(0);
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown static PetscErrorCode sol_true(PetscReal t,Vec U)
65c4762a1bSJed Brown {
66c4762a1bSJed Brown   PetscErrorCode ierr;
67c4762a1bSJed Brown   PetscScalar    *u;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   PetscFunctionBegin;
70c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
71c4762a1bSJed Brown   u[0] = PetscSqrtScalar(1.0+PetscCosScalar(t));
72c4762a1bSJed Brown   u[1] = PetscSqrtScalar(2.0+PetscCosScalar(5.0*t));
73c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
74c4762a1bSJed Brown   PetscFunctionReturn(0);
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown int main(int argc,char **argv)
78c4762a1bSJed Brown {
79c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
80c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
81c4762a1bSJed Brown   Vec            Utrue;
82c4762a1bSJed Brown   PetscErrorCode ierr;
83c4762a1bSJed Brown   PetscMPIInt    size;
84c4762a1bSJed Brown   AppCtx         ctx;
85c4762a1bSJed Brown   PetscScalar    *u;
86c4762a1bSJed Brown   IS             iss;
87c4762a1bSJed Brown   IS             isf;
88c4762a1bSJed Brown   PetscInt       *indicess;
89c4762a1bSJed Brown   PetscInt       *indicesf;
90c4762a1bSJed Brown   PetscInt       n=2;
91c4762a1bSJed Brown   PetscReal      error,tt;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94c4762a1bSJed Brown      Initialize program
95c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
97c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
98c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
99c4762a1bSJed Brown 
100c4762a1bSJed Brown    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101c4762a1bSJed Brown     Create index for slow part and fast part
102c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103c4762a1bSJed Brown   ierr = PetscMalloc1(1,&indicess);CHKERRQ(ierr);
104c4762a1bSJed Brown   indicess[0]=0;
105c4762a1bSJed Brown   ierr = PetscMalloc1(1,&indicesf);CHKERRQ(ierr);
106c4762a1bSJed Brown   indicesf[0]=1;
107c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF,1,indicess,PETSC_COPY_VALUES,&iss);CHKERRQ(ierr);
108c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF,1,indicesf,PETSC_COPY_VALUES,&isf);CHKERRQ(ierr);
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111c4762a1bSJed Brown     Create necesary vector
112c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&U);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = VecSetSizes(U,n,PETSC_DETERMINE);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = VecSetFromOptions(U);CHKERRQ(ierr);
116c4762a1bSJed Brown   ierr = VecDuplicate(U,&Utrue);
117c4762a1bSJed Brown   ierr = VecCopy(U,Utrue);
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120c4762a1bSJed Brown     Set initial condition
121c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
123c4762a1bSJed Brown   u[0] = PetscSqrtScalar(2.0);
124c4762a1bSJed Brown   u[1] = PetscSqrtScalar(3.0);
125c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128c4762a1bSJed Brown      Create timestepping solver context
129c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
131c4762a1bSJed Brown   ierr = TSSetType(ts,TSMPRK);CHKERRQ(ierr);
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr);
134c4762a1bSJed Brown   ierr = TSRHSSplitSetIS(ts,"slow",iss);CHKERRQ(ierr);
135c4762a1bSJed Brown   ierr = TSRHSSplitSetIS(ts,"fast",isf);CHKERRQ(ierr);
136109dc152SHong Zhang   ierr = TSRHSSplitSetRHSFunction(ts,"slow",NULL,(TSRHSFunction)RHSFunctionslow,&ctx);CHKERRQ(ierr);
137109dc152SHong Zhang   ierr = TSRHSSplitSetRHSFunction(ts,"fast",NULL,(TSRHSFunction)RHSFunctionfast,&ctx);CHKERRQ(ierr);
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140c4762a1bSJed Brown      Set initial conditions
141c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142c4762a1bSJed Brown   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown      Set solver options
146c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ODE options","");CHKERRQ(ierr);
148c4762a1bSJed Brown   {
149c4762a1bSJed Brown     ctx.Tf = 0.3;
150c4762a1bSJed Brown     ctx.dt = 0.01;
151c4762a1bSJed Brown     ierr   = PetscOptionsScalar("-Tf","","",ctx.Tf,&ctx.Tf,NULL);CHKERRQ(ierr);
152c4762a1bSJed Brown     ierr   = PetscOptionsScalar("-dt","","",ctx.dt,&ctx.dt,NULL);CHKERRQ(ierr);
153c4762a1bSJed Brown   }
154c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
155c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ctx.Tf);CHKERRQ(ierr);
156c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,ctx.dt);CHKERRQ(ierr);
157c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
158c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161c4762a1bSJed Brown      Solve linear system
162c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
165c4762a1bSJed Brown 
166c4762a1bSJed Brown  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167c4762a1bSJed Brown      Check the error of the Petsc solution
168c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169c4762a1bSJed Brown   ierr = TSGetTime(ts,&tt);CHKERRQ(ierr);
170c4762a1bSJed Brown   ierr = sol_true(tt,Utrue);CHKERRQ(ierr);
171c4762a1bSJed Brown   ierr = VecAXPY(Utrue,-1.0,U);CHKERRQ(ierr);
172c4762a1bSJed Brown   ierr = VecNorm(Utrue,NORM_2,&error);
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175c4762a1bSJed Brown      Print norm2 error
176c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"l2 error norm: %g\n", error);CHKERRQ(ierr);
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
181c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
183c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
184c4762a1bSJed Brown   ierr = VecDestroy(&Utrue);CHKERRQ(ierr);
185c4762a1bSJed Brown   ierr = ISDestroy(&iss);CHKERRQ(ierr);
186c4762a1bSJed Brown   ierr = ISDestroy(&isf);CHKERRQ(ierr);
187c4762a1bSJed Brown   ierr = PetscFree(indicess);CHKERRQ(ierr);
188c4762a1bSJed Brown   ierr = PetscFree(indicesf);CHKERRQ(ierr);
189c4762a1bSJed Brown   ierr = PetscFinalize();
190c4762a1bSJed Brown   return ierr;
191c4762a1bSJed Brown }
192109dc152SHong Zhang 
193109dc152SHong Zhang /*TEST
194109dc152SHong Zhang     build:
195*f56ea12dSJed Brown       requires: !complex
196109dc152SHong Zhang 
197109dc152SHong Zhang     test:
198109dc152SHong Zhang 
199109dc152SHong Zhang TEST*/
200