60 pressureIdx = Indices::pressureIdx,
61 switchIdx = Indices::switchIdx,
66 contiH2OEqIdx = Indices::conti0EqIdx + FluidSystem::H2OIdx,
67 contiO2EqIdx = Indices::conti0EqIdx + FluidSystem::O2Idx,
73 liquidPhaseIdx = FluidSystem::liquidPhaseIdx,
74 gasPhaseIdx = FluidSystem::gasPhaseIdx
77 using GridView =
typename GridGeometry::GridView;
79 using GlobalPosition =
typename Dune::FieldVector<typename GridView::ctype, GridView::dimensionworld>;
91 template<
class SourceValues>
93 const std::string& paramGroup =
"")
102 const auto lengthBox = gridYMax/nCellsY;
104 currentDensity *= 2.0/lengthBox;
106 currentDensity *= 1.0/lengthBox;
108 static Scalar transportNumberH2O =
getParam<Scalar>(
"ElectroChemistry.TransportNumberH20");
111 values[contiH2OEqIdx] = currentDensity/(2*
Constant::F);
112 values[contiH2OEqIdx] += currentDensity/
Constant::F*transportNumberH2O;
113 values[contiO2EqIdx] = -currentDensity/(4*
Constant::F);
121 template<
class VolumeVariables>
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");
131 Scalar currentDensity = reversibleVoltage - cellVoltage - 0.5;
132 Scalar increment = 1e-4;
133 Scalar deltaCurrentDensity = currentDensity*increment;
134 Scalar deltaVoltage = 1.0;
139 while (abs(deltaVoltage) > 1e-6)
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);
151 deltaVoltage = currentDensity*specificResistance - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
153 currentDensity = currentDensity - (deltaVoltage*deltaCurrentDensity)/(deltaCurrentDensity*specificResistance + activationLossesDiff);
155 activationLosses = calculateActivationLosses_(volVars, currentDensity);
157 deltaVoltage = currentDensity*specificResistance - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
165 deltaVoltage = currentDensity*specificResistance/(1-sw) - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
167 currentDensity = currentDensity - (deltaVoltage*deltaCurrentDensity)/(deltaCurrentDensity*specificResistance/(1-sw) + activationLossesDiff);
169 activationLosses = calculateActivationLosses_(volVars, currentDensity);
171 deltaVoltage = currentDensity*specificResistance/(1-sw) - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
176 if(iterations >= maxIter)
178 DUNE_THROW(
NumericalProblem,
"Newton solver for electrochemistry didn't converge");
183 return currentDensity*10000;
193 template<
class VolumeVariables>
194 static Scalar calculateActivationLosses_(
const VolumeVariables &volVars,
const Scalar currentDensity)
196 static Scalar refO2PartialPressure =
getParam<Scalar>(
"ElectroChemistry.RefO2PartialPressure");
197 static Scalar numElectrons =
getParam<Scalar>(
"ElectroChemistry.NumElectrons");
198 static Scalar transferCoefficient =
getParam<Scalar>(
"ElectroChemistry.TransferCoefficient");
201 Scalar sw = volVars.saturation(liquidPhaseIdx);
203 Scalar preFactor =
Constant::R*volVars.fluidState().temperature()/transferCoefficient/
Constant::F/numElectrons;
205 Scalar pO2 = volVars.pressure(gasPhaseIdx) * volVars.fluidState().moleFraction(gasPhaseIdx, FluidSystem::O2Idx);
214 *( log(abs(currentDensity)/abs(exchangeCurrentDensity_(volVars)))
215 - log(pO2/refO2PartialPressure)
222 *( log(abs(currentDensity)/abs(exchangeCurrentDensity_(volVars)))
223 - log(pO2/refO2PartialPressure)
234 template<
class VolumeVariables>
235 static Scalar calculateConcentrationLosses_(
const VolumeVariables &volVars)
238 static Scalar numElectrons =
getParam<Scalar>(
"ElectroChemistry.NumElectrons");
239 static Scalar transferCoefficient =
getParam<Scalar>(
"ElectroChemistry.TransferCoefficient");
244 Scalar pO2 = volVars.pressure(gasPhaseIdx) * volVars.fluidState().moleFraction(gasPhaseIdx, FluidSystem::O2Idx);
251 losses = -1.0*preFactor*(transferCoefficient/2)*log(pO2/pO2Inlet);
255 losses = -1.0*preFactor*(transferCoefficient/2+1)*log(pO2/pO2Inlet);
266 template<
class VolumeVariables>
267 static Scalar exchangeCurrentDensity_(
const VolumeVariables &volVars)
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");
275 Scalar T = volVars.fluidState().temperature();
276 Scalar refExchangeCurrentDensity = -1.0
278 * surfaceIncreasingFactor
279 * exp(-1.0 * activationBarrier /
Constant::R * (1/T-1/refTemperature));
281 return refExchangeCurrentDensity;
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:375
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition parameters.hh:365
static void reactionSource(SourceValues &values, Scalar currentDensity, const std::string ¶mGroup="")
Calculates reaction sources with an electrochemical model approach.
Definition electrochemistry.hh:92