https://mooseframework.inl.gov
Functions | Variables
EquilibriumConstantInterpolatorTest.C File Reference

Go to the source code of this file.

Functions

 TEST (EquilibriumConstantInterpolatorTest, constructor_except)
 exception tests for constructor More...
 
 TEST (EquilibriumConstantInterpolatorTest, constructor)
 
 TEST (EquilibriumConstantInterpolatorTest, linear)
 
 TEST (EquilibriumConstantInterpolatorTest, fourthOrder)
 
 TEST (EquilibriumConstantInterpolatorTest, maierKelly)
 
 TEST (EquilibriumConstantInterpolatorTest, linearNullValues)
 
 TEST (EquilibriumConstantInterpolatorTest, fourthOrderNullValues)
 
 TEST (EquilibriumConstantInterpolatorTest, fourthOrderUserDefinedNullValues)
 
 TEST (EquilibriumConstantInterpolatorTest, piecewiselinear)
 Piecewise-linear. More...
 
 TEST (EquilibriumConstantInterpolatorTest, piecewiselinear_oneval)
 Piecewise-linear with just one value. More...
 
 TEST (EquilibriumConstantInterpolatorTest, piecewiselinear_500)
 Piecewise-linear with 500 values. More...
 
 TEST (EquilibriumConstantInterpolatorTest, piecewiselinear_except)
 Piecewise-linear exception due to badly ordered T values. More...
 

Variables

const double tol = 1.0e-6
 

Function Documentation

◆ TEST() [1/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
constructor_except   
)

exception tests for constructor

Definition at line 17 of file EquilibriumConstantInterpolatorTest.C.

18 {
19  // Fourth-order polynomial
20  std::vector<double> T = {0.0, 25.0};
21  std::vector<double> k_bad = {-6.5570};
22  std::vector<double> k_good = {-6.5570, 1.0};
23 
24  try
25  {
26  EquilibriumConstantInterpolator logk(T, k_bad, "fourth-order");
27  FAIL() << "Missing expected exception.";
28  }
29  catch (const std::exception & e)
30  {
31  std::string msg(e.what());
32  ASSERT_TRUE(msg.find("The temperature and logk data sets must be equal in length in "
33  "EquilibriumConstantInterpolator") != std::string::npos)
34  << "Failed with unexpected error message: " << msg;
35  }
36 
37  try
38  {
39  EquilibriumConstantInterpolator logk(T, k_good, "bad-type");
40  FAIL() << "Missing expected exception.";
41  }
42  catch (const std::exception & e)
43  {
44  std::string msg(e.what());
45  ASSERT_TRUE(msg.find("Type bad-type is not supported in EquilibriumConstantInterpolator") !=
46  std::string::npos)
47  << "Failed with unexpected error message: " << msg;
48  }
49 
50  try
51  {
52  EquilibriumConstantInterpolator logk(T, k_good, "maier-kelly");
53  FAIL() << "Missing expected exception.";
54  }
55  catch (const std::exception & e)
56  {
57  std::string msg(e.what());
58  ASSERT_TRUE(
59  msg.find("A Maier-Kelly fit cannot be used when the temperature points include 0. Use a "
60  "fourth-order fit instead") != std::string::npos)
61  << "Failed with unexpected error message: " << msg;
62  }
63 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...

◆ TEST() [2/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
constructor   
)

Definition at line 65 of file EquilibriumConstantInterpolatorTest.C.

66 {
67  // Fourth-order polynomial
68  std::vector<double> T = {0.0, 25.0, 60.0, 100.0, 150.0, 200.0, 250.0, 300.0};
69  std::vector<double> k = {-6.5570, -6.3660, -6.3325, -6.4330, -6.7420, -7.1880, -7.7630, -8.4650};
70 
71  EquilibriumConstantInterpolator logk(T, k, "fourth-order");
72  EXPECT_EQ(logk.getSampleSize(), T.size());
73 
74  // Maier-Kelly fit
75  T = {0.01, 25.0, 60.0, 100.0, 150.0, 200.0, 250.0, 300.0};
76  k = {-6.5570, -6.3660, -6.3325, -6.4330, -6.7420, -7.1880, -7.7630, -8.4650};
77 
78  EquilibriumConstantInterpolator logk2(T, k, "maier-kelly");
79  EXPECT_EQ(logk2.getSampleSize(), T.size());
80 
81  // Piecewise-linear fit
82  EquilibriumConstantInterpolator logk3(T, k, "piecewise-linear");
83  EXPECT_EQ(logk3.getSampleSize(), T.size());
84 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
static const std::string k
Definition: NS.h:130

◆ TEST() [3/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
linear   
)

Definition at line 86 of file EquilibriumConstantInterpolatorTest.C.

87 {
88  const std::vector<double> T = {0.0, 25.0};
89  const std::vector<double> k = {1.2, 2.4};
90 
91  // Instantiate with fourth-order fit, but a linear fit will be constructed as there
92  // are not enought points to fit a fourth-order polynomial
93  EquilibriumConstantInterpolator logk(T, k, "fourth-order");
94  EXPECT_EQ(logk.getSampleSize(), T.size());
95 
96  logk.generate();
97 
98  EXPECT_NEAR(logk.sample(T[0]), 1.2, tol);
99  EXPECT_NEAR(logk.sample(T[1]), 2.4, tol);
100  EXPECT_NEAR(logk.sample(50.0), 3.6, tol);
101 
102  EXPECT_NEAR(logk.sampleDerivative(50.0), 0.048, tol);
103 
104  // Check sampleDerivative() with the AD version
105  ADReal Tad = 50.0;
106  Moose::derivInsert(Tad.derivatives(), 0, 1.0);
107  EXPECT_NEAR(logk.sampleDerivative(50.0), logk.sample(Tad).derivatives()[0], tol);
108 
109  // Should get identical results specifying a maier-kelly fit as only a linear
110  // fit can be used (note: shift both temperature points to keep linear fit the same
111  // as above)
112  const std::vector<double> T2 = {0.01, 25.01};
113  EquilibriumConstantInterpolator logk2(T2, k, "maier-kelly");
114  EXPECT_EQ(logk2.getSampleSize(), T2.size());
115 
116  logk2.generate();
117 
118  EXPECT_NEAR(logk2.sample(T2[0]), 1.2, tol);
119  EXPECT_NEAR(logk2.sample(T2[1]), 2.4, tol);
120  EXPECT_NEAR(logk2.sample(50.01), 3.6, tol);
121 
122  EXPECT_NEAR(logk2.sampleDerivative(50.0), 0.048, tol);
123 
124  // Check sampleDerivative() with the AD version
125  ADReal Tad2 = 50.01;
126  Moose::derivInsert(Tad2.derivatives(), 0, 1.0);
127  EXPECT_NEAR(logk.sampleDerivative(50.01), logk.sample(Tad2).derivatives()[0], tol);
128 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
DualNumber< Real, DNDerivativeType, true > ADReal
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
static const std::string k
Definition: NS.h:130

◆ TEST() [4/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
fourthOrder   
)

Compare with values calculated using scipy.optimize.curve_fit

from scipy.optimize import curve_fit

def logk(T, a0, a1, a2, a3, a4): return a0 + a1 * T + a2 * T**2 + a3 * T**3 + a4 * T**4

T = [0.0, 25.0, 60.0, 100.0, 150.0, 200.0, 250.0, 300.0] k = [-6.5570, -6.3660, -6.3325, -6.4330, -6.7420, -7.1880, -7.7630, -8.4650]

popt, pcov = curve_fit(logk, T, k)

print(logk(T[1], *popt)) print(logk(T[2], *popt)) print(logk(T[3], *popt)) print(logk(T[4], *popt)) print(logk(T[5], *popt)) print(logk(125, *popt))

Definition at line 130 of file EquilibriumConstantInterpolatorTest.C.

131 {
132  const std::vector<double> T = {0.0, 25.0, 60.0, 100.0, 150.0, 200.0, 250.0, 300.0};
133  const std::vector<double> k = {
134  -6.5570, -6.3660, -6.3325, -6.4330, -6.7420, -7.1880, -7.7630, -8.4650};
135 
136  EquilibriumConstantInterpolator logk(T, k, "fourth-order");
137  logk.generate();
138 
159  EXPECT_NEAR(logk.sample(T[1]), -6.38384537894225, tol);
160  EXPECT_NEAR(logk.sample(T[2]), -6.321712110081394, tol);
161  EXPECT_NEAR(logk.sample(T[3]), -6.429646151221264, tol);
162  EXPECT_NEAR(logk.sample(T[4]), -6.744733784802997, tol);
163  EXPECT_NEAR(logk.sample(T[5]), -7.193821631774174, tol);
164  EXPECT_NEAR(logk.sample(125.0), -6.567092778770316, tol);
165 
166  // Check sampleDerivative() with the AD version
167  ADReal T2 = 50.0;
168  Moose::derivInsert(T2.derivatives(), 0, 1.0);
169  EXPECT_NEAR(logk.sampleDerivative(50.0), logk.sample(T2).derivatives()[0], tol);
170 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
DualNumber< Real, DNDerivativeType, true > ADReal
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
static const std::string k
Definition: NS.h:130

◆ TEST() [5/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
maierKelly   
)

Compare with values calculated using scipy.optimize.curve_fit

from scipy.optimize import curve_fit import numpy as np

def logk(T, a0, a1, a2, a3, a4): return a0 * np.log(T) + a1 + a2 * T + a3 / T + a4 / T**2

T = [0.01, 25.0, 60.0, 100.0, 150.0, 200.0, 250.0, 300.0] k = [0.2081, 0.0579, -0.2746, -0.7311, -1.3659, -2.0618, -2.8403, -3.7681]

popt, pcov = curve_fit(logk, T, k)

print(logk(T[1], *popt)) print(logk(T[2], *popt)) print(logk(T[3], *popt)) print(logk(T[4], *popt)) print(logk(T[5], *popt))

Definition at line 172 of file EquilibriumConstantInterpolatorTest.C.

173 {
174  const std::vector<double> T = {0.01, 25.0, 60.0, 100.0, 150.0, 200.0, 250.0, 300.0};
175  const std::vector<double> k = {
176  0.2081, 0.0579, -0.2746, -0.7311, -1.3659, -2.0618, -2.8403, -3.7681};
177 
178  EquilibriumConstantInterpolator logk(T, k, "maier-kelly");
179  logk.generate();
180 
201  EXPECT_NEAR(logk.sample(T[1]), 0.061497433011205596, tol);
202  EXPECT_NEAR(logk.sample(T[2]), -0.30083706818706746, tol);
203  EXPECT_NEAR(logk.sample(T[3]), -0.7005921220757013, tol);
204  EXPECT_NEAR(logk.sample(T[4]), -1.3418011239865, tol);
205  EXPECT_NEAR(logk.sample(T[5]), -2.0827775371250423, tol);
206  EXPECT_NEAR(logk.sample(125.0), -1.0053697469019622, tol);
207 
208  // Check sampleDerivative() with the AD version
209  ADReal T2 = 50.0;
210  Moose::derivInsert(T2.derivatives(), 0, 1.0);
211  EXPECT_NEAR(logk.sampleDerivative(50.0), logk.sample(T2).derivatives()[0], tol);
212 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
DualNumber< Real, DNDerivativeType, true > ADReal
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
static const std::string k
Definition: NS.h:130

◆ TEST() [6/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
linearNullValues   
)

Definition at line 214 of file EquilibriumConstantInterpolatorTest.C.

215 {
216  const std::vector<double> T = {0.0, 25.0, 60.0, 100.0, 150.0, 200.0, 250.0, 300.0};
217  const std::vector<double> k = {17.8660, 15.9746, 500.0, 500.0, 500.0, 500.0, 500.0, 500.0};
218 
219  EquilibriumConstantInterpolator logk(T, k, "fourth-order");
220  logk.generate();
221 
222  // Should only be two samples in this case (ignoring values of 500)
223  EXPECT_EQ(logk.getSampleSize(), (std::size_t)2);
224 
225  // Compare with values calculated using linear fit and extrapolating
226  EXPECT_NEAR(logk.sample(T[0]), 17.866, tol);
227  EXPECT_NEAR(logk.sample(T[1]), 15.9746, tol);
228  EXPECT_NEAR(logk.sample(60.0), 13.32664, tol);
229  EXPECT_NEAR(logk.sample(125.0), 8.409, tol);
230  EXPECT_NEAR(logk.sample(300.0), -4.8308, tol);
231 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
static const std::string k
Definition: NS.h:130

◆ TEST() [7/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
fourthOrderNullValues   
)

Definition at line 233 of file EquilibriumConstantInterpolatorTest.C.

234 {
235  const std::vector<double> T = {0.0, 25.0, 60.0, 100.0, 150.0, 200.0, 250.0, 300.0};
236  const std::vector<double> k = {
237  -0.0162, -0.0264, -0.0340, -0.0459, -0.0618, -0.0803, 500.0, 500.0};
238 
239  EquilibriumConstantInterpolator logk(T, k, "fourth-order");
240  logk.generate();
241 
242  // Should only be six samples in this case (ignoring values of 500)
243  EXPECT_EQ(logk.getSampleSize(), (std::size_t)6);
244 
245  // Compare with values calculated using fourth order polynomial (python code above)
246  EXPECT_NEAR(logk.sample(T[0]), -0.016418224486741174, tol);
247  EXPECT_NEAR(logk.sample(T[1]), -0.025715932916029558, tol);
248  EXPECT_NEAR(logk.sample(T[2]), -0.03492782535210691, tol);
249  EXPECT_NEAR(logk.sample(T[3]), -0.04524532642418216, tol);
250  EXPECT_NEAR(logk.sample(T[4]), -0.0620327728353172, tol);
251  EXPECT_NEAR(logk.sample(250.0), -0.08739490878982198, tol);
252  EXPECT_NEAR(logk.sample(300.0), -0.060966149023565674, tol);
253 
254  // Check sampleDerivative() with the AD version
255  ADReal T2 = 50.0;
256  Moose::derivInsert(T2.derivatives(), 0, 1.0);
257  EXPECT_NEAR(logk.sampleDerivative(50.0), logk.sample(T2).derivatives()[0], tol);
258 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
DualNumber< Real, DNDerivativeType, true > ADReal
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
static const std::string k
Definition: NS.h:130

◆ TEST() [8/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
fourthOrderUserDefinedNullValues   
)

Definition at line 260 of file EquilibriumConstantInterpolatorTest.C.

261 {
262  const std::vector<double> T = {0.0, 25.0, 60.0, 100.0, 150.0, 200.0, 250.0, 300.0};
263  const std::vector<double> k = {
264  -0.0162, -0.0264, -0.0340, -0.0459, -0.0618, -0.0803, 999.999, 999.999};
265 
266  EquilibriumConstantInterpolator logk(T, k, "fourth-order", 999.999);
267  logk.generate();
268 
269  // Should only be six samples in this case (ignoring values of 500)
270  EXPECT_EQ(logk.getSampleSize(), (std::size_t)6);
271 
272  // Compare with values calculated using fourth order polynomial (python code above)
273  EXPECT_NEAR(logk.sample(T[0]), -0.016418224486741174, tol);
274  EXPECT_NEAR(logk.sample(T[1]), -0.025715932916029558, tol);
275  EXPECT_NEAR(logk.sample(T[2]), -0.03492782535210691, tol);
276  EXPECT_NEAR(logk.sample(T[3]), -0.04524532642418216, tol);
277  EXPECT_NEAR(logk.sample(T[4]), -0.0620327728353172, tol);
278  EXPECT_NEAR(logk.sample(250.0), -0.08739490878982198, tol);
279  EXPECT_NEAR(logk.sample(300.0), -0.060966149023565674, tol);
280 
281  // Check sampleDerivative() with the AD version
282  ADReal T2 = 50.0;
283  Moose::derivInsert(T2.derivatives(), 0, 1.0);
284  EXPECT_NEAR(logk.sampleDerivative(50.0), logk.sample(T2).derivatives()[0], tol);
285 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
DualNumber< Real, DNDerivativeType, true > ADReal
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
static const std::string k
Definition: NS.h:130

◆ TEST() [9/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
piecewiselinear   
)

Piecewise-linear.

Definition at line 288 of file EquilibriumConstantInterpolatorTest.C.

289 {
290  const std::vector<double> T = {0.0, 2.0, 6.0};
291  const std::vector<double> k = {10.0, 0.0, 40.0};
292 
293  EquilibriumConstantInterpolator logk(T, k, "piecewise-linear");
294  EXPECT_EQ(logk.getSampleSize(), T.size());
295 
296  for (unsigned i = 0; i < 3; ++i)
297  EXPECT_NEAR(logk.sample(T[i]), k[i], tol);
298  EXPECT_NEAR(logk.sample(1.0), 5.0, tol);
299  EXPECT_NEAR(logk.sample(3.0), 10.0, tol);
300  EXPECT_NEAR(logk.sample(-1.0), 10.0, tol);
301  EXPECT_NEAR(logk.sample(600), 40.0, tol);
302 
303  EXPECT_EQ(logk.sampleDerivative(-1.0), 0.0);
304  EXPECT_NEAR(logk.sampleDerivative(1.0), -5.0, tol);
305  EXPECT_NEAR(logk.sampleDerivative(3.0), 10.0, tol);
306  EXPECT_EQ(logk.sampleDerivative(600.0), 0.0);
307 
308  // Check sampleDerivative() with the AD version produces the correct error
309  ADReal T2 = 5.0;
310  Moose::derivInsert(T2.derivatives(), 0, 1.0);
311  try
312  {
313  EXPECT_NEAR(logk.sampleDerivative(5.0), logk.sample(T2).derivatives()[0], tol);
314  FAIL() << "Missing expected exception.";
315  }
316  catch (const std::exception & e)
317  {
318  std::string msg(e.what());
319  ASSERT_TRUE(
320  msg.find("Dual cannot be used for specified fit type in EquilibriumConstantInterpolator") !=
321  std::string::npos)
322  << "Failed with unexpected error message: " << msg;
323  }
324 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
DualNumber< Real, DNDerivativeType, true > ADReal
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
static const std::string k
Definition: NS.h:130

◆ TEST() [10/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
piecewiselinear_oneval   
)

Piecewise-linear with just one value.

Definition at line 327 of file EquilibriumConstantInterpolatorTest.C.

328 {
329  const std::vector<double> T = {2.0};
330  const std::vector<double> k = {10.0};
331 
332  EquilibriumConstantInterpolator logk(T, k, "piecewise-linear");
333  EXPECT_EQ(logk.getSampleSize(), T.size());
334 
335  EXPECT_EQ(logk.sample(1.0), 10.0);
336  EXPECT_EQ(logk.sample(3.0), 10.0);
337 
338  EXPECT_EQ(logk.sampleDerivative(-1.0), 0.0);
339  EXPECT_EQ(logk.sampleDerivative(600.0), 0.0);
340 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
static const std::string k
Definition: NS.h:130

◆ TEST() [11/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
piecewiselinear_500   
)

Piecewise-linear with 500 values.

Definition at line 343 of file EquilibriumConstantInterpolatorTest.C.

344 {
345  const std::vector<double> T = {0.0, 2.0, 4.0, 6.0, 8.0};
346  const std::vector<double> k = {10.0, 0.0, 500.0, 40.0, 500.0};
347 
348  EquilibriumConstantInterpolator logk(T, k, "piecewise-linear");
349  EXPECT_EQ(logk.getSampleSize(), (std::size_t)3);
350 
351  logk.generate(); // does nothing relevant because type = piecewise-linear
352 
353  EXPECT_NEAR(logk.sample(T[0]), k[0], tol);
354  EXPECT_NEAR(logk.sample(T[1]), k[1], tol);
355  EXPECT_NEAR(logk.sample(T[3]), k[3], tol);
356  EXPECT_NEAR(logk.sample(1.0), 5.0, tol);
357  EXPECT_NEAR(logk.sample(3.0), 10.0, tol);
358  EXPECT_NEAR(logk.sample(-1.0), 10.0, tol);
359  EXPECT_NEAR(logk.sample(600), 40.0, tol);
360 
361  EXPECT_EQ(logk.sampleDerivative(-1.0), 0.0);
362  EXPECT_NEAR(logk.sampleDerivative(1.0), -5.0, tol);
363  EXPECT_NEAR(logk.sampleDerivative(3.0), 10.0, tol);
364  EXPECT_EQ(logk.sampleDerivative(600.0), 0.0);
365 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
static const std::string k
Definition: NS.h:130

◆ TEST() [12/12]

TEST ( EquilibriumConstantInterpolatorTest  ,
piecewiselinear_except   
)

Piecewise-linear exception due to badly ordered T values.

Definition at line 368 of file EquilibriumConstantInterpolatorTest.C.

369 {
370  const std::vector<double> T = {2.0, 0.0, 6.0};
371  const std::vector<double> k = {10.0, 0.0, 40.0};
372 
373  try
374  {
375 
376  EquilibriumConstantInterpolator logk(T, k, "piecewise-linear");
377  FAIL() << "Missing expected exception.";
378  }
379  catch (const std::exception & e)
380  {
381  std::string msg(e.what());
382  ASSERT_TRUE(msg.find("EquilibriumConstantInterpolation: x-values are not strictly increasing: "
383  "x[0]: 2 x[1]: 0") != std::string::npos)
384  << "Failed with unexpected error message: " << msg;
385  }
386 }
Fit the equilibrium constant values read from the thermodynamic databse at specified temperature valu...
static const std::string k
Definition: NS.h:130

Variable Documentation

◆ tol

const double tol = 1.0e-6

Definition at line 14 of file EquilibriumConstantInterpolatorTest.C.

Referenced by TEST().