xref: /petsc/src/ts/tests/ex17.c (revision 4a2752e6099aa818be98754c4b44c7a908e87da5)
1*4a2752e6SStefano Zampini static char help[] = "Solves the ODE du/dt = poly(t), u(0) = 0. Tests TSResize for varying size.\n\n";
2*4a2752e6SStefano Zampini 
3*4a2752e6SStefano Zampini #include <petscts.h>
4*4a2752e6SStefano Zampini 
5*4a2752e6SStefano Zampini PetscScalar poly(PetscInt p, PetscReal t)
6*4a2752e6SStefano Zampini {
7*4a2752e6SStefano Zampini   return p ? t * poly(p - 1, t) : 1.0;
8*4a2752e6SStefano Zampini }
9*4a2752e6SStefano Zampini 
10*4a2752e6SStefano Zampini PetscScalar dpoly(PetscInt p, PetscReal t)
11*4a2752e6SStefano Zampini {
12*4a2752e6SStefano Zampini   return p > 0 ? (PetscReal)p * poly(p - 1, t) : 0.0;
13*4a2752e6SStefano Zampini }
14*4a2752e6SStefano Zampini 
15*4a2752e6SStefano Zampini PetscErrorCode CreateVec(PetscInt lsize, Vec *out)
16*4a2752e6SStefano Zampini {
17*4a2752e6SStefano Zampini   Vec x;
18*4a2752e6SStefano Zampini 
19*4a2752e6SStefano Zampini   PetscFunctionBeginUser;
20*4a2752e6SStefano Zampini   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
21*4a2752e6SStefano Zampini   PetscCall(VecSetSizes(x, lsize, PETSC_DECIDE));
22*4a2752e6SStefano Zampini   PetscCall(VecSetFromOptions(x));
23*4a2752e6SStefano Zampini   PetscCall(VecSetUp(x));
24*4a2752e6SStefano Zampini   *out = x;
25*4a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
26*4a2752e6SStefano Zampini }
27*4a2752e6SStefano Zampini 
28*4a2752e6SStefano Zampini PetscErrorCode CreateMat(PetscInt lsize, Mat *out)
29*4a2752e6SStefano Zampini {
30*4a2752e6SStefano Zampini   Mat A;
31*4a2752e6SStefano Zampini 
32*4a2752e6SStefano Zampini   PetscFunctionBeginUser;
33*4a2752e6SStefano Zampini   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
34*4a2752e6SStefano Zampini   PetscCall(MatSetSizes(A, lsize, lsize, PETSC_DECIDE, PETSC_DECIDE));
35*4a2752e6SStefano Zampini   PetscCall(MatSetFromOptions(A));
36*4a2752e6SStefano Zampini   PetscCall(MatSetUp(A));
37*4a2752e6SStefano Zampini   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
38*4a2752e6SStefano Zampini   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
39*4a2752e6SStefano Zampini   /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */
40*4a2752e6SStefano Zampini   PetscCall(MatShift(A, (PetscReal)1));
41*4a2752e6SStefano Zampini   PetscCall(MatShift(A, (PetscReal)-1));
42*4a2752e6SStefano Zampini   *out = A;
43*4a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
44*4a2752e6SStefano Zampini }
45*4a2752e6SStefano Zampini 
46*4a2752e6SStefano Zampini PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec x, Vec f, void *ctx)
47*4a2752e6SStefano Zampini {
48*4a2752e6SStefano Zampini   PetscInt *order = (PetscInt *)ctx;
49*4a2752e6SStefano Zampini 
50*4a2752e6SStefano Zampini   PetscFunctionBeginUser;
51*4a2752e6SStefano Zampini   PetscCall(VecSet(f, dpoly(*order, t)));
52*4a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
53*4a2752e6SStefano Zampini }
54*4a2752e6SStefano Zampini 
55*4a2752e6SStefano Zampini PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat B, void *ctx)
56*4a2752e6SStefano Zampini {
57*4a2752e6SStefano Zampini   PetscFunctionBeginUser;
58*4a2752e6SStefano Zampini   PetscCall(MatZeroEntries(B));
59*4a2752e6SStefano Zampini   if (B != A) PetscCall(MatZeroEntries(A));
60*4a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
61*4a2752e6SStefano Zampini }
62*4a2752e6SStefano Zampini 
63*4a2752e6SStefano Zampini PetscErrorCode Transfer(TS ts, PetscInt nv, Vec vecsin[], Vec vecsout[], void *ctx)
64*4a2752e6SStefano Zampini {
65*4a2752e6SStefano Zampini   PetscInt n, nnew;
66*4a2752e6SStefano Zampini 
67*4a2752e6SStefano Zampini   PetscFunctionBeginUser;
68*4a2752e6SStefano Zampini   PetscAssert(nv > 0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Zero vectors");
69*4a2752e6SStefano Zampini   PetscCall(VecGetLocalSize(vecsin[0], &n));
70*4a2752e6SStefano Zampini   nnew = n == 2 ? 1 : 2;
71*4a2752e6SStefano Zampini   for (PetscInt i = 0; i < nv; i++) {
72*4a2752e6SStefano Zampini     const PetscScalar *vals;
73*4a2752e6SStefano Zampini 
74*4a2752e6SStefano Zampini     PetscCall(CreateVec(nnew, &vecsout[i]));
75*4a2752e6SStefano Zampini     PetscCall(VecGetArrayRead(vecsin[i], &vals));
76*4a2752e6SStefano Zampini     PetscCall(VecSet(vecsout[i], vals[0]));
77*4a2752e6SStefano Zampini     PetscCall(VecRestoreArrayRead(vecsin[i], &vals));
78*4a2752e6SStefano Zampini   }
79*4a2752e6SStefano Zampini   Mat A;
80*4a2752e6SStefano Zampini   PetscCall(CreateMat(nnew, &A));
81*4a2752e6SStefano Zampini   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL));
82*4a2752e6SStefano Zampini   PetscCall(MatDestroy(&A));
83*4a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
84*4a2752e6SStefano Zampini }
85*4a2752e6SStefano Zampini 
86*4a2752e6SStefano Zampini PetscErrorCode TransferSetUp(TS ts, PetscInt step, PetscReal time, Vec sol, PetscBool *resize, void *ctx)
87*4a2752e6SStefano Zampini {
88*4a2752e6SStefano Zampini   PetscFunctionBeginUser;
89*4a2752e6SStefano Zampini   *resize = PETSC_TRUE;
90*4a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
91*4a2752e6SStefano Zampini }
92*4a2752e6SStefano Zampini 
93*4a2752e6SStefano Zampini PetscErrorCode Monitor(TS ts, PetscInt n, PetscReal t, Vec x, void *ctx)
94*4a2752e6SStefano Zampini {
95*4a2752e6SStefano Zampini   const PetscScalar *a;
96*4a2752e6SStefano Zampini   PetscScalar       *store = (PetscScalar *)ctx;
97*4a2752e6SStefano Zampini 
98*4a2752e6SStefano Zampini   PetscFunctionBeginUser;
99*4a2752e6SStefano Zampini   PetscCall(VecGetArrayRead(x, &a));
100*4a2752e6SStefano Zampini   if (n < 10) store[n] = a[0];
101*4a2752e6SStefano Zampini   PetscCall(VecRestoreArrayRead(x, &a));
102*4a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
103*4a2752e6SStefano Zampini }
104*4a2752e6SStefano Zampini 
105*4a2752e6SStefano Zampini int main(int argc, char **argv)
106*4a2752e6SStefano Zampini {
107*4a2752e6SStefano Zampini   TS          ts;
108*4a2752e6SStefano Zampini   Vec         x;
109*4a2752e6SStefano Zampini   Mat         A;
110*4a2752e6SStefano Zampini   PetscInt    order = 2;
111*4a2752e6SStefano Zampini   PetscScalar results[2][10];
112*4a2752e6SStefano Zampini   /* I would like to use 0 here, but linux-gcc-complex-opt-32bit
113*4a2752e6SStefano Zampini      errors with arkimex with 1.e-18 errors */
114*4a2752e6SStefano Zampini   PetscReal tol = PETSC_MACHINE_EPSILON;
115*4a2752e6SStefano Zampini 
116*4a2752e6SStefano Zampini   PetscFunctionBeginUser;
117*4a2752e6SStefano Zampini   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
118*4a2752e6SStefano Zampini   PetscCall(PetscOptionsGetInt(NULL, NULL, "-order", &order, NULL));
119*4a2752e6SStefano Zampini   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL));
120*4a2752e6SStefano Zampini 
121*4a2752e6SStefano Zampini   for (PetscInt i = 0; i < 2; i++) {
122*4a2752e6SStefano Zampini     PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
123*4a2752e6SStefano Zampini     PetscCall(TSSetProblemType(ts, TS_LINEAR));
124*4a2752e6SStefano Zampini 
125*4a2752e6SStefano Zampini     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &order));
126*4a2752e6SStefano Zampini 
127*4a2752e6SStefano Zampini     PetscCall(CreateMat(1, &A));
128*4a2752e6SStefano Zampini     PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL));
129*4a2752e6SStefano Zampini     PetscCall(MatDestroy(&A));
130*4a2752e6SStefano Zampini 
131*4a2752e6SStefano Zampini     PetscCall(CreateVec(1, &x));
132*4a2752e6SStefano Zampini     PetscCall(TSSetSolution(ts, x));
133*4a2752e6SStefano Zampini     PetscCall(VecDestroy(&x));
134*4a2752e6SStefano Zampini 
135*4a2752e6SStefano Zampini     for (PetscInt j = 0; j < 10; j++) results[i][j] = 0;
136*4a2752e6SStefano Zampini     PetscCall(TSMonitorSet(ts, Monitor, results[i], NULL));
137*4a2752e6SStefano Zampini     PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
138*4a2752e6SStefano Zampini     if (i) PetscCall(TSSetResize(ts, TransferSetUp, Transfer, NULL));
139*4a2752e6SStefano Zampini     PetscCall(TSSetTime(ts, 0));
140*4a2752e6SStefano Zampini     PetscCall(TSSetTimeStep(ts, 1. / 4.));
141*4a2752e6SStefano Zampini     PetscCall(TSSetMaxSteps(ts, 10));
142*4a2752e6SStefano Zampini     PetscCall(TSSetFromOptions(ts));
143*4a2752e6SStefano Zampini 
144*4a2752e6SStefano Zampini     PetscCall(TSSolve(ts, NULL));
145*4a2752e6SStefano Zampini 
146*4a2752e6SStefano Zampini     PetscCall(TSDestroy(&ts));
147*4a2752e6SStefano Zampini   }
148*4a2752e6SStefano Zampini 
149*4a2752e6SStefano Zampini   /* Dump errors if any */
150*4a2752e6SStefano Zampini   PetscBool flg = PETSC_FALSE;
151*4a2752e6SStefano Zampini   for (PetscInt i = 0; i < 10; i++) {
152*4a2752e6SStefano Zampini     PetscReal err = PetscAbsScalar(results[0][i] - results[1][i]);
153*4a2752e6SStefano Zampini     if (err > tol) {
154*4a2752e6SStefano Zampini       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error step %" PetscInt_FMT ": %g\n", i, (double)err));
155*4a2752e6SStefano Zampini       flg = PETSC_TRUE;
156*4a2752e6SStefano Zampini     }
157*4a2752e6SStefano Zampini   }
158*4a2752e6SStefano Zampini   if (flg) {
159*4a2752e6SStefano Zampini     PetscCall(PetscScalarView(10, results[0], PETSC_VIEWER_STDOUT_WORLD));
160*4a2752e6SStefano Zampini     PetscCall(PetscScalarView(10, results[1], PETSC_VIEWER_STDOUT_WORLD));
161*4a2752e6SStefano Zampini   }
162*4a2752e6SStefano Zampini   PetscCall(PetscFinalize());
163*4a2752e6SStefano Zampini   return 0;
164*4a2752e6SStefano Zampini }
165*4a2752e6SStefano Zampini 
166*4a2752e6SStefano Zampini /*TEST
167*4a2752e6SStefano Zampini 
168*4a2752e6SStefano Zampini     test:
169*4a2752e6SStefano Zampini       suffix: bdf
170*4a2752e6SStefano Zampini       args: -ts_adapt_wnormtype infinity -ts_type bdf -ts_bdf_order {{2 3 4 5 6}} -order 6 -ts_adapt_type {{none basic dsp}} -ksp_type preonly -pc_type lu
171*4a2752e6SStefano Zampini       output_file: output/ex17.out
172*4a2752e6SStefano Zampini 
173*4a2752e6SStefano Zampini     test:
174*4a2752e6SStefano Zampini       suffix: expl
175*4a2752e6SStefano Zampini       args: -ts_adapt_wnormtype infinity -ts_type {{euler rk ssp}} -order 6 -ts_adapt_type {{none basic dsp}}
176*4a2752e6SStefano Zampini       output_file: output/ex17.out
177*4a2752e6SStefano Zampini 
178*4a2752e6SStefano Zampini     test:
179*4a2752e6SStefano Zampini       suffix: impl
180*4a2752e6SStefano Zampini       args: -ts_adapt_wnormtype infinity -ts_type {{rosw beuler cn alpha theta arkimex}} -order 6 -ts_adapt_type {{none basic dsp}} -ksp_type preonly -pc_type lu
181*4a2752e6SStefano Zampini       output_file: output/ex17.out
182*4a2752e6SStefano Zampini 
183*4a2752e6SStefano Zampini TEST*/
184