xref: /petsc/src/ts/tests/ex18.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 
main(int argc,char ** argv)19d27334e2SStefano Zampini int main(int argc, char **argv)
20d27334e2SStefano Zampini {
21d27334e2SStefano Zampini   TS        ts;
22d27334e2SStefano Zampini   Vec       x;
230467964bSStefano Zampini   PetscBool dae = PETSC_TRUE, random = PETSC_FALSE;
24d27334e2SStefano Zampini 
25d27334e2SStefano Zampini   PetscFunctionBeginUser;
26c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
27d27334e2SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dae", &dae, NULL));
280467964bSStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-random", &random, NULL));
29d27334e2SStefano Zampini 
30d27334e2SStefano Zampini   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
31d27334e2SStefano Zampini   PetscCall(TSSetIFunction(ts, NULL, IFunction, &dae));
32d27334e2SStefano Zampini   PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobian, &dae));
33d27334e2SStefano Zampini 
34d27334e2SStefano Zampini   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
35d27334e2SStefano Zampini   PetscCall(VecSetSizes(x, 2, PETSC_DECIDE));
36d27334e2SStefano Zampini   PetscCall(VecSetFromOptions(x));
37d27334e2SStefano Zampini   PetscCall(VecSetUp(x));
380467964bSStefano Zampini   if (random) {
390467964bSStefano Zampini     PetscCall(VecSetRandom(x, NULL));
400467964bSStefano Zampini     PetscCall(VecRealPart(x));
410467964bSStefano Zampini   } else PetscCall(VecSet(x, 0.5));
42d27334e2SStefano Zampini   PetscCall(TSSetSolution(ts, x));
43d27334e2SStefano Zampini   PetscCall(VecDestroy(&x));
44d27334e2SStefano Zampini 
45d27334e2SStefano Zampini   PetscCall(TSSetTimeStep(ts, 1.0 / 16.0));
46d27334e2SStefano Zampini   PetscCall(TSSetMaxTime(ts, 1.0));
47d27334e2SStefano Zampini   PetscCall(TSSetFromOptions(ts));
48d27334e2SStefano Zampini   PetscCall(TSSolve(ts, NULL));
49d27334e2SStefano Zampini 
50d27334e2SStefano Zampini   PetscCall(TSDestroy(&ts));
51d27334e2SStefano Zampini   PetscCall(PetscFinalize());
52d27334e2SStefano Zampini   return 0;
53d27334e2SStefano Zampini }
54d27334e2SStefano Zampini 
IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,PetscCtx ctx)55*2a8381b2SBarry Smith PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx)
56d27334e2SStefano Zampini {
57d27334e2SStefano Zampini   const PetscScalar *xdot, *x;
58d27334e2SStefano Zampini   PetscScalar       *f;
59d27334e2SStefano Zampini   PetscBool          dae = *(PetscBool *)(ctx);
60d27334e2SStefano Zampini 
61d27334e2SStefano Zampini   PetscFunctionBeginUser;
62d27334e2SStefano Zampini   PetscCall(VecGetArrayRead(Xdot, &xdot));
63d27334e2SStefano Zampini   PetscCall(VecGetArrayRead(X, &x));
64d27334e2SStefano Zampini   PetscCall(VecGetArrayWrite(F, &f));
65d27334e2SStefano Zampini   if (dae) {
66d27334e2SStefano Zampini     f[0] = x[0] * xdot[0] - x[0] * x[1];
67d27334e2SStefano Zampini     f[1] = x[0] - x[1];
68d27334e2SStefano Zampini   } else {
69d27334e2SStefano Zampini     f[0] = xdot[0] - x[0];
70d27334e2SStefano Zampini     f[1] = xdot[1] - x[1];
71d27334e2SStefano Zampini   }
72d27334e2SStefano Zampini   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
73d27334e2SStefano Zampini   PetscCall(VecRestoreArrayRead(X, &x));
74d27334e2SStefano Zampini   PetscCall(VecRestoreArrayWrite(F, &f));
75d27334e2SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
76d27334e2SStefano Zampini }
77d27334e2SStefano Zampini 
IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat A,Mat B,PetscCtx ctx)78*2a8381b2SBarry Smith PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal shift, Mat A, Mat B, PetscCtx ctx)
79d27334e2SStefano Zampini {
80d27334e2SStefano Zampini   const PetscScalar *xdot, *x;
81d27334e2SStefano Zampini   PetscBool          dae = *(PetscBool *)(ctx);
82d27334e2SStefano Zampini 
83d27334e2SStefano Zampini   PetscFunctionBeginUser;
84d27334e2SStefano Zampini   PetscCall(VecGetArrayRead(Xdot, &xdot));
85d27334e2SStefano Zampini   PetscCall(VecGetArrayRead(X, &x));
86d27334e2SStefano Zampini   if (dae) {
87d27334e2SStefano Zampini     PetscCall(MatSetValue(B, 0, 0, shift * x[0] + xdot[0] - x[1], INSERT_VALUES));
88d27334e2SStefano Zampini     PetscCall(MatSetValue(B, 0, 1, -x[0], INSERT_VALUES));
89d27334e2SStefano Zampini     PetscCall(MatSetValue(B, 1, 0, 1.0, INSERT_VALUES));
90d27334e2SStefano Zampini     PetscCall(MatSetValue(B, 1, 1, -1.0, INSERT_VALUES));
91d27334e2SStefano Zampini   } else {
92d27334e2SStefano Zampini     PetscCall(MatZeroEntries(B));
93d27334e2SStefano Zampini     PetscCall(MatShift(B, shift));
94d27334e2SStefano Zampini   }
95d27334e2SStefano Zampini   PetscCall(VecRestoreArrayRead(X, &x));
96d27334e2SStefano Zampini   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
97d27334e2SStefano Zampini   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
98d27334e2SStefano Zampini   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
99d27334e2SStefano Zampini   if (A != B) {
100d27334e2SStefano Zampini     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
101d27334e2SStefano Zampini     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
102d27334e2SStefano Zampini   }
103d27334e2SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
104d27334e2SStefano Zampini }
105d27334e2SStefano Zampini 
106d27334e2SStefano Zampini /*TEST
107d27334e2SStefano Zampini 
108d27334e2SStefano Zampini   testset:
109188af4bfSBarry Smith     args: -ts_view_solution -ts_max_steps 10 -ts_time_step 0.1 -ts_view_solution -ts_adapt_type {{none basic}} -ts_exact_final_time matchstep -snes_error_if_not_converged
110d27334e2SStefano Zampini 
111d27334e2SStefano Zampini     test:
112d27334e2SStefano Zampini       output_file: output/ex18_1.out
113d27334e2SStefano Zampini       suffix: bdf
114d27334e2SStefano Zampini       args: -ts_type bdf
115d27334e2SStefano Zampini 
116d27334e2SStefano Zampini     test:
117d27334e2SStefano Zampini       output_file: output/ex18_1.out
118d27334e2SStefano Zampini       suffix: dirk
1193405e92cSStefano 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}}
120d27334e2SStefano Zampini 
1210467964bSStefano Zampini     test:
1220467964bSStefano Zampini       output_file: output/ex18_1.out
1230467964bSStefano Zampini       suffix: dirk_explicit_first_random_dae
124188af4bfSBarry Smith       args: -dae -ts_type dirk -ts_dirk_type es122sal -random 1 -ts_max_step_rejections -1
1250467964bSStefano Zampini 
126d27334e2SStefano Zampini TEST*/
127