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
369371c9d4SSatish Balay typedef enum {
379371c9d4SSatish Balay VAR_CONSERVATIVE,
389371c9d4SSatish Balay VAR_NONCONSERVATIVE,
399371c9d4SSatish Balay VAR_TRANSIENTVAR
409371c9d4SSatish Balay } VarMode;
41c4762a1bSJed Brown static const char *const VarModes[] = {"CONSERVATIVE", "NONCONSERVATIVE", "TRANSIENTVAR", "VarMode", "VAR_", NULL};
42c4762a1bSJed Brown
IFunction_Conservative(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,PetscCtx ctx)43*2a8381b2SBarry Smith static PetscErrorCode IFunction_Conservative(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx)
44d71ae5a4SJacob Faibussowitsch {
45c4762a1bSJed Brown const PetscScalar *u, *udot;
46c4762a1bSJed Brown PetscScalar *f;
47c4762a1bSJed Brown
487510d9b0SBarry Smith PetscFunctionBeginUser;
499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot));
519566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
52c4762a1bSJed Brown
53c4762a1bSJed Brown f[0] = udot[0] + u[0];
54c4762a1bSJed Brown f[1] = udot[1] - u[0];
55c4762a1bSJed Brown
569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot));
589566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
60c4762a1bSJed Brown }
61c4762a1bSJed Brown
IFunction_Nonconservative(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,PetscCtx ctx)62*2a8381b2SBarry Smith static PetscErrorCode IFunction_Nonconservative(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx)
63d71ae5a4SJacob Faibussowitsch {
64c4762a1bSJed Brown const PetscScalar *u, *udot;
65c4762a1bSJed Brown PetscScalar *f;
66c4762a1bSJed Brown
677510d9b0SBarry Smith PetscFunctionBeginUser;
689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot));
709566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
71c4762a1bSJed Brown
72c4762a1bSJed Brown f[0] = PetscExpScalar(u[0]) * udot[0] + PetscExpScalar(u[0]);
73c4762a1bSJed Brown f[1] = PetscExpScalar(u[1]) * udot[1] - PetscExpScalar(u[0]);
74c4762a1bSJed Brown
759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot));
779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
79c4762a1bSJed Brown }
80c4762a1bSJed Brown
IFunction_TransientVar(TS ts,PetscReal t,Vec U,Vec Cdot,Vec F,PetscCtx ctx)81*2a8381b2SBarry Smith static PetscErrorCode IFunction_TransientVar(TS ts, PetscReal t, Vec U, Vec Cdot, Vec F, PetscCtx ctx)
82d71ae5a4SJacob Faibussowitsch {
83c4762a1bSJed Brown const PetscScalar *u, *cdot;
84c4762a1bSJed Brown PetscScalar *f;
85c4762a1bSJed Brown
867510d9b0SBarry Smith PetscFunctionBeginUser;
879566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Cdot, &cdot));
899566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
90c4762a1bSJed Brown
91c4762a1bSJed Brown f[0] = cdot[0] + PetscExpScalar(u[0]);
92c4762a1bSJed Brown f[1] = cdot[1] - PetscExpScalar(u[0]);
93c4762a1bSJed Brown
949566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Cdot, &cdot));
969566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
98c4762a1bSJed Brown }
99c4762a1bSJed Brown
TransientVar(TS ts,Vec U,Vec C,PetscCtx ctx)100*2a8381b2SBarry Smith static PetscErrorCode TransientVar(TS ts, Vec U, Vec C, PetscCtx ctx)
101d71ae5a4SJacob Faibussowitsch {
1027510d9b0SBarry Smith PetscFunctionBeginUser;
1039566063dSJacob Faibussowitsch PetscCall(VecCopy(U, C));
1049566063dSJacob Faibussowitsch PetscCall(VecExp(C));
1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
106c4762a1bSJed Brown }
107c4762a1bSJed Brown
main(int argc,char * argv[])108d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
109d71ae5a4SJacob Faibussowitsch {
110c4762a1bSJed Brown TS ts;
111c4762a1bSJed Brown DM dm;
112c4762a1bSJed Brown Vec U;
113c4762a1bSJed Brown VarMode var = VAR_CONSERVATIVE;
114c4762a1bSJed Brown PetscScalar sum;
115c4762a1bSJed Brown
116327415f7SBarry Smith PetscFunctionBeginUser;
1179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
118d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "TS conservation example", "");
1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-var", "Variable formulation", NULL, VarModes, (PetscEnum)var, (PetscEnum *)&var, NULL));
120d0609cedSBarry Smith PetscOptionsEnd();
121c4762a1bSJed Brown
1229566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1239566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBDF));
1249566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
1259566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 2, &U));
1269566063dSJacob Faibussowitsch PetscCall(VecSetValue(U, 0, 2., INSERT_VALUES));
1279566063dSJacob Faibussowitsch PetscCall(VecSetValue(U, 1, 1., INSERT_VALUES));
128c4762a1bSJed Brown switch (var) {
129d71ae5a4SJacob Faibussowitsch case VAR_CONSERVATIVE:
130d71ae5a4SJacob Faibussowitsch PetscCall(DMTSSetIFunction(dm, IFunction_Conservative, NULL));
131d71ae5a4SJacob Faibussowitsch break;
132c4762a1bSJed Brown case VAR_NONCONSERVATIVE:
1339566063dSJacob Faibussowitsch PetscCall(VecLog(U));
1349566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunction(dm, IFunction_Nonconservative, NULL));
135c4762a1bSJed Brown break;
136c4762a1bSJed Brown case VAR_TRANSIENTVAR:
1379566063dSJacob Faibussowitsch PetscCall(VecLog(U));
1389566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunction(dm, IFunction_TransientVar, NULL));
1399566063dSJacob Faibussowitsch PetscCall(DMTSSetTransientVariable(dm, TransientVar, NULL));
140c4762a1bSJed Brown }
1419566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1.));
1429566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
143c4762a1bSJed Brown
1449566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
145c4762a1bSJed Brown switch (var) {
146d71ae5a4SJacob Faibussowitsch case VAR_CONSERVATIVE:
147d71ae5a4SJacob Faibussowitsch break;
148c4762a1bSJed Brown case VAR_NONCONSERVATIVE:
149d71ae5a4SJacob Faibussowitsch case VAR_TRANSIENTVAR:
150d71ae5a4SJacob Faibussowitsch PetscCall(VecExp(U));
151d71ae5a4SJacob Faibussowitsch break;
152c4762a1bSJed Brown }
1539566063dSJacob Faibussowitsch PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
1549566063dSJacob Faibussowitsch PetscCall(VecSum(U, &sum));
155300f1712SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Conservation error %g\n", (double)PetscRealPart(sum - 3)));
156c4762a1bSJed Brown
1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U));
1589566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1599566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
160b122ec5aSJacob Faibussowitsch return 0;
161c4762a1bSJed Brown }
162c4762a1bSJed Brown
163c4762a1bSJed Brown /*TEST
164c4762a1bSJed Brown
165c4762a1bSJed Brown test:
166c4762a1bSJed Brown suffix: conservative
167c4762a1bSJed Brown args: -snes_fd -var conservative
168c4762a1bSJed Brown test:
169c4762a1bSJed Brown suffix: nonconservative
170c4762a1bSJed Brown args: -snes_fd -var nonconservative
171c4762a1bSJed Brown test:
172c4762a1bSJed Brown suffix: transientvar
173c4762a1bSJed Brown args: -snes_fd -var transientvar
174c4762a1bSJed Brown
175c4762a1bSJed Brown TEST*/
176