xref: /petsc/src/ts/tests/ex18.c (revision 3405e92ce0a429885711dcb2a39f6cd8565d5879)
1d27334e2SStefano Zampini static char help[] = "Solves a DAE with a non-trivial mass matrix. \n\n";
2d27334e2SStefano Zampini /*
3d27334e2SStefano Zampini    Solves:
4d27334e2SStefano Zampini         U * dU/dt = U*V
5d27334e2SStefano Zampini         U - V = 0
6d27334e2SStefano Zampini 
7d27334e2SStefano Zampini    that can be rewritten in implicit form F = 0, with F as
8d27334e2SStefano Zampini                  x[0] * xdot[0] - x[0] * x[1]
9d27334e2SStefano Zampini    F(t,x,xdot) =
10d27334e2SStefano Zampini                  x[0] - x[1]
11d27334e2SStefano Zampini    It is equivalent to solve dU/dt = U, U = U0 with solution U = U0 * exp(tfinal)
12d27334e2SStefano Zampini */
13d27334e2SStefano Zampini 
14d27334e2SStefano Zampini #include <petscts.h>
15d27334e2SStefano Zampini 
16d27334e2SStefano Zampini PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *);
17d27334e2SStefano Zampini PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
18d27334e2SStefano Zampini 
19d27334e2SStefano Zampini int main(int argc, char **argv)
20d27334e2SStefano Zampini {
21d27334e2SStefano Zampini   TS        ts;
22d27334e2SStefano Zampini   Vec       x;
23d27334e2SStefano Zampini   PetscBool dae = PETSC_TRUE;
24d27334e2SStefano Zampini 
25d27334e2SStefano Zampini   PetscFunctionBeginUser;
26d27334e2SStefano Zampini   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
27d27334e2SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dae", &dae, NULL));
28d27334e2SStefano Zampini 
29d27334e2SStefano Zampini   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
30d27334e2SStefano Zampini   PetscCall(TSSetIFunction(ts, NULL, IFunction, &dae));
31d27334e2SStefano Zampini   PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobian, &dae));
32d27334e2SStefano Zampini 
33d27334e2SStefano Zampini   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
34d27334e2SStefano Zampini   PetscCall(VecSetSizes(x, 2, PETSC_DECIDE));
35d27334e2SStefano Zampini   PetscCall(VecSetFromOptions(x));
36d27334e2SStefano Zampini   PetscCall(VecSetUp(x));
37d27334e2SStefano Zampini   PetscCall(VecSet(x, 0.5));
38d27334e2SStefano Zampini   PetscCall(TSSetSolution(ts, x));
39d27334e2SStefano Zampini   PetscCall(VecDestroy(&x));
40d27334e2SStefano Zampini 
41d27334e2SStefano Zampini   PetscCall(TSSetTimeStep(ts, 1.0 / 16.0));
42d27334e2SStefano Zampini   PetscCall(TSSetMaxTime(ts, 1.0));
43d27334e2SStefano Zampini   PetscCall(TSSetFromOptions(ts));
44d27334e2SStefano Zampini   PetscCall(TSSolve(ts, NULL));
45d27334e2SStefano Zampini 
46d27334e2SStefano Zampini   PetscCall(TSDestroy(&ts));
47d27334e2SStefano Zampini   PetscCall(PetscFinalize());
48d27334e2SStefano Zampini   return 0;
49d27334e2SStefano Zampini }
50d27334e2SStefano Zampini 
51d27334e2SStefano Zampini PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
52d27334e2SStefano Zampini {
53d27334e2SStefano Zampini   const PetscScalar *xdot, *x;
54d27334e2SStefano Zampini   PetscScalar       *f;
55d27334e2SStefano Zampini   PetscBool          dae = *(PetscBool *)(ctx);
56d27334e2SStefano Zampini 
57d27334e2SStefano Zampini   PetscFunctionBeginUser;
58d27334e2SStefano Zampini   PetscCall(VecGetArrayRead(Xdot, &xdot));
59d27334e2SStefano Zampini   PetscCall(VecGetArrayRead(X, &x));
60d27334e2SStefano Zampini   PetscCall(VecGetArrayWrite(F, &f));
61d27334e2SStefano Zampini   if (dae) {
62d27334e2SStefano Zampini     f[0] = x[0] * xdot[0] - x[0] * x[1];
63d27334e2SStefano Zampini     f[1] = x[0] - x[1];
64d27334e2SStefano Zampini   } else {
65d27334e2SStefano Zampini     f[0] = xdot[0] - x[0];
66d27334e2SStefano Zampini     f[1] = xdot[1] - x[1];
67d27334e2SStefano Zampini   }
68d27334e2SStefano Zampini   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
69d27334e2SStefano Zampini   PetscCall(VecRestoreArrayRead(X, &x));
70d27334e2SStefano Zampini   PetscCall(VecRestoreArrayWrite(F, &f));
71d27334e2SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
72d27334e2SStefano Zampini }
73d27334e2SStefano Zampini 
74d27334e2SStefano Zampini PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal shift, Mat A, Mat B, void *ctx)
75d27334e2SStefano Zampini {
76d27334e2SStefano Zampini   const PetscScalar *xdot, *x;
77d27334e2SStefano Zampini   PetscBool          dae = *(PetscBool *)(ctx);
78d27334e2SStefano Zampini 
79d27334e2SStefano Zampini   PetscFunctionBeginUser;
80d27334e2SStefano Zampini   PetscCall(VecGetArrayRead(Xdot, &xdot));
81d27334e2SStefano Zampini   PetscCall(VecGetArrayRead(X, &x));
82d27334e2SStefano Zampini   if (dae) {
83d27334e2SStefano Zampini     PetscCall(MatSetValue(B, 0, 0, shift * x[0] + xdot[0] - x[1], INSERT_VALUES));
84d27334e2SStefano Zampini     PetscCall(MatSetValue(B, 0, 1, -x[0], INSERT_VALUES));
85d27334e2SStefano Zampini     PetscCall(MatSetValue(B, 1, 0, 1.0, INSERT_VALUES));
86d27334e2SStefano Zampini     PetscCall(MatSetValue(B, 1, 1, -1.0, INSERT_VALUES));
87d27334e2SStefano Zampini   } else {
88d27334e2SStefano Zampini     PetscCall(MatZeroEntries(B));
89d27334e2SStefano Zampini     PetscCall(MatShift(B, shift));
90d27334e2SStefano Zampini   }
91d27334e2SStefano Zampini   PetscCall(VecRestoreArrayRead(X, &x));
92d27334e2SStefano Zampini   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
93d27334e2SStefano Zampini   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
94d27334e2SStefano Zampini   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
95d27334e2SStefano Zampini   if (A != B) {
96d27334e2SStefano Zampini     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
97d27334e2SStefano Zampini     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
98d27334e2SStefano Zampini   }
99d27334e2SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
100d27334e2SStefano Zampini }
101d27334e2SStefano Zampini 
102d27334e2SStefano Zampini /*TEST
103d27334e2SStefano Zampini 
104d27334e2SStefano Zampini   testset:
105d27334e2SStefano Zampini     args: -ts_view_solution -ts_max_steps 10 -ts_dt 0.1 -ts_view_solution -ts_adapt_type {{none basic}} -ts_exact_final_time matchstep
106d27334e2SStefano Zampini 
107d27334e2SStefano Zampini     test:
108d27334e2SStefano Zampini       output_file: output/ex18_1.out
109d27334e2SStefano Zampini       suffix: bdf
110d27334e2SStefano Zampini       args: -ts_type bdf
111d27334e2SStefano Zampini 
112d27334e2SStefano Zampini     test:
113d27334e2SStefano Zampini       output_file: output/ex18_1.out
114d27334e2SStefano Zampini       suffix: dirk
115*3405e92cSStefano Zampini       args: -dae {{0 1}} -ts_type dirk -ts_dirk_type {{s212 es122sal es213sal es324sal es325sal 657a es648sa 658a s659a 7510sal es7510sa 759a s7511sal 8614a 8616sal es8516sal}}
116d27334e2SStefano Zampini 
117d27334e2SStefano Zampini TEST*/
118