xref: /petsc/src/ts/tests/ex15.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
36c4762a1bSJed Brown typedef enum {VAR_CONSERVATIVE, VAR_NONCONSERVATIVE, VAR_TRANSIENTVAR} VarMode;
37c4762a1bSJed Brown static const char *const VarModes[] = {"CONSERVATIVE", "NONCONSERVATIVE", "TRANSIENTVAR", "VarMode", "VAR_", NULL};
38c4762a1bSJed Brown 
39c4762a1bSJed Brown static PetscErrorCode IFunction_Conservative(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
40c4762a1bSJed Brown {
41c4762a1bSJed Brown   const PetscScalar *u,*udot;
42c4762a1bSJed Brown   PetscScalar       *f;
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   PetscFunctionBegin;
459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
469566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   f[0] = udot[0] + u[0];
50c4762a1bSJed Brown   f[1] = udot[1] - u[0];
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
55c4762a1bSJed Brown   PetscFunctionReturn(0);
56c4762a1bSJed Brown }
57c4762a1bSJed Brown 
58c4762a1bSJed Brown static PetscErrorCode IFunction_Nonconservative(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
59c4762a1bSJed Brown {
60c4762a1bSJed Brown   const PetscScalar *u,*udot;
61c4762a1bSJed Brown   PetscScalar       *f;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   PetscFunctionBegin;
649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
669566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   f[0] = PetscExpScalar(u[0])*udot[0] + PetscExpScalar(u[0]);
69c4762a1bSJed Brown   f[1] = PetscExpScalar(u[1])*udot[1] - PetscExpScalar(u[0]);
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
74c4762a1bSJed Brown   PetscFunctionReturn(0);
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown static PetscErrorCode IFunction_TransientVar(TS ts,PetscReal t,Vec U,Vec Cdot,Vec F,void *ctx)
78c4762a1bSJed Brown {
79c4762a1bSJed Brown   const PetscScalar *u,*cdot;
80c4762a1bSJed Brown   PetscScalar       *f;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   PetscFunctionBegin;
839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Cdot,&cdot));
859566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   f[0] = cdot[0] + PetscExpScalar(u[0]);
88c4762a1bSJed Brown   f[1] = cdot[1] - PetscExpScalar(u[0]);
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Cdot,&cdot));
929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown static PetscErrorCode TransientVar(TS ts,Vec U,Vec C,void *ctx)
97c4762a1bSJed Brown {
98c4762a1bSJed Brown   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(VecCopy(U,C));
1009566063dSJacob Faibussowitsch   PetscCall(VecExp(C));
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown int main(int argc, char *argv[])
105c4762a1bSJed Brown {
106c4762a1bSJed Brown   TS             ts;
107c4762a1bSJed Brown   DM             dm;
108c4762a1bSJed Brown   Vec            U;
109c4762a1bSJed Brown   VarMode        var = VAR_CONSERVATIVE;
110c4762a1bSJed Brown   PetscScalar    sum;
111c4762a1bSJed Brown 
112*327415f7SBarry Smith   PetscFunctionBeginUser;
1139566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
114d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"TS conservation example","");
1159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-var","Variable formulation",NULL,VarModes,(PetscEnum)var,(PetscEnum*)&var,NULL));
116d0609cedSBarry Smith   PetscOptionsEnd();
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
1199566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSBDF));
1209566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&dm));
1219566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,2,&U));
1229566063dSJacob Faibussowitsch   PetscCall(VecSetValue(U,0,2.,INSERT_VALUES));
1239566063dSJacob Faibussowitsch   PetscCall(VecSetValue(U,1,1.,INSERT_VALUES));
124c4762a1bSJed Brown   switch (var) {
125c4762a1bSJed Brown   case VAR_CONSERVATIVE:
1269566063dSJacob Faibussowitsch     PetscCall(DMTSSetIFunction(dm,IFunction_Conservative,NULL));
127c4762a1bSJed Brown     break;
128c4762a1bSJed Brown   case VAR_NONCONSERVATIVE:
1299566063dSJacob Faibussowitsch     PetscCall(VecLog(U));
1309566063dSJacob Faibussowitsch     PetscCall(DMTSSetIFunction(dm,IFunction_Nonconservative,NULL));
131c4762a1bSJed Brown     break;
132c4762a1bSJed Brown   case VAR_TRANSIENTVAR:
1339566063dSJacob Faibussowitsch     PetscCall(VecLog(U));
1349566063dSJacob Faibussowitsch     PetscCall(DMTSSetIFunction(dm,IFunction_TransientVar,NULL));
1359566063dSJacob Faibussowitsch     PetscCall(DMTSSetTransientVariable(dm,TransientVar,NULL));
136c4762a1bSJed Brown   }
1379566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,1.));
1389566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
139c4762a1bSJed Brown 
1409566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,U));
141c4762a1bSJed Brown   switch (var) {
142c4762a1bSJed Brown   case VAR_CONSERVATIVE:
143c4762a1bSJed Brown     break;
144c4762a1bSJed Brown   case VAR_NONCONSERVATIVE:
145c4762a1bSJed Brown   case VAR_TRANSIENTVAR:
1469566063dSJacob Faibussowitsch     PetscCall(VecExp(U));
147c4762a1bSJed Brown     break;
148c4762a1bSJed Brown   }
1499566063dSJacob Faibussowitsch   PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD));
1509566063dSJacob Faibussowitsch   PetscCall(VecSum(U,&sum));
15163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Conservation error %g\n", (double)PetscRealPart(sum - 3.)));
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
1549566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
156b122ec5aSJacob Faibussowitsch   return 0;
157c4762a1bSJed Brown }
158c4762a1bSJed Brown 
159c4762a1bSJed Brown /*TEST
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   test:
162c4762a1bSJed Brown     suffix: conservative
163c4762a1bSJed Brown     args: -snes_fd -var conservative
164c4762a1bSJed Brown   test:
165c4762a1bSJed Brown     suffix: nonconservative
166c4762a1bSJed Brown     args: -snes_fd -var nonconservative
167c4762a1bSJed Brown   test:
168c4762a1bSJed Brown     suffix: transientvar
169c4762a1bSJed Brown     args: -snes_fd -var transientvar
170c4762a1bSJed Brown 
171c4762a1bSJed Brown TEST*/
172