Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 3 : //* 4 : //* All rights reserved, see COPYRIGHT for full restrictions 5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT 6 : //* 7 : //* Licensed under LGPL 2.1, please see LICENSE for details 8 : //* https://www.gnu.org/licenses/lgpl-2.1.html 9 : 10 : #include "EquilibriumConstantAux.h" 11 : 12 : registerMooseObject("ChemicalReactionsApp", EquilibriumConstantAux); 13 : 14 : InputParameters 15 123 : EquilibriumConstantAux::validParams() 16 : { 17 123 : InputParameters params = AuxKernel::validParams(); 18 246 : params.addCoupledVar( 19 : "temperature", 298.15, "The temperature of the aqueous phase (K). Default is 298.15K"); 20 246 : params.addRequiredParam<std::vector<Real>>( 21 : "temperature_points", "Temperature points where log(Keq) data is evaluated (K)"); 22 246 : params.addRequiredParam<std::vector<Real>>( 23 : "logk_points", "log(Keq) data evaluated at each value of temperature_points"); 24 123 : params.addClassDescription( 25 : "Equilibrium constant for a given equilibrium species (in form log10(Keq))"); 26 123 : return params; 27 0 : } 28 : 29 66 : EquilibriumConstantAux::EquilibriumConstantAux(const InputParameters & parameters) 30 : : AuxKernel(parameters), 31 66 : _temperature(coupledValue("temperature")), 32 132 : _temperature_points(getParam<std::vector<Real>>("temperature_points")), 33 198 : _logk_points(getParam<std::vector<Real>>("logk_points")) 34 : { 35 : // Check that the number of temperature_points and logk_points are equal 36 66 : if (_temperature_points.size() != _logk_points.size()) 37 0 : mooseError("The number of temperature_points and logk_points must be equal in ", _name); 38 : 39 66 : if (_temperature_points.size() >= 5) 40 : { 41 : // If there at least 5 values, then use the Maier-Kelley fit 42 44 : _logk = std::make_unique<EquilibriumConstantFit>(_temperature_points, _logk_points); 43 22 : _logk->generate(); 44 : } 45 44 : else if ((_temperature_points.size() >= 2) && (_temperature_points.size() <= 4)) 46 : { 47 : // If between 2 and 4 values are provided, use a linear fit 48 44 : _linear_logk = std::make_unique<PolynomialFit>(_temperature_points, _logk_points, 1); 49 22 : _linear_logk->generate(); 50 : } 51 66 : } 52 : 53 : Real 54 6732 : EquilibriumConstantAux::computeValue() 55 : { 56 6732 : if (_temperature_points.size() == 1) 57 2244 : return -_logk_points[0]; 58 : 59 4488 : if (_temperature_points.size() > 5) 60 2244 : return -_logk->sample(_temperature[_qp]); 61 : 62 2244 : return -_linear_logk->sample(_temperature[_qp]); 63 : }