xref: /petsc/src/ts/tutorials/multirate/ex3fastslowsplit.c (revision 3a2a065f407e7da3bb6fc7eaf619dbab56034fde)
1*3a2a065fSHong Zhang static char help[] = "A fast-slow system for testing ARKIMEX.\n";
2*3a2a065fSHong Zhang 
3*3a2a065fSHong Zhang /*F
4*3a2a065fSHong Zhang 
5*3a2a065fSHong Zhang \begin{eqnarray}
6*3a2a065fSHong Zhang                  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*3a2a065fSHong Zhang                  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*3a2a065fSHong Zhang \end{eqnarray}
9*3a2a065fSHong Zhang 
10*3a2a065fSHong Zhang F*/
11*3a2a065fSHong Zhang 
12*3a2a065fSHong Zhang /*
13*3a2a065fSHong Zhang   This example demonstrates how to use ARKIMEX for solving a fast-slow system. The system is partitioned additively and component-wise at the same time.
14*3a2a065fSHong Zhang   ys stands for the slow component and yf stands for the fast component. On the RHS for yf, only the term -\frac{-2.0+yf^2-\cos(5.0*t)}{2.0*yf} is treated implicitly while the rest is treated explicitly.
15*3a2a065fSHong Zhang */
16*3a2a065fSHong Zhang 
17*3a2a065fSHong Zhang #include <petscts.h>
18*3a2a065fSHong Zhang 
19*3a2a065fSHong Zhang typedef struct {
20*3a2a065fSHong Zhang   PetscReal Tf, dt;
21*3a2a065fSHong Zhang } AppCtx;
22*3a2a065fSHong Zhang 
23*3a2a065fSHong Zhang static PetscErrorCode RHSFunctionslow(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
24*3a2a065fSHong Zhang {
25*3a2a065fSHong Zhang   const PetscScalar *u;
26*3a2a065fSHong Zhang   PetscScalar       *f;
27*3a2a065fSHong Zhang 
28*3a2a065fSHong Zhang   PetscFunctionBegin;
29*3a2a065fSHong Zhang   PetscCall(VecGetArrayRead(U, &u));
30*3a2a065fSHong Zhang   PetscCall(VecGetArray(F, &f));
31*3a2a065fSHong Zhang   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]);
32*3a2a065fSHong Zhang   PetscCall(VecRestoreArrayRead(U, &u));
33*3a2a065fSHong Zhang   PetscCall(VecRestoreArray(F, &f));
34*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
35*3a2a065fSHong Zhang }
36*3a2a065fSHong Zhang 
37*3a2a065fSHong Zhang static PetscErrorCode RHSFunctionfast(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
38*3a2a065fSHong Zhang {
39*3a2a065fSHong Zhang   const PetscScalar *u;
40*3a2a065fSHong Zhang   PetscScalar       *f;
41*3a2a065fSHong Zhang 
42*3a2a065fSHong Zhang   PetscFunctionBegin;
43*3a2a065fSHong Zhang   PetscCall(VecGetArrayRead(U, &u));
44*3a2a065fSHong Zhang   PetscCall(VecGetArray(F, &f));
45*3a2a065fSHong Zhang   f[0] = 0.05 * (-1.0 + u[0] * u[0] - PetscCosScalar(t)) / (2.0 * u[0]) - 5.0 * PetscSinScalar(5.0 * t) / (2.0 * u[1]);
46*3a2a065fSHong Zhang   PetscCall(VecRestoreArrayRead(U, &u));
47*3a2a065fSHong Zhang   PetscCall(VecRestoreArray(F, &f));
48*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
49*3a2a065fSHong Zhang }
50*3a2a065fSHong Zhang 
51*3a2a065fSHong Zhang static PetscErrorCode IFunctionfast(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
52*3a2a065fSHong Zhang {
53*3a2a065fSHong Zhang   const PetscScalar *u, *udot;
54*3a2a065fSHong Zhang   PetscScalar       *f;
55*3a2a065fSHong Zhang 
56*3a2a065fSHong Zhang   PetscFunctionBegin;
57*3a2a065fSHong Zhang   PetscCall(VecGetArrayRead(U, &u));
58*3a2a065fSHong Zhang   PetscCall(VecGetArrayRead(Udot, &udot));
59*3a2a065fSHong Zhang   PetscCall(VecGetArray(F, &f));
60*3a2a065fSHong Zhang   f[0] = udot[0] + (-2.0 + u[1] * u[1] - PetscCosScalar(5.0 * t)) / (2.0 * u[1]);
61*3a2a065fSHong Zhang   PetscCall(VecRestoreArrayRead(Udot, &udot));
62*3a2a065fSHong Zhang   PetscCall(VecRestoreArrayRead(U, &u));
63*3a2a065fSHong Zhang   PetscCall(VecRestoreArray(F, &f));
64*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
65*3a2a065fSHong Zhang }
66*3a2a065fSHong Zhang 
67*3a2a065fSHong Zhang static PetscErrorCode IJacobianfast(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx)
68*3a2a065fSHong Zhang {
69*3a2a065fSHong Zhang   PetscInt           rowcol[] = {0};
70*3a2a065fSHong Zhang   const PetscScalar *u;
71*3a2a065fSHong Zhang   PetscScalar        J[1][1];
72*3a2a065fSHong Zhang 
73*3a2a065fSHong Zhang   PetscFunctionBeginUser;
74*3a2a065fSHong Zhang   PetscCall(VecGetArrayRead(U, &u));
75*3a2a065fSHong Zhang   J[0][0] = a + (2.0 + PetscCosScalar(5.0 * t)) / (2.0 * u[1] * u[1]) + 0.5;
76*3a2a065fSHong Zhang   PetscCall(MatSetValues(B, 1, rowcol, 1, rowcol, &J[0][0], INSERT_VALUES));
77*3a2a065fSHong Zhang   PetscCall(VecRestoreArrayRead(U, &u));
78*3a2a065fSHong Zhang 
79*3a2a065fSHong Zhang   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
80*3a2a065fSHong Zhang   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
81*3a2a065fSHong Zhang   if (A != B) {
82*3a2a065fSHong Zhang     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
83*3a2a065fSHong Zhang     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
84*3a2a065fSHong Zhang   }
85*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
86*3a2a065fSHong Zhang }
87*3a2a065fSHong Zhang 
88*3a2a065fSHong Zhang /*
89*3a2a065fSHong Zhang   Define the analytic solution for check method easily
90*3a2a065fSHong Zhang */
91*3a2a065fSHong Zhang static PetscErrorCode sol_true(PetscReal t, Vec U)
92*3a2a065fSHong Zhang {
93*3a2a065fSHong Zhang   PetscScalar *u;
94*3a2a065fSHong Zhang 
95*3a2a065fSHong Zhang   PetscFunctionBegin;
96*3a2a065fSHong Zhang   PetscCall(VecGetArray(U, &u));
97*3a2a065fSHong Zhang   u[0] = PetscSqrtScalar(1.0 + PetscCosScalar(t));
98*3a2a065fSHong Zhang   u[1] = PetscSqrtScalar(2.0 + PetscCosScalar(5.0 * t));
99*3a2a065fSHong Zhang   PetscCall(VecRestoreArray(U, &u));
100*3a2a065fSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
101*3a2a065fSHong Zhang }
102*3a2a065fSHong Zhang 
103*3a2a065fSHong Zhang int main(int argc, char **argv)
104*3a2a065fSHong Zhang {
105*3a2a065fSHong Zhang   TS           ts; /* ODE integrator */
106*3a2a065fSHong Zhang   Vec          U;  /* solution will be stored here */
107*3a2a065fSHong Zhang   Vec          Utrue;
108*3a2a065fSHong Zhang   Mat          A;
109*3a2a065fSHong Zhang   PetscMPIInt  size;
110*3a2a065fSHong Zhang   AppCtx       ctx;
111*3a2a065fSHong Zhang   PetscScalar *u;
112*3a2a065fSHong Zhang   IS           iss;
113*3a2a065fSHong Zhang   IS           isf;
114*3a2a065fSHong Zhang   PetscInt    *indicess;
115*3a2a065fSHong Zhang   PetscInt    *indicesf;
116*3a2a065fSHong Zhang   PetscInt     n = 2;
117*3a2a065fSHong Zhang   PetscScalar  error, tt;
118*3a2a065fSHong Zhang 
119*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120*3a2a065fSHong Zhang      Initialize program
121*3a2a065fSHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122*3a2a065fSHong Zhang   PetscFunctionBeginUser;
123*3a2a065fSHong Zhang   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
124*3a2a065fSHong Zhang   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
125*3a2a065fSHong Zhang   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
126*3a2a065fSHong Zhang 
127*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128*3a2a065fSHong Zhang     Create index for slow part and fast part
129*3a2a065fSHong Zhang     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130*3a2a065fSHong Zhang   PetscCall(PetscMalloc1(1, &indicess));
131*3a2a065fSHong Zhang   indicess[0] = 0;
132*3a2a065fSHong Zhang   PetscCall(PetscMalloc1(1, &indicesf));
133*3a2a065fSHong Zhang   indicesf[0] = 1;
134*3a2a065fSHong Zhang   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, indicess, PETSC_COPY_VALUES, &iss));
135*3a2a065fSHong Zhang   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, indicesf, PETSC_COPY_VALUES, &isf));
136*3a2a065fSHong Zhang 
137*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138*3a2a065fSHong Zhang     Create necessary vector
139*3a2a065fSHong Zhang     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140*3a2a065fSHong Zhang   PetscCall(VecCreate(PETSC_COMM_WORLD, &U));
141*3a2a065fSHong Zhang   PetscCall(VecSetSizes(U, n, PETSC_DETERMINE));
142*3a2a065fSHong Zhang   PetscCall(VecSetFromOptions(U));
143*3a2a065fSHong Zhang   PetscCall(VecDuplicate(U, &Utrue));
144*3a2a065fSHong Zhang   PetscCall(VecCopy(U, Utrue));
145*3a2a065fSHong Zhang 
146*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147*3a2a065fSHong Zhang     Set runtime options
148*3a2a065fSHong Zhang     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149*3a2a065fSHong Zhang   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ODE options", "");
150*3a2a065fSHong Zhang   ctx.Tf = 0.3;
151*3a2a065fSHong Zhang   ctx.dt = 0.01;
152*3a2a065fSHong Zhang   PetscCall(PetscOptionsScalar("-Tf", "", "", ctx.Tf, &ctx.Tf, NULL));
153*3a2a065fSHong Zhang   PetscCall(PetscOptionsScalar("-dt", "", "", ctx.dt, &ctx.dt, NULL));
154*3a2a065fSHong Zhang   PetscOptionsEnd();
155*3a2a065fSHong Zhang 
156*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157*3a2a065fSHong Zhang      Initialize the solution
158*3a2a065fSHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159*3a2a065fSHong Zhang   PetscCall(VecGetArray(U, &u));
160*3a2a065fSHong Zhang   u[0] = PetscSqrtScalar(2.0);
161*3a2a065fSHong Zhang   u[1] = PetscSqrtScalar(3.0);
162*3a2a065fSHong Zhang   PetscCall(VecRestoreArray(U, &u));
163*3a2a065fSHong Zhang 
164*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165*3a2a065fSHong Zhang      Create timestepping solver context
166*3a2a065fSHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167*3a2a065fSHong Zhang   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
168*3a2a065fSHong Zhang   PetscCall(TSSetType(ts, TSARKIMEX));
169*3a2a065fSHong Zhang   PetscCall(TSARKIMEXSetFastSlowSplit(ts, PETSC_TRUE));
170*3a2a065fSHong Zhang 
171*3a2a065fSHong Zhang   PetscCall(TSRHSSplitSetIS(ts, "slow", iss));
172*3a2a065fSHong Zhang   PetscCall(TSRHSSplitSetIS(ts, "fast", isf));
173*3a2a065fSHong Zhang   PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, (TSRHSFunctionFn *)RHSFunctionslow, &ctx));
174*3a2a065fSHong Zhang   PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, (TSRHSFunctionFn *)RHSFunctionfast, &ctx));
175*3a2a065fSHong Zhang   PetscCall(TSRHSSplitSetIFunction(ts, "fast", NULL, (TSIFunctionFn *)IFunctionfast, &ctx));
176*3a2a065fSHong Zhang   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
177*3a2a065fSHong Zhang   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
178*3a2a065fSHong Zhang   PetscCall(MatSetFromOptions(A));
179*3a2a065fSHong Zhang   PetscCall(MatSetUp(A));
180*3a2a065fSHong Zhang   PetscCall(TSRHSSplitSetIJacobian(ts, "fast", A, A, (TSIJacobianFn *)IJacobianfast, &ctx));
181*3a2a065fSHong Zhang 
182*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183*3a2a065fSHong Zhang      Set initial conditions
184*3a2a065fSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185*3a2a065fSHong Zhang   PetscCall(TSSetSolution(ts, U));
186*3a2a065fSHong Zhang 
187*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188*3a2a065fSHong Zhang      Set solver options
189*3a2a065fSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
190*3a2a065fSHong Zhang   PetscCall(TSSetMaxTime(ts, ctx.Tf));
191*3a2a065fSHong Zhang   PetscCall(TSSetTimeStep(ts, ctx.dt));
192*3a2a065fSHong Zhang   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
193*3a2a065fSHong Zhang   PetscCall(TSSetFromOptions(ts));
194*3a2a065fSHong Zhang 
195*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196*3a2a065fSHong Zhang      Solve linear system
197*3a2a065fSHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198*3a2a065fSHong Zhang   PetscCall(TSSolve(ts, U));
199*3a2a065fSHong Zhang   PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
200*3a2a065fSHong Zhang 
201*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202*3a2a065fSHong Zhang      Check the error of the Petsc solution
203*3a2a065fSHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
204*3a2a065fSHong Zhang   PetscCall(TSGetTime(ts, &tt));
205*3a2a065fSHong Zhang   PetscCall(sol_true(tt, Utrue));
206*3a2a065fSHong Zhang   PetscCall(VecAXPY(Utrue, -1.0, U));
207*3a2a065fSHong Zhang   PetscCall(VecNorm(Utrue, NORM_2, &error));
208*3a2a065fSHong Zhang 
209*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210*3a2a065fSHong Zhang      Print norm2 error
211*3a2a065fSHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212*3a2a065fSHong Zhang   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error norm: %g\n", (double)PetscAbsScalar(error)));
213*3a2a065fSHong Zhang 
214*3a2a065fSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215*3a2a065fSHong Zhang      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
216*3a2a065fSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217*3a2a065fSHong Zhang   PetscCall(VecDestroy(&U));
218*3a2a065fSHong Zhang   PetscCall(TSDestroy(&ts));
219*3a2a065fSHong Zhang   PetscCall(VecDestroy(&Utrue));
220*3a2a065fSHong Zhang   PetscCall(ISDestroy(&iss));
221*3a2a065fSHong Zhang   PetscCall(ISDestroy(&isf));
222*3a2a065fSHong Zhang   PetscCall(PetscFree(indicess));
223*3a2a065fSHong Zhang   PetscCall(PetscFree(indicesf));
224*3a2a065fSHong Zhang   PetscCall(MatDestroy(&A));
225*3a2a065fSHong Zhang   PetscCall(PetscFinalize());
226*3a2a065fSHong Zhang   return 0;
227*3a2a065fSHong Zhang }
228*3a2a065fSHong Zhang 
229*3a2a065fSHong Zhang /*TEST
230*3a2a065fSHong Zhang     build:
231*3a2a065fSHong Zhang       requires: !complex
232*3a2a065fSHong Zhang 
233*3a2a065fSHong Zhang     test:
234*3a2a065fSHong Zhang 
235*3a2a065fSHong Zhang     test:
236*3a2a065fSHong Zhang       suffix: 2
237*3a2a065fSHong Zhang       args: -ts_arkimex_initial_guess_extrapolate 1
238*3a2a065fSHong Zhang       output_file: output/ex3fastslowsplit_1.out
239*3a2a065fSHong Zhang 
240*3a2a065fSHong Zhang TEST*/
241