xref: /petsc/src/ts/tests/ex15.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Test conservation properties for 2-variable system\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*F
4c4762a1bSJed Brown We consider a linear reaction system with two concentrations
5c4762a1bSJed Brown 
6c4762a1bSJed Brown \begin{align}
7c4762a1bSJed Brown   \frac{\partial c_0}{\partial t} &= -c_0 \\
8c4762a1bSJed Brown   \frac{\partial c_1}{\partial t} &= c_0,
9c4762a1bSJed Brown \end{align}
10c4762a1bSJed Brown 
11c4762a1bSJed Brown wherethe sum $c_0 + c_1$ is conserved, as can be seen by adding the two equations.
12c4762a1bSJed Brown 
13c4762a1bSJed Brown We now consider a different set of variables, defined implicitly by $c(u) = e^u$.  This type of transformation is
14c4762a1bSJed Brown sometimes used to ensure positivity, and related transformations are sometimes used to develop a well-conditioned
15c4762a1bSJed Brown formulation in limits such as zero Mach number.  In this instance, the relation is explicitly invertible, but that is
16c4762a1bSJed Brown not always the case.  We can rewrite the differential equation in terms of non-conservative variables u,
17c4762a1bSJed Brown 
18c4762a1bSJed Brown \begin{align}
19c4762a1bSJed Brown   \frac{\partial c_0}{\partial u_0} \frac{\partial u_0}{\partial t} &= -c_0(u_0) \\
20c4762a1bSJed Brown   \frac{\partial c_1}{\partial u_1} \frac{\partial u_1}{\partial t} &= c_0(u_0).
21c4762a1bSJed Brown \end{align}
22c4762a1bSJed Brown 
23c4762a1bSJed Brown We'll consider this three ways, each using an IFunction
24c4762a1bSJed Brown 
25c4762a1bSJed Brown 1. CONSERVATIVE: standard integration in conservative variables: F(C, Cdot) = 0
26c4762a1bSJed Brown 2. NONCONSERVATIVE: chain rule formulation entirely in primitive variables: F(U, Udot) = 0
27c4762a1bSJed Brown 3. TRANSIENTVAR: Provide function C(U) and solve F(U, Cdot) = 0, where the time integrators handles the transformation
28c4762a1bSJed Brown 
29c4762a1bSJed Brown We will see that 1 and 3 are conservative (up to machine precision/solver tolerance, independent of temporal
30c4762a1bSJed Brown discretization error) while 2 is not conservative (i.e., scales with temporal discretization error).
31c4762a1bSJed Brown 
32c4762a1bSJed Brown F*/
33c4762a1bSJed Brown 
34c4762a1bSJed Brown #include <petscts.h>
35c4762a1bSJed Brown 
36*9371c9d4SSatish Balay typedef enum {
37*9371c9d4SSatish Balay   VAR_CONSERVATIVE,
38*9371c9d4SSatish Balay   VAR_NONCONSERVATIVE,
39*9371c9d4SSatish Balay   VAR_TRANSIENTVAR
40*9371c9d4SSatish Balay } VarMode;
41c4762a1bSJed Brown static const char *const VarModes[] = {"CONSERVATIVE", "NONCONSERVATIVE", "TRANSIENTVAR", "VarMode", "VAR_", NULL};
42c4762a1bSJed Brown 
43*9371c9d4SSatish Balay static PetscErrorCode IFunction_Conservative(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
44c4762a1bSJed Brown   const PetscScalar *u, *udot;
45c4762a1bSJed Brown   PetscScalar       *f;
46c4762a1bSJed Brown 
477510d9b0SBarry Smith   PetscFunctionBeginUser;
489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   f[0] = udot[0] + u[0];
53c4762a1bSJed Brown   f[1] = udot[1] - u[0];
54c4762a1bSJed Brown 
559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
58c4762a1bSJed Brown   PetscFunctionReturn(0);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61*9371c9d4SSatish Balay static PetscErrorCode IFunction_Nonconservative(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
62c4762a1bSJed Brown   const PetscScalar *u, *udot;
63c4762a1bSJed Brown   PetscScalar       *f;
64c4762a1bSJed Brown 
657510d9b0SBarry Smith   PetscFunctionBeginUser;
669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   f[0] = PetscExpScalar(u[0]) * udot[0] + PetscExpScalar(u[0]);
71c4762a1bSJed Brown   f[1] = PetscExpScalar(u[1]) * udot[1] - PetscExpScalar(u[0]);
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
76c4762a1bSJed Brown   PetscFunctionReturn(0);
77c4762a1bSJed Brown }
78c4762a1bSJed Brown 
79*9371c9d4SSatish Balay static PetscErrorCode IFunction_TransientVar(TS ts, PetscReal t, Vec U, Vec Cdot, Vec F, void *ctx) {
80c4762a1bSJed Brown   const PetscScalar *u, *cdot;
81c4762a1bSJed Brown   PetscScalar       *f;
82c4762a1bSJed Brown 
837510d9b0SBarry Smith   PetscFunctionBeginUser;
849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Cdot, &cdot));
869566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   f[0] = cdot[0] + PetscExpScalar(u[0]);
89c4762a1bSJed Brown   f[1] = cdot[1] - PetscExpScalar(u[0]);
90c4762a1bSJed Brown 
919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Cdot, &cdot));
939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
94c4762a1bSJed Brown   PetscFunctionReturn(0);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97*9371c9d4SSatish Balay static PetscErrorCode TransientVar(TS ts, Vec U, Vec C, void *ctx) {
987510d9b0SBarry Smith   PetscFunctionBeginUser;
999566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, C));
1009566063dSJacob Faibussowitsch   PetscCall(VecExp(C));
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104*9371c9d4SSatish Balay int main(int argc, char *argv[]) {
105c4762a1bSJed Brown   TS          ts;
106c4762a1bSJed Brown   DM          dm;
107c4762a1bSJed Brown   Vec         U;
108c4762a1bSJed Brown   VarMode     var = VAR_CONSERVATIVE;
109c4762a1bSJed Brown   PetscScalar sum;
110c4762a1bSJed Brown 
111327415f7SBarry Smith   PetscFunctionBeginUser;
1129566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
113d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "TS conservation example", "");
1149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-var", "Variable formulation", NULL, VarModes, (PetscEnum)var, (PetscEnum *)&var, NULL));
115d0609cedSBarry Smith   PetscOptionsEnd();
116c4762a1bSJed Brown 
1179566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1189566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBDF));
1199566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
1209566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 2, &U));
1219566063dSJacob Faibussowitsch   PetscCall(VecSetValue(U, 0, 2., INSERT_VALUES));
1229566063dSJacob Faibussowitsch   PetscCall(VecSetValue(U, 1, 1., INSERT_VALUES));
123c4762a1bSJed Brown   switch (var) {
124*9371c9d4SSatish Balay   case VAR_CONSERVATIVE: PetscCall(DMTSSetIFunction(dm, IFunction_Conservative, NULL)); break;
125c4762a1bSJed Brown   case VAR_NONCONSERVATIVE:
1269566063dSJacob Faibussowitsch     PetscCall(VecLog(U));
1279566063dSJacob Faibussowitsch     PetscCall(DMTSSetIFunction(dm, IFunction_Nonconservative, NULL));
128c4762a1bSJed Brown     break;
129c4762a1bSJed Brown   case VAR_TRANSIENTVAR:
1309566063dSJacob Faibussowitsch     PetscCall(VecLog(U));
1319566063dSJacob Faibussowitsch     PetscCall(DMTSSetIFunction(dm, IFunction_TransientVar, NULL));
1329566063dSJacob Faibussowitsch     PetscCall(DMTSSetTransientVariable(dm, TransientVar, NULL));
133c4762a1bSJed Brown   }
1349566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 1.));
1359566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
138c4762a1bSJed Brown   switch (var) {
139*9371c9d4SSatish Balay   case VAR_CONSERVATIVE: break;
140c4762a1bSJed Brown   case VAR_NONCONSERVATIVE:
141*9371c9d4SSatish Balay   case VAR_TRANSIENTVAR: PetscCall(VecExp(U)); break;
142c4762a1bSJed Brown   }
1439566063dSJacob Faibussowitsch   PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
1449566063dSJacob Faibussowitsch   PetscCall(VecSum(U, &sum));
14563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Conservation error %g\n", (double)PetscRealPart(sum - 3.)));
146c4762a1bSJed Brown 
1479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
1489566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1499566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
150b122ec5aSJacob Faibussowitsch   return 0;
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown /*TEST
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   test:
156c4762a1bSJed Brown     suffix: conservative
157c4762a1bSJed Brown     args: -snes_fd -var conservative
158c4762a1bSJed Brown   test:
159c4762a1bSJed Brown     suffix: nonconservative
160c4762a1bSJed Brown     args: -snes_fd -var nonconservative
161c4762a1bSJed Brown   test:
162c4762a1bSJed Brown     suffix: transientvar
163c4762a1bSJed Brown     args: -snes_fd -var transientvar
164c4762a1bSJed Brown 
165c4762a1bSJed Brown TEST*/
166