48        pressureIdx = Indices::pressureIdx, 
 
   49        switchIdx = Indices::switchIdx,     
 
   54        contiH2OEqIdx = Indices::conti0EqIdx + FluidSystem::H2OIdx,
 
   55        contiO2EqIdx = Indices::conti0EqIdx + FluidSystem::O2Idx,
 
   61        liquidPhaseIdx = FluidSystem::liquidPhaseIdx,
 
   62        gasPhaseIdx = FluidSystem::gasPhaseIdx
 
   65    using GridView = 
typename GridGeometry::GridView;
 
   67    using GlobalPosition = 
typename Dune::FieldVector<typename GridView::ctype, GridView::dimensionworld>;
 
   79    template<
class SourceValues>
 
   81                               const std::string& paramGroup = 
"")
 
   90        const auto lengthBox = gridYMax/nCellsY;
 
   92            currentDensity *= 2.0/lengthBox;
 
   94            currentDensity *= 1.0/lengthBox;
 
   96        static Scalar transportNumberH2O = 
getParam<Scalar>(
"ElectroChemistry.TransportNumberH20");
 
   99        values[contiH2OEqIdx] = currentDensity/(2*
Constant::F);                  
 
  100        values[contiH2OEqIdx] += currentDensity/
Constant::F*transportNumberH2O;  
 
  101        values[contiO2EqIdx]  = -currentDensity/(4*
Constant::F);                 
 
 
  109    template<
class VolumeVariables>
 
  113        static int maxIter = 
getParam<int>(
"ElectroChemistry.MaxIterations");
 
  114        static Scalar specificResistance = 
getParam<Scalar>(
"ElectroChemistry.SpecificResistance");
 
  115        static Scalar reversibleVoltage = 
getParam<Scalar>(
"ElectroChemistry.ReversibleVoltage");
 
  116        static Scalar cellVoltage = 
getParam<Scalar>(
"ElectroChemistry.CellVoltage");
 
  119        Scalar currentDensity = reversibleVoltage - cellVoltage - 0.5;
 
  120        Scalar increment = 1e-4;
 
  121        Scalar deltaCurrentDensity = currentDensity*increment;
 
  122        Scalar deltaVoltage = 1.0;
 
  127        while (abs(deltaVoltage) > 1e-6)
 
  130            Scalar activationLosses        = calculateActivationLosses_(volVars, currentDensity);
 
  131            Scalar activationLossesNext    = calculateActivationLosses_(volVars, currentDensity+deltaCurrentDensity);
 
  132            Scalar concentrationLosses     = calculateConcentrationLosses_(volVars);
 
  133            Scalar activationLossesDiff    = activationLossesNext - activationLosses;
 
  134            Scalar sw                      = volVars.saturation(liquidPhaseIdx);
 
  139                deltaVoltage = currentDensity*specificResistance - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
 
  141                currentDensity = currentDensity - (deltaVoltage*deltaCurrentDensity)/(deltaCurrentDensity*specificResistance + activationLossesDiff);
 
  143                activationLosses = calculateActivationLosses_(volVars, currentDensity);
 
  145                deltaVoltage = currentDensity*specificResistance - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
 
  153                deltaVoltage = currentDensity*specificResistance/(1-sw) - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
 
  155                currentDensity = currentDensity - (deltaVoltage*deltaCurrentDensity)/(deltaCurrentDensity*specificResistance/(1-sw) + activationLossesDiff);
 
  157                activationLosses = calculateActivationLosses_(volVars, currentDensity);
 
  159                deltaVoltage = currentDensity*specificResistance/(1-sw) - reversibleVoltage + cellVoltage + activationLosses + concentrationLosses;
 
  164            if(iterations >= maxIter)
 
  166                DUNE_THROW(
NumericalProblem, 
"Newton solver for electrochemistry didn't converge");
 
  171        return currentDensity*10000;
 
 
  181    template<
class VolumeVariables>
 
  182    static Scalar calculateActivationLosses_(
const VolumeVariables &volVars, 
const Scalar currentDensity)
 
  184        static Scalar refO2PartialPressure = 
getParam<Scalar>(
"ElectroChemistry.RefO2PartialPressure");
 
  185        static Scalar numElectrons = 
getParam<Scalar>(
"ElectroChemistry.NumElectrons");
 
  186        static Scalar transferCoefficient = 
getParam<Scalar>(
"ElectroChemistry.TransferCoefficient");
 
  189        Scalar sw = volVars.saturation(liquidPhaseIdx);
 
  191        Scalar preFactor = 
Constant::R*volVars.fluidState().temperature()/transferCoefficient/
Constant::F/numElectrons;
 
  193        Scalar pO2 = volVars.pressure(gasPhaseIdx) * volVars.fluidState().moleFraction(gasPhaseIdx, FluidSystem::O2Idx);
 
  202                            *(  log(abs(currentDensity)/abs(exchangeCurrentDensity_(volVars)))
 
  203                            - log(pO2/refO2PartialPressure)
 
  210            *(  log(abs(currentDensity)/abs(exchangeCurrentDensity_(volVars)))
 
  211                - log(pO2/refO2PartialPressure)
 
  222    template<
class VolumeVariables>
 
  223    static Scalar calculateConcentrationLosses_(
const VolumeVariables &volVars)
 
  226        static Scalar numElectrons = 
getParam<Scalar>(
"ElectroChemistry.NumElectrons");
 
  227        static Scalar transferCoefficient = 
getParam<Scalar>(
"ElectroChemistry.TransferCoefficient");
 
  232        Scalar pO2 = volVars.pressure(gasPhaseIdx) * volVars.fluidState().moleFraction(gasPhaseIdx, FluidSystem::O2Idx);
 
  239            losses = -1.0*preFactor*(transferCoefficient/2)*log(pO2/pO2Inlet);
 
  243            losses = -1.0*preFactor*(transferCoefficient/2+1)*log(pO2/pO2Inlet);
 
  254    template<
class VolumeVariables>
 
  255    static Scalar exchangeCurrentDensity_(
const VolumeVariables &volVars)
 
  258        static Scalar activationBarrier = 
getParam<Scalar>(
"ElectroChemistry.ActivationBarrier");
 
  259        static Scalar surfaceIncreasingFactor = 
getParam<Scalar>(
"ElectroChemistry.SurfaceIncreasingFactor");
 
  260        static Scalar refTemperature = 
getParam<Scalar>(
"ElectroChemistry.RefTemperature");
 
  261        static Scalar refCurrentDensity = 
getParam<Scalar>(
"ElectroChemistry.RefCurrentDensity");
 
  263        Scalar T = volVars.fluidState().temperature();
 
  264        Scalar refExchangeCurrentDensity = -1.0
 
  266                            * surfaceIncreasingFactor
 
  267                            * exp(-1.0 * activationBarrier / 
Constant::R * (1/T-1/refTemperature));
 
  269        return refExchangeCurrentDensity;
 
 
static void reactionSource(SourceValues &values, Scalar currentDensity, const std::string ¶mGroup="")
Calculates reaction sources with an electrochemical model approach.
Definition electrochemistry.hh:80
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition parameters.hh:139