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