1fa9584fbSIlya Fursov static char help[] = "Solves the trivial ODE du/dt = 1, u(0) = 0. \n\n";
2fa9584fbSIlya Fursov
3fa9584fbSIlya Fursov #include <petscts.h>
4fa9584fbSIlya Fursov #include <petscpc.h>
5fa9584fbSIlya Fursov
6fa9584fbSIlya Fursov static PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
7fa9584fbSIlya Fursov static PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
8fa9584fbSIlya Fursov
9fa9584fbSIlya Fursov static PetscErrorCode PreStep(TS);
10fa9584fbSIlya Fursov static PetscErrorCode PostStep(TS);
11fa9584fbSIlya Fursov static PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
12fa9584fbSIlya Fursov static PetscErrorCode Event(TS, PetscReal, Vec, PetscReal *, void *);
13fa9584fbSIlya Fursov static PetscErrorCode PostEvent(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *);
14fa9584fbSIlya Fursov static PetscErrorCode TransferSetUp(TS, PetscInt, PetscReal, Vec, PetscBool *, void *);
15fa9584fbSIlya Fursov static PetscErrorCode Transfer(TS, PetscInt, Vec[], Vec[], void *);
16fa9584fbSIlya Fursov
main(int argc,char ** argv)17fa9584fbSIlya Fursov int main(int argc, char **argv)
18fa9584fbSIlya Fursov {
19fa9584fbSIlya Fursov TS ts;
20ecc87898SStefano Zampini PetscInt n, ntransfer[] = {2, 2};
21fa9584fbSIlya Fursov const PetscInt n_end = 11;
22fa9584fbSIlya Fursov PetscReal t;
23fa9584fbSIlya Fursov const PetscReal t_end = 11;
24fa9584fbSIlya Fursov Vec x;
25fa9584fbSIlya Fursov Vec f;
26fa9584fbSIlya Fursov Mat A;
27fa9584fbSIlya Fursov
28fa9584fbSIlya Fursov PetscFunctionBeginUser;
29c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
30fa9584fbSIlya Fursov
31fa9584fbSIlya Fursov PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
32fa9584fbSIlya Fursov
33fa9584fbSIlya Fursov PetscCall(VecCreate(PETSC_COMM_WORLD, &f));
34fa9584fbSIlya Fursov PetscCall(VecSetSizes(f, 1, PETSC_DECIDE));
35fa9584fbSIlya Fursov PetscCall(VecSetFromOptions(f));
36fa9584fbSIlya Fursov PetscCall(VecSetUp(f));
37fa9584fbSIlya Fursov PetscCall(TSSetRHSFunction(ts, f, RHSFunction, NULL));
38fa9584fbSIlya Fursov PetscCall(VecDestroy(&f));
39fa9584fbSIlya Fursov
40fa9584fbSIlya Fursov PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
41fa9584fbSIlya Fursov PetscCall(MatSetSizes(A, 1, 1, PETSC_DECIDE, PETSC_DECIDE));
42fa9584fbSIlya Fursov PetscCall(MatSetFromOptions(A));
43fa9584fbSIlya Fursov PetscCall(MatSetUp(A));
44fa9584fbSIlya Fursov PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
45fa9584fbSIlya Fursov PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
46fa9584fbSIlya Fursov /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */
47fa9584fbSIlya Fursov PetscCall(MatShift(A, (PetscReal)1));
48fa9584fbSIlya Fursov PetscCall(MatShift(A, (PetscReal)-1));
49fa9584fbSIlya Fursov PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL));
50fa9584fbSIlya Fursov PetscCall(MatDestroy(&A));
51fa9584fbSIlya Fursov
52fa9584fbSIlya Fursov PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
53fa9584fbSIlya Fursov PetscCall(VecSetSizes(x, 1, PETSC_DECIDE));
54fa9584fbSIlya Fursov PetscCall(VecSetFromOptions(x));
55fa9584fbSIlya Fursov PetscCall(VecSetUp(x));
56fa9584fbSIlya Fursov PetscCall(TSSetSolution(ts, x));
57fa9584fbSIlya Fursov PetscCall(VecDestroy(&x));
58fa9584fbSIlya Fursov
59fa9584fbSIlya Fursov PetscCall(TSMonitorSet(ts, Monitor, NULL, NULL));
60fa9584fbSIlya Fursov PetscCall(TSSetPreStep(ts, PreStep));
61fa9584fbSIlya Fursov PetscCall(TSSetPostStep(ts, PostStep));
62fa9584fbSIlya Fursov
63fa9584fbSIlya Fursov {
64fa9584fbSIlya Fursov TSAdapt adapt;
65fa9584fbSIlya Fursov PetscCall(TSGetAdapt(ts, &adapt));
66fa9584fbSIlya Fursov PetscCall(TSAdaptSetType(adapt, TSADAPTNONE));
67fa9584fbSIlya Fursov }
68fa9584fbSIlya Fursov {
69fa9584fbSIlya Fursov PetscInt direction[3];
70fa9584fbSIlya Fursov PetscBool terminate[3];
71fa9584fbSIlya Fursov direction[0] = +1;
72fa9584fbSIlya Fursov terminate[0] = PETSC_FALSE;
73fa9584fbSIlya Fursov direction[1] = -1;
74fa9584fbSIlya Fursov terminate[1] = PETSC_FALSE;
75fa9584fbSIlya Fursov direction[2] = 0;
76fa9584fbSIlya Fursov terminate[2] = PETSC_FALSE;
77fa9584fbSIlya Fursov PetscCall(TSSetTimeStep(ts, 1));
78fa9584fbSIlya Fursov PetscCall(TSSetEventHandler(ts, 3, direction, terminate, Event, PostEvent, NULL));
79fa9584fbSIlya Fursov }
80fa9584fbSIlya Fursov PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
81ecc87898SStefano Zampini PetscCall(TSSetResize(ts, PETSC_TRUE, TransferSetUp, Transfer, ntransfer));
82fa9584fbSIlya Fursov PetscCall(TSSetFromOptions(ts));
83fa9584fbSIlya Fursov
84fa9584fbSIlya Fursov /* --- First Solve --- */
85fa9584fbSIlya Fursov
86fa9584fbSIlya Fursov PetscCall(TSSetStepNumber(ts, 0));
87fa9584fbSIlya Fursov PetscCall(TSSetTimeStep(ts, 1));
88fa9584fbSIlya Fursov PetscCall(TSSetTime(ts, 0));
89fa9584fbSIlya Fursov PetscCall(TSSetMaxTime(ts, PETSC_MAX_REAL));
90fa9584fbSIlya Fursov PetscCall(TSSetMaxSteps(ts, 3));
91fa9584fbSIlya Fursov
92fa9584fbSIlya Fursov PetscCall(TSGetTime(ts, &t));
93fa9584fbSIlya Fursov PetscCall(TSGetSolution(ts, &x));
94fa9584fbSIlya Fursov PetscCall(VecSet(x, t));
95fa9584fbSIlya Fursov while (t < t_end) {
96fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "TSSolve: Begin\n"));
97fa9584fbSIlya Fursov PetscCall(TSSolve(ts, NULL));
98fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "TSSolve: End\n\n"));
99fa9584fbSIlya Fursov PetscCall(TSGetTime(ts, &t));
100fa9584fbSIlya Fursov PetscCall(TSGetStepNumber(ts, &n));
101fa9584fbSIlya Fursov PetscCall(TSSetMaxSteps(ts, PetscMin(n + 3, n_end)));
102fa9584fbSIlya Fursov }
103fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "TSSolve: Begin\n"));
104fa9584fbSIlya Fursov PetscCall(TSSolve(ts, NULL));
105fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "TSSolve: End\n\n"));
106fa9584fbSIlya Fursov
107fa9584fbSIlya Fursov /* --- Second Solve --- */
108fa9584fbSIlya Fursov
109fa9584fbSIlya Fursov PetscCall(TSSetStepNumber(ts, 0));
110fa9584fbSIlya Fursov PetscCall(TSSetTimeStep(ts, 1));
111fa9584fbSIlya Fursov PetscCall(TSSetTime(ts, 0));
112fa9584fbSIlya Fursov PetscCall(TSSetMaxTime(ts, 3));
1131690c2aeSBarry Smith PetscCall(TSSetMaxSteps(ts, PETSC_INT_MAX));
114fa9584fbSIlya Fursov
115fa9584fbSIlya Fursov PetscCall(TSGetTime(ts, &t));
116fa9584fbSIlya Fursov PetscCall(TSGetSolution(ts, &x));
117fa9584fbSIlya Fursov PetscCall(VecSet(x, t));
118fa9584fbSIlya Fursov while (t < t_end) {
119fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "TSSolve: Begin\n"));
120fa9584fbSIlya Fursov PetscCall(TSSolve(ts, NULL));
121fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "TSSolve: End\n\n"));
122fa9584fbSIlya Fursov PetscCall(TSGetTime(ts, &t));
123fa9584fbSIlya Fursov PetscCall(TSSetMaxTime(ts, PetscMin(t + 3, t_end)));
124fa9584fbSIlya Fursov }
125fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "TSSolve: Begin\n"));
126fa9584fbSIlya Fursov PetscCall(TSSolve(ts, NULL));
127fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "TSSolve: End\n\n"));
128fa9584fbSIlya Fursov
129fa9584fbSIlya Fursov /* --- */
130fa9584fbSIlya Fursov
131fa9584fbSIlya Fursov PetscCall(TSDestroy(&ts));
132fa9584fbSIlya Fursov
133fa9584fbSIlya Fursov PetscCall(PetscFinalize());
134fa9584fbSIlya Fursov return 0;
135fa9584fbSIlya Fursov }
136fa9584fbSIlya Fursov
137fa9584fbSIlya Fursov /* -------------------------------------------------------------------*/
138fa9584fbSIlya Fursov
RHSFunction(TS ts,PetscReal t,Vec x,Vec f,PetscCtx ctx)139*2a8381b2SBarry Smith PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec x, Vec f, PetscCtx ctx)
140fa9584fbSIlya Fursov {
141fa9584fbSIlya Fursov PetscFunctionBeginUser;
142fa9584fbSIlya Fursov PetscCall(VecSet(f, (PetscReal)1));
143fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
144fa9584fbSIlya Fursov }
145fa9584fbSIlya Fursov
RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat B,PetscCtx ctx)146*2a8381b2SBarry Smith PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat B, PetscCtx ctx)
147fa9584fbSIlya Fursov {
148fa9584fbSIlya Fursov PetscFunctionBeginUser;
149fa9584fbSIlya Fursov PetscCall(MatZeroEntries(B));
150fa9584fbSIlya Fursov if (B != A) PetscCall(MatZeroEntries(A));
151fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
152fa9584fbSIlya Fursov }
153fa9584fbSIlya Fursov
PreStep(TS ts)154fa9584fbSIlya Fursov PetscErrorCode PreStep(TS ts)
155fa9584fbSIlya Fursov {
156fa9584fbSIlya Fursov PetscInt n;
157fa9584fbSIlya Fursov PetscReal t;
158fa9584fbSIlya Fursov Vec x;
159fa9584fbSIlya Fursov const PetscScalar *a;
160c6bf8827SStefano Zampini PetscBool flg;
161fa9584fbSIlya Fursov
162fa9584fbSIlya Fursov PetscFunctionBeginUser;
163fa9584fbSIlya Fursov PetscCall(TSGetStepNumber(ts, &n));
164fa9584fbSIlya Fursov PetscCall(TSGetTime(ts, &t));
165fa9584fbSIlya Fursov PetscCall(TSGetSolution(ts, &x));
166fa9584fbSIlya Fursov PetscCall(VecGetArrayRead(x, &a));
167c6bf8827SStefano Zampini PetscCall(TSGetStepResize(ts, &flg));
168c6bf8827SStefano Zampini PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%-10s-> step %" PetscInt_FMT " time %g value %g%s\n", PETSC_FUNCTION_NAME, n, (double)t, (double)PetscRealPart(a[0]), flg ? " resized" : ""));
169fa9584fbSIlya Fursov PetscCall(VecRestoreArrayRead(x, &a));
170fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
171fa9584fbSIlya Fursov }
172fa9584fbSIlya Fursov
PostStep(TS ts)173fa9584fbSIlya Fursov PetscErrorCode PostStep(TS ts)
174fa9584fbSIlya Fursov {
175fa9584fbSIlya Fursov PetscInt n;
176fa9584fbSIlya Fursov PetscReal t;
177fa9584fbSIlya Fursov Vec x;
178fa9584fbSIlya Fursov const PetscScalar *a;
179fa9584fbSIlya Fursov
180fa9584fbSIlya Fursov PetscFunctionBeginUser;
181fa9584fbSIlya Fursov PetscCall(TSGetStepNumber(ts, &n));
182fa9584fbSIlya Fursov PetscCall(TSGetTime(ts, &t));
183fa9584fbSIlya Fursov PetscCall(TSGetSolution(ts, &x));
184fa9584fbSIlya Fursov PetscCall(VecGetArrayRead(x, &a));
185fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%-10s-> step %" PetscInt_FMT " time %g value %g\n", PETSC_FUNCTION_NAME, n, (double)t, (double)PetscRealPart(a[0])));
186fa9584fbSIlya Fursov PetscCall(VecRestoreArrayRead(x, &a));
187fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
188fa9584fbSIlya Fursov }
189fa9584fbSIlya Fursov
Monitor(TS ts,PetscInt n,PetscReal t,Vec x,PetscCtx ctx)190*2a8381b2SBarry Smith PetscErrorCode Monitor(TS ts, PetscInt n, PetscReal t, Vec x, PetscCtx ctx)
191fa9584fbSIlya Fursov {
192fa9584fbSIlya Fursov const PetscScalar *a;
193fa9584fbSIlya Fursov
194fa9584fbSIlya Fursov PetscFunctionBeginUser;
195fa9584fbSIlya Fursov PetscCall(VecGetArrayRead(x, &a));
196fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%-10s-> step %" PetscInt_FMT " time %g value %g\n", PETSC_FUNCTION_NAME, n, (double)t, (double)PetscRealPart(a[0])));
197fa9584fbSIlya Fursov PetscCall(VecRestoreArrayRead(x, &a));
198fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
199fa9584fbSIlya Fursov }
200fa9584fbSIlya Fursov
Event(TS ts,PetscReal t,Vec x,PetscReal * fvalue,PetscCtx ctx)201*2a8381b2SBarry Smith PetscErrorCode Event(TS ts, PetscReal t, Vec x, PetscReal *fvalue, PetscCtx ctx)
202fa9584fbSIlya Fursov {
203fa9584fbSIlya Fursov PetscFunctionBeginUser;
204fa9584fbSIlya Fursov fvalue[0] = t - 5;
205fa9584fbSIlya Fursov fvalue[1] = 7 - t;
206fa9584fbSIlya Fursov fvalue[2] = t - 9;
207fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
208fa9584fbSIlya Fursov }
209fa9584fbSIlya Fursov
PostEvent(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec x,PetscBool forwardsolve,PetscCtx ctx)210*2a8381b2SBarry Smith PetscErrorCode PostEvent(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec x, PetscBool forwardsolve, PetscCtx ctx)
211fa9584fbSIlya Fursov {
212fa9584fbSIlya Fursov PetscInt i;
213fa9584fbSIlya Fursov const PetscScalar *a;
214fa9584fbSIlya Fursov
215fa9584fbSIlya Fursov PetscFunctionBeginUser;
216fa9584fbSIlya Fursov PetscCall(TSGetStepNumber(ts, &i));
217fa9584fbSIlya Fursov PetscCall(VecGetArrayRead(x, &a));
218fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%-10s-> step %" PetscInt_FMT " time %g value %g\n", PETSC_FUNCTION_NAME, i, (double)t, (double)PetscRealPart(a[0])));
219fa9584fbSIlya Fursov PetscCall(VecRestoreArrayRead(x, &a));
220fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
221fa9584fbSIlya Fursov }
222fa9584fbSIlya Fursov
TransferSetUp(TS ts,PetscInt n,PetscReal t,Vec x,PetscBool * flg,PetscCtx ctx)223*2a8381b2SBarry Smith PetscErrorCode TransferSetUp(TS ts, PetscInt n, PetscReal t, Vec x, PetscBool *flg, PetscCtx ctx)
224fa9584fbSIlya Fursov {
225fa9584fbSIlya Fursov PetscInt *nt = (PetscInt *)ctx;
226fa9584fbSIlya Fursov
227fa9584fbSIlya Fursov PetscFunctionBeginUser;
228ecc87898SStefano Zampini if (n == 1) {
229ecc87898SStefano Zampini nt[0] = 2;
230ecc87898SStefano Zampini nt[1] = 2;
231ecc87898SStefano Zampini }
232ecc87898SStefano Zampini *flg = (PetscBool)(nt[0] && n && n % (nt[0]) == 0);
233ecc87898SStefano Zampini if (*flg) nt[0] += nt[1];
234fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
235fa9584fbSIlya Fursov }
236fa9584fbSIlya Fursov
Transfer(TS ts,PetscInt nv,Vec vin[],Vec vout[],PetscCtx ctx)237*2a8381b2SBarry Smith PetscErrorCode Transfer(TS ts, PetscInt nv, Vec vin[], Vec vout[], PetscCtx ctx)
238fa9584fbSIlya Fursov {
239fa9584fbSIlya Fursov PetscInt i;
240fa9584fbSIlya Fursov
241fa9584fbSIlya Fursov PetscFunctionBeginUser;
242fa9584fbSIlya Fursov PetscCall(TSGetStepNumber(ts, &i));
243fa9584fbSIlya Fursov PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%-10s-> step %" PetscInt_FMT " nv %" PetscInt_FMT "\n", PETSC_FUNCTION_NAME, i, nv));
244fa9584fbSIlya Fursov for (i = 0; i < nv; i++) {
245fa9584fbSIlya Fursov PetscCall(VecDuplicate(vin[i], &vout[i]));
246fa9584fbSIlya Fursov PetscCall(VecCopy(vin[i], vout[i]));
247fa9584fbSIlya Fursov }
248fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS);
249fa9584fbSIlya Fursov }
250fa9584fbSIlya Fursov
251fa9584fbSIlya Fursov /*TEST
252fa9584fbSIlya Fursov
253fa9584fbSIlya Fursov test:
254fa9584fbSIlya Fursov suffix: euler
255fa9584fbSIlya Fursov diff_args: -j
256fa9584fbSIlya Fursov args: -ts_type euler
257fa9584fbSIlya Fursov output_file: output/ex1.out
258fa9584fbSIlya Fursov
259fa9584fbSIlya Fursov test:
260fa9584fbSIlya Fursov suffix: ssp
261fa9584fbSIlya Fursov diff_args: -j
262fa9584fbSIlya Fursov args: -ts_type ssp
263fa9584fbSIlya Fursov output_file: output/ex1.out
264fa9584fbSIlya Fursov
265fa9584fbSIlya Fursov test:
266fa9584fbSIlya Fursov suffix: rk
267fa9584fbSIlya Fursov diff_args: -j
268fa9584fbSIlya Fursov args: -ts_type rk
269fa9584fbSIlya Fursov output_file: output/ex1.out
270fa9584fbSIlya Fursov
271fa9584fbSIlya Fursov test:
272fa9584fbSIlya Fursov suffix: beuler
273fa9584fbSIlya Fursov diff_args: -j
274fa9584fbSIlya Fursov args: -ts_type beuler
275ecc87898SStefano Zampini output_file: output/ex1_theta.out
276fa9584fbSIlya Fursov
277fa9584fbSIlya Fursov test:
278fa9584fbSIlya Fursov suffix: cn
279fa9584fbSIlya Fursov diff_args: -j
280fa9584fbSIlya Fursov args: -ts_type cn
281ecc87898SStefano Zampini output_file: output/ex1_theta.out
282fa9584fbSIlya Fursov
283fa9584fbSIlya Fursov test:
284fa9584fbSIlya Fursov suffix: theta
285fa9584fbSIlya Fursov args: -ts_type theta
286fa9584fbSIlya Fursov diff_args: -j
287ecc87898SStefano Zampini output_file: output/ex1_theta.out
288fa9584fbSIlya Fursov
289fa9584fbSIlya Fursov test:
290fa9584fbSIlya Fursov suffix: bdf
291fa9584fbSIlya Fursov diff_args: -j
292fa9584fbSIlya Fursov args: -ts_type bdf
293fa9584fbSIlya Fursov output_file: output/ex1_bdf.out
294fa9584fbSIlya Fursov
295fa9584fbSIlya Fursov test:
296fa9584fbSIlya Fursov suffix: alpha
297fa9584fbSIlya Fursov diff_args: -j
298fa9584fbSIlya Fursov args: -ts_type alpha
299ecc87898SStefano Zampini output_file: output/ex1_alpha.out
300fa9584fbSIlya Fursov
301fa9584fbSIlya Fursov test:
302fa9584fbSIlya Fursov suffix: rosw
303fa9584fbSIlya Fursov diff_args: -j
304fa9584fbSIlya Fursov args: -ts_type rosw
305fa9584fbSIlya Fursov output_file: output/ex1.out
306fa9584fbSIlya Fursov
307fa9584fbSIlya Fursov test:
308fa9584fbSIlya Fursov suffix: arkimex
309fa9584fbSIlya Fursov diff_args: -j
310fa9584fbSIlya Fursov args: -ts_type arkimex
311ecc87898SStefano Zampini output_file: output/ex1_arkimex.out
312fa9584fbSIlya Fursov
313fa9584fbSIlya Fursov TEST*/
314