3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
electrochemistry.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_ELECTROCHEMISTRY_HH
25#define DUMUX_ELECTROCHEMISTRY_HH
26
27#include <cmath>
28
33
35
36namespace Dumux {
37
43
51template <class Scalar, class Indices, class FluidSystem, class GridGeometry, ElectroChemistryModel electroChemistryModel>
53{
54
56
57 enum
58 {
59 //indices of the primary variables
60 pressureIdx = Indices::pressureIdx, // gas-phase pressure
61 switchIdx = Indices::switchIdx, // liquid saturation or mole fraction
62 };
63 enum
64 {
65 //equation indices
66 contiH2OEqIdx = Indices::conti0EqIdx + FluidSystem::H2OIdx,
67 contiO2EqIdx = Indices::conti0EqIdx + FluidSystem::O2Idx,
68 };
69
70 enum
71 {
72 // phase indices
73 liquidPhaseIdx = FluidSystem::liquidPhaseIdx,
74 gasPhaseIdx = FluidSystem::gasPhaseIdx
75 };
76
77 using GridView = typename GridGeometry::GridView;
78 static constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethod::box;
79 using GlobalPosition = typename Dune::FieldVector<typename GridView::ctype, GridView::dimensionworld>;
80
81public:
91 template<class SourceValues>
92 static void reactionSource(SourceValues &values, Scalar currentDensity,
93 const std::string& paramGroup = "")
94 {
95 // correction to account for actually relevant reaction area
96 // current density has to be divided by the half length of the box
97 // \todo Do we have multiply with the electrochemically active surface area (ECSA) here instead?
98 static Scalar gridYMax = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.UpperRight")[1];
99 static Scalar nCellsY = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.Cells")[1];
100
101 // Warning: This assumes the reaction layer is always just one cell (cell-centered) or half a box (box) thick
102 const auto lengthBox = gridYMax/nCellsY;
103 if (isBox)
104 currentDensity *= 2.0/lengthBox;
105 else
106 currentDensity *= 1.0/lengthBox;
107
108 static Scalar transportNumberH2O = getParam<Scalar>("ElectroChemistry.TransportNumberH20");
109
110 //calculation of flux terms with Faraday equation
111 values[contiH2OEqIdx] = currentDensity/(2*Constant::F); //reaction term in reaction layer
112 values[contiH2OEqIdx] += currentDensity/Constant::F*transportNumberH2O; //osmotic term in membrane
113 values[contiO2EqIdx] = -currentDensity/(4*Constant::F); //O2-equation
114 }
115
121 template<class VolumeVariables>
122 static Scalar calculateCurrentDensity(const VolumeVariables &volVars)
123 {
124
125 static int maxIter = getParam<int>("ElectroChemistry.MaxIterations");
126 static Scalar specificResistance = getParam<Scalar>("ElectroChemistry.SpecificResistance");
127 static Scalar reversibleVoltage = getParam<Scalar>("ElectroChemistry.ReversibleVoltage");
128 static Scalar cellVoltage = getParam<Scalar>("ElectroChemistry.CellVoltage");
129
130 //initial guess for the current density and initial newton solver parameters
131 Scalar currentDensity = reversibleVoltage - cellVoltage - 0.5;
132 Scalar increment = 1e-4;
133 Scalar deltaCurrentDensity = currentDensity*increment;
134 Scalar deltaVoltage = 1.0;
135 int iterations = 0;
136
137 //Newton Solver for current Density
138 using std::abs;
139 while (abs(deltaVoltage) > 1e-6)
140 {
141
142 Scalar activationLosses = calculateActivationLosses_(volVars, currentDensity);
143 Scalar activationLossesNext = calculateActivationLosses_(volVars, currentDensity+deltaCurrentDensity);
144 Scalar concentrationLosses = calculateConcentrationLosses_(volVars);
145 Scalar activationLossesDiff = activationLossesNext - activationLosses;
146 Scalar sw = volVars.saturation(liquidPhaseIdx);
147
148 if(electroChemistryModel == ElectroChemistryModel::Acosta)
149 {
150 // Acosta calculation
151 deltaVoltage = currentDensity*specificResistance - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
152
153 currentDensity = currentDensity - (deltaVoltage*deltaCurrentDensity)/(deltaCurrentDensity*specificResistance + activationLossesDiff);
154
155 activationLosses = calculateActivationLosses_(volVars, currentDensity);
156
157 deltaVoltage = currentDensity*specificResistance - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
158
159 iterations++;
160
161 }
162 else
163 {
164 // Ochs & other calculation
165 deltaVoltage = currentDensity*specificResistance/(1-sw) - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
166
167 currentDensity = currentDensity - (deltaVoltage*deltaCurrentDensity)/(deltaCurrentDensity*specificResistance/(1-sw) + activationLossesDiff);
168
169 activationLosses = calculateActivationLosses_(volVars, currentDensity);
170
171 deltaVoltage = currentDensity*specificResistance/(1-sw) - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
172
173 iterations++;
174 }
175
176 if(iterations >= maxIter)
177 {
178 DUNE_THROW(NumericalProblem, "Newton solver for electrochemistry didn't converge");
179 }
180 }
181
182 //conversion from [A/cm^2] to [A/m^2]
183 return currentDensity*10000;
184 }
185
186private:
187
193 template<class VolumeVariables>
194 static Scalar calculateActivationLosses_(const VolumeVariables &volVars, const Scalar currentDensity)
195 {
196 static Scalar refO2PartialPressure = getParam<Scalar>("ElectroChemistry.RefO2PartialPressure");
197 static Scalar numElectrons = getParam<Scalar>("ElectroChemistry.NumElectrons");
198 static Scalar transferCoefficient = getParam<Scalar>("ElectroChemistry.TransferCoefficient");
199
200 //Saturation sw for Acosta calculation
201 Scalar sw = volVars.saturation(liquidPhaseIdx);
202 //Calculate prefactor
203 Scalar preFactor = Constant::R*volVars.fluidState().temperature()/transferCoefficient/Constant::F/numElectrons;
204 //Get partial pressure of O2 in the gas phase
205 Scalar pO2 = volVars.pressure(gasPhaseIdx) * volVars.fluidState().moleFraction(gasPhaseIdx, FluidSystem::O2Idx);
206
207 Scalar losses = 0.0;
208 //Calculate activation losses
209 using std::log;
210 using std::abs;
211 if(electroChemistryModel == ElectroChemistryModel::Acosta)
212 {
213 losses = preFactor
214 *( log(abs(currentDensity)/abs(exchangeCurrentDensity_(volVars)))
215 - log(pO2/refO2PartialPressure)
216 - log(1 - sw)
217 );
218 }
219 else
220 {
221 losses = preFactor
222 *( log(abs(currentDensity)/abs(exchangeCurrentDensity_(volVars)))
223 - log(pO2/refO2PartialPressure)
224 );
225 }
226 return losses;
227 }
228
229
234 template<class VolumeVariables>
235 static Scalar calculateConcentrationLosses_(const VolumeVariables &volVars)
236 {
237 static Scalar pO2Inlet = getParam<Scalar>("ElectroChemistry.pO2Inlet");
238 static Scalar numElectrons = getParam<Scalar>("ElectroChemistry.NumElectrons");
239 static Scalar transferCoefficient = getParam<Scalar>("ElectroChemistry.TransferCoefficient");
240
241 //Calculate preFactor
242 Scalar preFactor = Constant::R*volVars.temperature()/transferCoefficient/Constant::F/numElectrons;
243 //Get partial pressure of O2 in the gas phase
244 Scalar pO2 = volVars.pressure(gasPhaseIdx) * volVars.fluidState().moleFraction(gasPhaseIdx, FluidSystem::O2Idx);
245
246 Scalar losses = 0.0;
247 //Calculate concentration losses
248 using std::log;
249 if(electroChemistryModel == ElectroChemistryModel::Acosta)
250 {
251 losses = -1.0*preFactor*(transferCoefficient/2)*log(pO2/pO2Inlet);
252 }else
253 {
254 // +1 is the Nernst part of the equation
255 losses = -1.0*preFactor*(transferCoefficient/2+1)*log(pO2/pO2Inlet);
256 }
257
258 return losses;
259 }
260
261
266 template<class VolumeVariables>
267 static Scalar exchangeCurrentDensity_(const VolumeVariables &volVars)
268 {
269 using std::exp;
270 static Scalar activationBarrier = getParam<Scalar>("ElectroChemistry.ActivationBarrier");
271 static Scalar surfaceIncreasingFactor = getParam<Scalar>("ElectroChemistry.SurfaceIncreasingFactor");
272 static Scalar refTemperature = getParam<Scalar>("ElectroChemistry.RefTemperature");
273 static Scalar refCurrentDensity = getParam<Scalar>("ElectroChemistry.RefCurrentDensity");
274
275 Scalar T = volVars.fluidState().temperature();
276 Scalar refExchangeCurrentDensity = -1.0
277 * refCurrentDensity
278 * surfaceIncreasingFactor
279 * exp(-1.0 * activationBarrier / Constant::R * (1/T-1/refTemperature));
280
281 return refExchangeCurrentDensity;
282 }
283};
284
285} // end namespace Dumux
286
287#endif
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
A two-phase (water and air) fluid system with water, nitrogen and oxygen as components.
A central place for various physical constants occuring in some equations.
ElectroChemistryModel
The type of electrochemistry models implemented here (Ochs (2008) or Acosta et al....
Definition: electrochemistry.hh:42
@ Acosta
Definition: electrochemistry.hh:42
@ Ochs
Definition: electrochemistry.hh:42
Definition: adapt.hh:29
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
This class calculates source terms and current densities for fuel cells with the electrochemical mode...
Definition: electrochemistry.hh:53
static void reactionSource(SourceValues &values, Scalar currentDensity, const std::string &paramGroup="")
Calculates reaction sources with an electrochemical model approach.
Definition: electrochemistry.hh:92
static Scalar calculateCurrentDensity(const VolumeVariables &volVars)
Newton solver for calculation of the current density.
Definition: electrochemistry.hh:122
A central place for various physical constants occuring in some equations.
Definition: constants.hh:39
static constexpr Scalar F
Faraday constant .
Definition: constants.hh:66
static constexpr Scalar R
The ideal gas constant .
Definition: constants.hh:44