![]() |
SLIDE
3.0.0
A simulator for lithium-ion battery pack degradation
|
#include <Cell_SPM.hpp>
Public Types | |
using | sigma_type = std::array< double, settings::nch+2 > |
Public Member Functions | |
Cell_SPM (OCVcurves OCV_curves_) | |
< Constructor More... | |
Cell_SPM (std::string IDi, const DEG_ID °id, double capf, double resf, double degfsei, double degflam) | |
Cell_SPM () | |
Default constructor. More... | |
Cell_SPM (Model_SPM *M_ptr) | |
getters More... | |
double | T () noexcept override |
returns the uniform battery temperature in [K] More... | |
double | getTenv () const noexcept |
get the environmental temperature [K] More... | |
Status | setCurrent (double Inew, bool checkV=true, bool print=true) override |
Status | setSOC (double SOCnew, bool checkV=true, bool print=true) override |
auto & | getStateObj () |
auto | setStateObj (State_SPM &st_new) |
std::array< double, 4 > | getVariations () const noexcept override |
void | getTemperatures (double *Tenv, double *Tref) noexcept |
< get the environmental and reference temperature More... | |
double | getThermalSurface () override |
return the 'A' for the thermal model of this SU (Q = hA*dT) More... | |
double | thermalModel (int Nneighb, double Tneighb[], double Kneighb[], double Aneighb[], double tim) override |
double | thermal_getTotalHeat () |
function for unit testing More... | |
double | getRdc () noexcept |
calculate the resistance from every component More... | |
double | getRtot () override |
double | getAnodeSurface () noexcept |
get the anode pure surface area (without cracks) product of the effective surface area (an) with the electrode volume More... | |
double | I () const override |
get the cell current [A] positive for discharging More... | |
double | V () override |
print is an optional argument More... | |
bool | getCSurf (double &cps, double &cns, bool print) |
get the surface concentrations More... | |
void | getC (double cp[], double cn[]) noexcept |
get the concentrations at all nodes More... | |
int | getVoltage (bool print, double *V, double *OCVp, double *OCVn, double *etap, double *etan, double *Rdrop, double *Temp) |
get the cell's voltage More... | |
int | getVoltage_ne (bool print, double *V, double *OCVp, double *OCVn, double *etap, double *etan, double *Rdrop, double *Temp) |
get the cell's voltage noexcept More... | |
void | getDaiStress (double *sigma_p, double *sigma_n, sigma_type &sigma_r_p, sigma_type &sigma_r_n, sigma_type &sigma_t_p, sigma_type &sigma_t_n, sigma_type &sigma_h_p, sigma_type &sigma_h_n) noexcept |
get the stresses at all nodes according to Dai's stress model More... | |
void | updateDaiStress () noexcept |
updated the stored stress values for Dai's stress model More... | |
void | getLaresgoitiStress (bool print, double *sigma_n) |
get the stresses at all nodes according to Laresgoiti's stress model More... | |
void | updateLaresgoitiStress (bool print) |
update the stored stress values for Laresgoiti's stress model More... | |
void | setT (double T) override |
set the cell's temperature More... | |
void | setTenv (double Tenv) |
void setStates(const State_SPM &si, double I); //!< set the cell's states to the states in the State object and the cell current to the given value More... | |
void | setC (double cp0, double cn0) |
void setCurrent(bool critical, bool check, double I); //!< set the cell's current to the specified value -> From old slide. More... | |
void | peekVoltage (double I) |
Peeks voltage for state for given I. More... | |
ThroughputData | getThroughputs () override |
void | overwriteCharacterisationStates (double Dpi, double Dni, double ri) |
void | overwriteGeometricStates (double thickpi, double thickni, double epi, double eni, double api, double ani) |
void ETI_electr(bool print, double I, double dti, bool blockDegradation, bool pos); //!< step forward with only one electrode using forward Euler time integration More... | |
State_SPM::states_type | dState (bool print, State_SPM &d_state) |
Utility. More... | |
void | checkModelparam () |
check if the inputs to the MATLAB code are the same as the ones here in the C++ code More... | |
void | getStates (getStates_t s) override |
returns the states of the cell collectively. More... | |
std::span< double > | viewStates () override |
Only for cells to see individual states. More... | |
double | getOCV () override |
Status | setStates (setStates_t sSpan, bool checkV, bool print) override |
opposite of getStates, check the states are valid? More... | |
bool | validStates (bool print=true) override |
checks if a state array is valid More... | |
double | SOC () override |
void | timeStep_CC (double dt, int steps=1) override |
setStates(std::move(states)); //!< store new states, checks if they are illegal (throws an error in that case) More... | |
Cell_SPM * | copy () override |
Obsolete functions (do not use): More... | |
void | ETI (bool print, double dti, bool blockDegradation) |
step forward in time using forward Eurler time integration More... | |
void | setOCVcurve (const std::string &namepos, const std::string &nameneg) |
sets the OCV curve of the cell to the given value More... | |
void | setInitialConcentration (double cmaxp, double cmaxn, double lifracp, double lifracn) |
sets the initial concentration More... | |
void | setGeometricParameters (double capnom, double elec_surf, double ep, double en, double thickp, double thickn) |
void setRamping(double Istep, double tstep); //!< sets the ramping parameters More... | |
void | setCharacterisationParam (double Dp, double Dn, double kp, double kn, double Rdc) |
sets the parameters related to the characterisation of the cell More... | |
![]() | |
Cell () | |
Cell (const std::string &ID_) | |
virtual | ~Cell ()=default |
double | Cap () const final override |
void | setCapacity (double capacity) |
constexpr double | Vmin () const override |
lower voltage limit of the cell More... | |
constexpr double | VMIN () const override |
safety cut off More... | |
constexpr double | VMAX () const override |
safety cut off More... | |
constexpr double | Vmax () const override |
upper voltage limit of the cell More... | |
constexpr double | Tmax () |
constexpr double | Tmin () |
double | getVhigh () final |
return the voltage of the cell with the highest voltage More... | |
double | getVlow () final |
return the voltage of the cell with the lowest voltage More... | |
virtual Status | setSOC (double SOCnew, bool checkV=true, bool print=true)=0 |
virtual double | SOC ()=0 |
virtual double | getThotSpot () override |
the T of the hottest element in the SU More... | |
size_t | getNcells () override final |
this is a single cell More... | |
virtual Status | checkCurrent (bool checkV, bool print) noexcept |
virtual Status | checkVoltage (double &v, bool print) noexcept override |
Check the voltage status of the cell. More... | |
virtual double | thermalModel (int Nneighb, double Tneighb[], double Kneighb[], double Aneighb[], double tim) override |
Calculate the thermal model of the cell. More... | |
virtual void | storeData () override |
Add another data point in the array. More... | |
virtual void | writeData (const std::string &prefix) override |
virtual ThroughputData | getThroughputs () |
virtual std::array< double, 4 > | getVariations () const noexcept |
![]() | |
StorageUnit ()=default | |
< basic getters and setters More... | |
StorageUnit (std::string_view ID_) | |
StorageUnit (std::string_view ID_, StorageUnit *parent_, bool blockDegAndTherm_) | |
virtual | ~StorageUnit ()=default |
const std::string & | getID () |
void | setID (std::string IDi) |
Return the full ID string, including the ID of the parent module. More... | |
virtual std::string | getFullID () |
virtual double | Cap () const =0 |
auto * | getParent () |
virtual double | I () const =0 |
virtual double | getRtot ()=0 |
virtual size_t | getNcells ()=0 |
return the number of single cells connected to this SU More... | |
bool | isCharging () |
negative means charge. More... | |
bool | isDischarging () |
positive means discharge. More... | |
virtual void | getStates (getStates_t s)=0 |
returns one long array with the states More... | |
virtual viewStates_t | viewStates () |
Only for cells to see individual states. More... | |
void | setBlockDegAndTherm (bool block) |
virtual void | setParent (StorageUnit *p) |
set the parent More... | |
virtual Status | setCurrent (double Inew, bool checkV=true, bool print=true)=0 |
virtual Status | setVoltage (double Vnew, bool checkI=true, bool print=true) |
virtual Status | setStates (setStates_t s, bool checkStates=true, bool print=true)=0 |
opposite of getStates, check the states are valid? More... | |
virtual void | backupStates () |
Back-up states. More... | |
virtual void | restoreStates () |
restore backed-up states. More... | |
virtual double | getOCV ()=0 |
virtual double | V ()=0 |
print is an optional argument More... | |
virtual Status | checkVoltage (double &v, bool print) noexcept=0 |
get the voltage and check if it is valid More... | |
virtual double | getVhigh ()=0 |
return the voltage of the cell with the highest voltage More... | |
virtual double | getVlow ()=0 |
return the voltage of the cell with the lowest voltage More... | |
virtual double | Vmin () const =0 |
lower voltage limit of the cell More... | |
virtual double | VMIN () const =0 |
safety cut off More... | |
virtual double | Vmax () const =0 |
upper voltage limit of the cell More... | |
virtual double | VMAX () const =0 |
safety cut off More... | |
virtual double | T ()=0 |
virtual double | getThotSpot ()=0 |
the T of the hottest element in the SU More... | |
virtual double | getThermalSurface ()=0 |
return the 'A' for the thermal model of this SU (Q = hA*dT) More... | |
virtual double | thermalModel (int Nneighb, double Tneighb[], double Kneighb[], double Aneighb[], double tim)=0 |
calculate the thermal model of this SU More... | |
virtual void | setT (double Tnew)=0 |
functionality More... | |
virtual bool | validStates (bool print=true)=0 |
checks if a state array is valid More... | |
virtual StorageUnit * | copy ()=0 |
copy this SU to a new object More... | |
virtual void | timeStep_CC (double dt, int steps=1)=0 |
take a number of time steps More... | |
virtual void | storeData ()=0 |
virtual void | writeData (const std::string &prefix)=0 |
Public Attributes | |
DEG_ID | deg_id |
structure with the identification of which degradation model(s) to use #TODO may be protected. More... | |
Protected Member Functions | |
std::pair< double, double > | calcSurfaceConcentration (double jp, double jn, double Dpt, double Dnt) |
std::pair< double, double > | calcOverPotential (double cps, double cns, double i_app) |
Should not throw normally, except divide by zero? More... | |
double | calcArrheniusCoeff () |
Calculates Arrhenius coefficient. More... | |
std::array< double, 2 > | calcDiffConstant () |
Calculate the diffusion constant at the battery temperature using an Arrhenius relation. More... | |
std::array< double, 3 > | calcMolarFlux () |
Calculate molar flux. More... | |
void | SEI (double OCVnt, double etan, double *isei, double *den) |
calculate the effect of SEI growth More... | |
void | CS (double OCVnt, double etan, double *isei_multiplyer, double *dCS, double *dDn) |
calculate the effect of surface crack growth More... | |
void | LAM (bool critical, double zp_surf, double etap, double *dthickp, double *dthickn, double *dap, double *dan, double *dep, double *den) |
calculate the effect of LAM More... | |
double | LiPlating (double OCVnt, double etan) |
state space model More... | |
void | dState_diffusion (bool print, State_SPM &d_state) |
just diffusion PDE More... | |
void | dState_thermal (bool print, double &dQgen) |
calculate the heat generation More... | |
void | dState_degradation (bool print, State_SPM &d_state) |
calculate the effect of lithium plating More... | |
void | dState_all (bool print, State_SPM &d_state) |
individual functions are combined in one function to gain speed. More... | |
double | thermalModel_cell () |
double | thermalModel_coupled (int Nneighb, double Tneighb[], double Kneighb[], double Aneighb[], double tim) |
cell to cell variations // #TODO why do we store variations? More... | |
![]() | |
virtual size_t | calculateNcells () |
Protected Attributes | |
State_SPM | st {} |
< protected such that child classes can access the class variables More... | |
State_SPM | s_ini {} |
the battery current/initial state, grouping all parameter which change over the battery's lifetime (see State_SPM.hpp) More... | |
double | Cmaxpos { 51385 } |
maximum lithium concentration in the cathode [mol m-3] value for NMC More... | |
double | Cmaxneg { 30555 } |
maximum lithium concentration in the anode [mol m-3] value for C More... | |
double | C_elec { 1000 } |
Li- concentration in electrolyte [mol m-3] standard concentration of 1 molar. More... | |
double | n { 1 } |
number of electrons involved in the main reaction [-] #TODO if really constant? More... | |
double | kp { 5e-11 } |
rate constant of main reaction at positive electrode at reference temperature More... | |
double | kp_T { 58000 } |
activation energy for the Arrhenius relation of kp More... | |
double | kn { 1.7640e-11 } |
rate constant of main reaction at negative electrode at reference temperature More... | |
double | kn_T { 20000 } |
The diffusion constants at reference temperature are part of State because they can change over the battery's lifetime. More... | |
double | Dp_T { 29000 } |
activation energy for the Arrhenius relation of Dp More... | |
double | Dn_T { 35000.0 / 5.0 } |
activation energy for the Arrhenius relation of Dn More... | |
double | Therm_Qgen {} |
total heat generation since the last update [J] More... | |
double | Therm_Qgentot {} |
variable for unit testing, total heat generation since the beginning of this cell's life [J] More... | |
double | Therm_time {} |
time since the last update of the thermal model More... | |
double | T_env { settings::T_ENV } |
environment temperature [K] More... | |
double | T_ref { 25.0_degC } |
reference temperature [K] More... | |
double | Qch { 45 } |
convective heat transfer coefficient per volume [W K-1 m-3] More... | |
double | rho { 1626 } |
density of the battery More... | |
double | Cp { 750 } |
thermal capacity of the battery #TODO = units missing. More... | |
param::Geometry_SPM | geo {} |
other geometric parameters are part of State because they can change over the battery's lifetime More... | |
param::StressParam | sparam { param::def::StressParam_Kokam } |
Stress parameters. More... | |
double | rsei { 2037.4 * 50 } |
specific resistance times real surface area of the SEI film [Ohm m] ? #TODO if unit is correct. aging fit. More... | |
double | nsei { 1 } |
number of electrons involved in the SEI reaction [-] More... | |
double | alphasei { 1 } |
charge transfer coefficient of the SEI reaction [-] More... | |
double | OCVsei { 0.4 } |
equilibrium potential of the SEI side reaction [V] More... | |
double | rhosei { 100e3 } |
partial molar volume of the SEI layer [m3 mol-1] More... | |
double | c_elec0 { 4.541e-3 } |
bulk concentration of the electrolyte molecule participating in the SEI growth (e.g. EC) [mol m-3] More... | |
double | Vmain { 13 } |
partial molar volume of the main reaction, see Ashwin et al, 2016 More... | |
double | Vsei { 64.39 } |
partial molar volume of the SEI side reaction, see Ashwin et al., 2016 More... | |
param::SEIparam | sei_p { param::def::SEIparam_Kokam } |
structure with the fitting parameters of the different SEI growth models More... | |
param::CSparam | csparam {} |
structure with the fitting parameters of the different crack growth models More... | |
double | OCVnmc { 4.1 } |
equilibrium potential of the NMC dissolution side reaction [V] More... | |
param::LAMparam | lam_p { param::def::LAMparam_Kokam } |
structure with the fitting parameters of the different LAM models More... | |
double | npl { 1 } |
number of electrons involved in the plating reaction [-] More... | |
double | alphapl { 1 } |
charge transfer constant for the plating reaction [-] More... | |
double | OCVpl { 0 } |
OCV of the plating reaction [V]. More... | |
double | rhopl { 10e6 } |
density of the plated lithium layer More... | |
param::PLparam | pl_p |
structure with the fitting parameters of the different plating models More... | |
Model_SPM * | M { Model_SPM::makeModel() } |
OCV curves. More... | |
OCVcurves | OCV_curves |
bool | Vcell_valid { false } |
Functions. More... | |
double | var_cap { 1 } |
relative factor increasing the capacity of the cell More... | |
double | var_R { 1 } |
relative factor increasing the DC resistance More... | |
double | var_degSEI { 1 } |
relative factor to speed up or slow down the rate of SEI growth More... | |
double | var_degLAM { 1 } |
relative factor to speed up or slow down the rate of LAM More... | |
![]() | |
double | capNom { 16 } |
capacity [Ah]. More... | |
CellData< settings::DATASTORE_CELL > | cellData |
Cell data storage. More... | |
![]() | |
std::string | ID { "StorageUnit" } |
identification string More... | |
StorageUnit * | parent { nullptr } |
pointer to the SU 'above' this one [e.g. the module to which a cell is connected] More... | |
bool | blockDegAndTherm { false } |
if true, degradation and the thermal ODE are ignored More... | |
Additional Inherited Members | |
![]() | |
static constexpr CellLimits | limits { defaultCellLimits } |
![]() | |
using | setStates_t = std::span< double > & |
To pass states to read, non-expandable container. More... | |
using | getStates_t = std::vector< double > & |
To pass states to save, expandable container. More... | |
using | viewStates_t = std::span< double > |
using slide::Cell_SPM::sigma_type = std::array<double, settings::nch + 2> |
|
inline |
< Constructor
slide::Cell_SPM::Cell_SPM | ( | std::string | IDi, |
const DEG_ID & | degid, | ||
double | capf, | ||
double | resf, | ||
double | degfsei, | ||
double | degflam | ||
) |
< Set the resistance
< current collector
< cathode
< anode #TODO if you memoize Rdc then getRdc here.
< set the capacity
< nominal capacity
< surface area of the electrodes (current to current density)
< set the degradation factor
< set the degradation ID and related settings
< Check if we will have to calculate the stress according to Dai's stress model
< check if we need to calculate the stress according to Laresgoiti's stress model
slide::Cell_SPM::Cell_SPM | ( | ) |
Default constructor.
< ID string
< Parameters are given for 16 Ah high-power, prismatic KokamNMC cell. (SLPB78205130H)
< Set initial state:
< SEI thickness. Start with a fresh cell, which has undergone some formation cycles so it has an initial SEI layer. never start with a value of 0, because some equations have a term 1/delta, which would give nan or inf so this will give errors in the code
< lost lithium. Start with 0 so we can keep track of how much li we lose while cycling the cell
< diffusion constant of the cathode at reference temperature
< diffusion constant of the anode at reference temperature
< thickness of the positive electrode
< thickness of the negative electrode
< volume fraction of active material in the cathode
< volume fraction of active material in the anode
< effective surface area of the cathode, the 'real' surface area is the product of the effective surface area (a) with the electrode volume (elec_surf * thick)
< effective surface area of the anode
< initial crack surface. Start with 1% of the real surface area
< R = Rdc * (thickp * ap * elec_surf + thickn * an * elec_surf) / 2; //!< specific resistance of the combined electrodes, see State::iniStates
< note: this is divided by geometric surface (elec_surf) instead of effective surf (a*thick*elec_surf)
< thickness of the plated lithium layer. You can start with 0 here
< 0.689332 lithium fraction in the cathode at 50% soc (3.68136 V) [-]
< 0.479283 lithium fraction in the anode at 50% soc (3.68136 V) [-]
< Since fp and fn are set at 50%.
< set the states, with a random value for the concentration
< Default values for not defined other param:
< assume the maximum crack surface is 5 times the initial anode surface
|
inline |
getters
|
inlineprotected |
Calculates Arrhenius coefficient.
|
protected |
Calculate the diffusion constant at the battery temperature using an Arrhenius relation.
< Calculate the diffusion constant at the battery temperature using an Arrhenius relation
< diffusion constant of the positive particle [m s-1]
< diffusion constant of the negative particle [m s-1]
|
protected |
Calculate molar flux.
void setStates(State_SPM &&states); //!< set the cell's states to the states in the array degradation models
< current density on the electrode [A m-2]
< molar flux on the positive particle [mol m-2 s-1]
< molar flux on the negative particle [mol m-2 s-1]
|
protected |
Should not throw normally, except divide by zero?
< Calculate the rate constants at the cell's temperature using an Arrhenius relation
< Rate constant at the positive electrode at the cell's temperature [m s-1]
< Rate constant at the negative electrode at the cell's temperature [m s-1]
< Calculate the overpotential using the Bulter-Volmer equation if alpha is 0.5, the Bulter-Volmer relation can be inverted to eta = 2RT / (nF) asinh(x) and asinh(x) = ln(x + sqrt(1+x^2) -> to asinh(x) function.
< exchange current density of the positive electrode
< exchange current density of the negative electrode
< x for the cathode
< x for the anode
< cathode overpotential [V], < 0 on discharge
< anode overpotential [V], > 0 on discharge
i_app | Should not throw normally, except divide by zero? |
|
protected |
< Calculate the surface concentration at the positive particle cp_surf = M->Cp[0][:] * zp[:] + M->Dp*jp/Dpt
< Calculate the surface concentration at the negative particle cn_surf = M->Cn[0][:] * zn[:] + M->Dn*jn/Dnt
Dnt | Should not throw normally, except divide by zero? |
void slide::Cell_SPM::checkModelparam | ( | ) |
check if the inputs to the MATLAB code are the same as the ones here in the C++ code
— slide-pack functions — //
< check if the inputs to the MATLAB code are the same as the ones here in the C++ code input: M->Input[0] has to be the same as nch (defined in State.hpp) M->Input[1] has to be the same as Rp (defined earlier in this constructor) M->Input[2] has to be the same as Rn (defined earlier in this constructor) M->input[3] has to give the location of the 0 eigenvalue
< allow a relative difference of e-10 due to numerical errors
< allow a relative difference of e-10 due to numerical errors
< allow a relative difference of e-10 due to numerical errors
< allow a relative difference of e-10 due to numerical errors
|
inlineoverridevirtual |
Obsolete functions (do not use):
Implements slide::StorageUnit.
|
protected |
calculate the effect of surface crack growth
< output parameters
< isei multiplier from all models to be considered
< increase in surface area from all models to be considered
< active surface area of the anode [m2] this active area is used to translate absolute values (such as currents) to relative values (such as current densities)
< Loop for each model we want to use
< a switch to calculate the effect according to model i
< no surface cracks
< Laresgoiti's stress and crack growth model (Laresgoiti, Kablitz, Ecker, Sauer, Journal of Power Sources 300, 2015) this model calculates crack growth due to temporal variations in the li-concentration
< check the calculated stress values are up to date
< if you see this error, you have to call Cell_SPM::updateLaresgoitiStress(), which calculates the stress and stores the values, before you call this function
< equations (22)+ (27) from the paper current density on particle = I /(elec_surf * thick * a) isei also acts on this scale since it is an extra boundary condition ( itot = (jn + isei) = (surface gradient)/nF ) crack growth -> increase isei_on_particle -> (jn + isei*(initial+crack_surface)/initial)*nF such that if the crack surface area is the same as the initial electrode surface area, we double isei
< increase SEI growth proportionally the crack surface
< Laresgoiti's crack growth model (Laresgoiti, Kablitz, Ecker, Sauer, Journal of Power Sources 300, 2015) but with stress model from Dai, Cai, White, Journal of Power sources 247, 2014 instead of Laresgoiti's stress from figure 5 this model calculates crack growth due to spatial (Dai) and temporal (Laresgoiti) variations in the li-concentration
< ensure the stress values are up to date
< if you see this error, you have to call Cell_SPM::updateDaiStress(), which calculates the stress and stores the values before you call this function
< equations (22)+ (27) from the paper but with Dai's stress
< increase SEI growth proportionally the crack surface
< model by Deshpande & Bernardi,Journal of the Electrochemical Society 164 (2), 2017 this model is adapted to calculate crack growth due to spatial variations in the li-concentration
< get concentrations
< Add the effects of this model
< equations (8) + (21) Note that eqn (8) refers to the change with respect to time while here we use the spatial variation This makes the model capture more of the spatial variation in stress Laresgoiti's model already accounted for temporal variation, so simply use that if you are interested in temporal rather than spatial variation
< increase SEI growth proportionally the crack surface
< model from Barai, Smith, Chen, Kim, Mukherjee, Journal of the Electrochemical Society 162 (9), 2015 equation (1a): CS = Amax(1 - exp(-m * Ah)) with m a fitting parameter and Ah the charge throughput up to now dCS/dAh = m Amax exp(-m Ah) = m Amax - m Amax + m Amax exp(-m Ah) = m Amax - m CS = m(Amax - CS) dCS/dt = dCS/dAh * dAh/dt = dCS/dAh * abs(I) = m(Amax - CS)*abs(I) where we use the absolute value of I because 'Ah' is the total charge throughput, i.e. int ( abs(I) dt )
< 'maximum crack surface area', a fitting parameters avoid negative crack growth if the crack surface becomes slightly larger than Amax this is possible due to discrete time steps: CS(t) is just smaller, but CS (t+1) = CS(t) + dCS*dt is just larger
< Add the effects of this model
< see above, with m = csparam.CS4
< increase SEI growth proportionally the crack surface
< model from Ekstrom and Lindbergh, Journal of the Electrochemical Society 162 (6), 2015 overpotential for the crack-side-reaction = overpotential for the SEI reaction
< overpotential [V], equation (6)
< get surface concentration
< get the surface lithium concentration //!< #TODO only cns is used.
< rate constant for the side reaction
< Calculate the rate constant, equation (11) with an Arrhenius relation for the temperature (which wasn't considered by Ekstrom)
< Add the effects of this model
< equation (9)
< increase SEI growth proportionally the crack surface
< unknown degradation model
< Decrease the negative diffusion constant if needed
< don't decrease the negative diffusion constant
< decrease it according to Barai, Smith, Chen, Kim, Mukherjee, Journal of the Electrochemical Society 162 (9), 2015
< equation (2) D(t) = D0 (1 - CS)^gamma but this can become negative is CS is larger than 1, we assume there should be a term /Amax in eqn (2): D(t) = D0 (1 - (CS/Amax))^gamma, which becomes 0 if CS = Amax so dD/dt = - gamma D0 (1-CS/Amax)^(gamma-1) 1/Amax dCS/dt
< 'maximum crack surface area', a fitting parameters avoid increasing diffusion coefficient if the crack surface becomes larger than Amax this is possible if the user chooses a different CS growth model, which does give larger crack surfaces (i.e. not CS4 which is Barai's crack growth model)
< cap the decrease rate at a maximum value of 2e-7 (2e-5% = kill the battery in about 1000h)
< unknown degradation model
State_SPM::states_type slide::Cell_SPM::dState | ( | bool | print, |
State_SPM & | d_state | ||
) |
Utility.
|
protected |
individual functions are combined in one function to gain speed.
thermal model
|
protected |
calculate the effect of lithium plating
< Calculcate the lithium fractions at the surface of the particles
< current density, molar flux on the pos/neg particle
< check if the surface concentration is within the allowed range 0 < cp < Cmaxpos 0 < cn < Cmaxneg
< Do not delete if you cannot ensure zp/zn between 0-1
< lithium fraction (0 to 1)
< Calculate the overpotentials if needed
< calculate the anode potential (needed for various degradation models)
< entropic coefficient of the anode potential [V K-1]
< anode potential [V]
< anode potential at the cell's temperature [V]
< SEI growth
< current density of the SEI growth side reaction [A m-2]
< decrease in volume fraction due to SEI growth [s-1]
< additional diffusion in the anode due to isei
< Throws but not wrapped in try-catch since only appears here.
< Subtract Li from negative electrode (like an extra current density -> add it in the boundary conditions: dCn/dx = jn + isei/nF)
< crack growth leading to additional exposed surface area
< relative increase in isei due to additional SEI growth on the extra exposed surface area [-]
< increase in crack surface area [m2 s-1]
< change in negative diffusion constant [m s-1 s-1]
< additional diffusion in the anode due to extra SEI growth on the crack surface
< Throws but not wrapped in try-catch since only appears here.
< crack surface leads to extra SEI growth because the exposed surface area increases. (like an extra current density -> add it in the boundary conditions: dCn/dx = jn + isei/nF + isei_CS/nF)
< extra SEI side reaction current density due to the crack surface area [A m-2]
< loss of active material LAM
< change in geometric parameters describing the amount of active material
< Throws but not wrapped in try-catch since only appears here.
< lithium plating
< current density of the plating side reaction [A m-2]
< additional diffusion in the anode due to ipl
< Subtract Li from negative electrode (like an extra current density -> add it in the boundary conditions: dCn/dx = jn + ipl/nF)
< #TODO (npl * F) division is unnecessary it is already multiple in the function.
< output
< dzp += 0 //!< dzp should be added from diffusion function
< dzn jtot = jn + isei/nF + isei_CS/nF + ipl/nF
< ddelta SEI thickness
< dLLI lost lithium
< dthickp electrode thickness
< dthickn
< dep volume fraction of active material
< den
< dap effective surface are, a = 3 e/R
< dan
< dCS surface area of the cracks
< dDp diffusion constant
< dDn
< drdc_p cathode resistance
< drdc_n anode resistance
< drdc_cc current collector resistance
< ddelta_pl thickness of the plated lithium
|
protected |
just diffusion PDE
< current density, molar flux on the pos/neg particle
< Calculate the effect of the main li-reaction on the (transformed) concentration
< dz/dt = D * A * z + B * j
< loop for each row of the matrix-vector product A * z
< A is diagonal, so the array M->A has only the diagonal elements
< dz/dt = D * A * z + B * j
< dSOC state of charge
|
protected |
calculate the heat generation
< Calculcate the lithium fractions at the surface of the particles
< current density, molar flux on the pos/neg particle
< check if the surface concentration is within the allowed range 0 < cp < Cmaxpos 0 < cn < Cmaxneg
< Do not delete if you cannot ensure zp/zn between 0-1
< lithium fraction (0 to 1)
< const double zn_surf = (cns / Cmaxneg);
< Calculate the overpotentials if needed
< #TODO why don't we use in the new one?
< only consider negative electrode, ignore the positive electrode
< only consider positive electrode, ignore the negative electrode
< in linear interpolation, throw an error if you are outside of the allowed range of the data
< entropic coefficient of the entire cell OCV [V K-1]
< temperature model Calculate the thermal sources/sinks/transfers per unit of volume of the battery The battery volume is given by the product of the cell thickness and the electrode surface
< reversible heat due to entropy changes [W]
< reaction heat due to the kinetics [W]
< Ohmic heat due to electrode resistance [W]
< const double Qc = -Qch * SAV * (st.T() - T_env); //!< cooling with the environment [W m-2]
< total heat generation and cooling in W
< dQcool = (Qc) * (L*elec_surf);
< dT = 1/(rho*Cp)*(Qrev+Qrea+Qohm+Qc); //!< dT cell temperature
void slide::Cell_SPM::ETI | ( | bool | print, |
double | dti, | ||
bool | blockDegradation | ||
) |
step forward in time using forward Eurler time integration
Functions for fitting:
|
inlinenoexcept |
get the anode pure surface area (without cracks) product of the effective surface area (an) with the electrode volume
|
noexcept |
get the concentrations at all nodes
< Calculate the diffusion constant at the battery temperature using an Arrhenius relation
< Calculate the molar flux on the surfaces
< current density, molar flux on the pos/neg particle
< Calculate concentration at the surface and inner nodes using the matrices from the spatial discretisation of the solid diffusion PDE cp = M->Cp[:][:] * zp[:] + M->Dp*jp/Dpt cn = M->Cn[:][:] * zn[:] + M->Dn*jn/Dnt
< Problem here!!!!!!! #TODO
< loop to calculate at each surface + inner node
bool slide::Cell_SPM::getCSurf | ( | double & | cps, |
double & | cns, | ||
bool | |||
) |
get the surface concentrations
< current density, molar flux on the pos/neg particle
< check if the surface concentration is within the allowed range 0 < cp < Cmaxpos && 0 < cn < Cmaxneg
|
noexcept |
get the stresses at all nodes according to Dai's stress model
< Get the locations of the Chebyshev nodes
< location (x-value) of the positive Chebyshev nodes
< location (x-value) of the positive and negative Chebyshev nodes [-surface .. centre .. +surface]
< centre node
< negative nodes
< positive nodes
< positive and negative nodes, [+surface .. inner .. centre]
< concentrations at all nodes, [-surface .. inner .. centre .. inner .. +surface]
< cathode centre node
< anode centre node
< cathode negative points
< anode negative points
< cathode positive points
< anode positive points
< Calculate the matrix-vector product of the required functions as given in the paper by Dai et al. (concentration * radius^2) we need to remember the transformation of variables from x (-1<x<1) to r (-R<r<R) int( c * r^2 dr) = int(c * (x*R)^2 * d(R*x)) = int(c x^2 R^3 dx) so the matrix-vector product we need is F = Q * (concentration * x^2 * R^3)
< Note: since Fp (Fn) is multiplied by Rp^3 (Rn^3) then divided by them they are eliminated to remove numerical problems.
< arrays with the product for the positive and negative electrode
< loop for each row (one row = one node)
< calculate the matrix-vector product for row i as you would do it by hand: F(i) = sum(Q(i,j)*C(j)*x(j)^2*R^3, j=0..2*nch+2)
< loop through the columns to calculate the sum
< Q(i,j)*C(j)*x(j)^2
< int( cp*r^2, r=0..Rp ) // Note r^3 is eliminated!
< int( cn*r^2, r=0..Rn )
< Reset the stress values.
< loop for the positive nodes
< r(i) = R * x(i) radius of positive node i in the positive particle radius of positive node i in the negative particle
< integral from the centre to positive node i int(cp*zp^2, zp=0..rp(i)) = F[nch+1+i] - F[nch+1]
< integral from the centre to positive node i int(cn*zn^2, zn=0..rn(i))
< Implement the equations from Dai et al. Flip all arrays to get the opposite order (now it is [centre .. +surface] and we want [+surface .. centre] and store in the output arrays
< centre node -> special formula (31 & 33) in Dai, Cai, White
< other nodes -> equation 13 in Dai, Cai, White
< ap = int (c x^2, x=0..R), bp = int (c x^2 , x=0..r)
< calculate hydrostatic stress
< find the maximum (in absolute value) of the stresses
< node with the maximum hydrostatic stress in the positive/negative particle
void slide::Cell_SPM::getLaresgoitiStress | ( | bool | print, |
double * | sigma_n | ||
) |
get the stresses at all nodes according to Laresgoiti's stress model
< Arrays with the stress from Laresgoiti, Kablitz, Ecker, Sauer, Journal of Power Sources 300, 2015 (figure 5)
< li-fraction in the graphite
< stress [MPa]
< Change this if xx is not fixed time step.
< Get the surface concentration
< get the surface lithium concentration [mol m-3]
< lithium fraction on negative surface [0, 1]
< check if the surface concentration is within the allowed range 0 < cp < Cmaxpos 0 < cn < Cmaxneg
|
overridevirtual |
< print if the (global) verbose-setting is above the threshold
< Get the surface concentrations
< Calculate the li-fraction (instead of the li-concentration)
< in linear interpolation, throw an error if you are out of the allowed range
< entropic coefficient of the total cell voltage [V/K]
< anode potential [V]
< cathode potential [V]
Implements slide::StorageUnit.
|
noexcept |
calculate the resistance from every component
< real surface area of the sei-layer
< real surface area of the positive electrode
< real surface area of the negative electrode
< calculate the resistance from every component
|
inlineoverridevirtual |
|
inline |
|
inlineoverridevirtual |
returns the states of the cell collectively.
Implements slide::StorageUnit.
|
inlinenoexcept |
< get the environmental and reference temperature
thermal model
|
inlinenoexcept |
get the environmental temperature [K]
|
overridevirtual |
return the 'A' for the thermal model of this SU (Q = hA*dT)
Implements slide::StorageUnit.
|
inlineoverridevirtual |
|
inlineoverridevirtualnoexcept |
Reimplemented from slide::Cell.
int slide::Cell_SPM::getVoltage | ( | bool | print, |
double * | V, | ||
double * | OCVp, | ||
double * | OCVn, | ||
double * | etap, | ||
double * | etan, | ||
double * | Rdrop, | ||
double * | Temp | ||
) |
get the cell's voltage
int slide::Cell_SPM::getVoltage_ne | ( | bool | print, |
double * | V, | ||
double * | OCVp, | ||
double * | OCVn, | ||
double * | etap, | ||
double * | etan, | ||
double * | Rdrop, | ||
double * | Temp | ||
) |
get the cell's voltage noexcept
|
inlineoverridevirtual |
get the cell current [A] positive for discharging
Implements slide::StorageUnit.
|
protected |
calculate the effect of LAM
< output parameters
< loop for each model to use
< calculate the effect of this model
< no LAM
< Stress model from Dai, Cai, White, Journal of Power sources 247, 2014 LAM equation similar to CS equation from Laresgoiti (Laresgoiti, Kablitz, Ecker, Sauer, Journal of Power Sources 300, 2015)
< ensure the stress values are up to date
< if you see this error, you have to call Cell_SPM::updateDaiStress(), which calculates the stress and stores the values before calling this function
< with Dai's stress model (values stored in s_dai_p)
< #TODO sparam.s_dt was 2.0 in slide, why?
< ageing fit you need to divide by the time step to counter the effect of the time step larger dt has double effect: increase the stress difference, and time integrate by larger time period assuming stress changes are constant at ds per second, it would be ds (time_now - time_prev), or ds * dt and to cover a period of T (so we need T/dt steps of dt each), the total effect of stress is (ds*dt) * dt * T/dt = ds*dt*T so the effect of stress increases linearly with the size of the time setp so divide by total time (dt*nstep) to cancel this out, i.e. ds is per second (and not per total time period) assume the other effects are 0
< Model by Delacourt & Safari, Journal of the Electrochemical Society 159 (8), 2012
< Get the molar flux on each particle
< current density, molar flux on the pos/neg particle
< Use Arrhenius relations to update the fitting parameters for the cell temperature
< Add the effects of this model
< equation (5) from the paper
< assume the other effects are 0
< Model by Kindermann, Keil, Frank, Jossen, Journal of the Electrochemical Society 164 (12), 2017
< cathode potential
< get OCV of positive electrode, throw error if out of bounds this should be updated for the cell's temperature using the entropic coefficient of the cathode but I couldn't find any data on this, so I have ignored the effect
< std::cout << "Throw test: " << 40 << '
';
< equation (9) from the paper
< temperature dependent rate constant
< Arrhenius law
< current density of the NMC dissolution reaction
< equation (8) from the paper
< cap the effect at 5e-6 to avoid a very fast drop of capacity (which could cause an error) a value of 5e-6 gives dap = -3.5. The initial value is about 17000, so the cell is dead in 5,000 seconds so this cap is quite high
< Add the effects of this model
< assume the other effects are 0
< Model by Narayanrao, Joglekar, Inguva, Journal of the Electrochemical Society 160 (1), 2012
< Add the effects of this model
< equation (7) from the paper
< assume the other effects are 0
< unknown degradation model
|
protected |
state space model
< Arrhenius relation for temperature-dependent plating parameters
< Rate constant
< Yang, Leng, Zhang, Ge, Wang, Journal of Power Sources 360, 2017
|
inline |
< Overwrite both current and initial states.
|
inline |
void ETI_electr(bool print, double I, double dti, bool blockDegradation, bool pos); //!< step forward with only one electrode using forward Euler time integration
time integration void integratorStep(bool print, double dti, bool blockDegradation); //!< step forward in time using specified numerical integrator Base integrators: State_SPM::states_type Int_FWEuler(bool print, double dti, bool blockDegradation); //!< Forward euler integrator. State_SPM::states_type Int_RK4(bool print, double dti, bool blockDegradation); //!< Runge-Kutta 4 integrator. Calculate the time derivatives of the states at the actual cell current (state-space model)
< Overwrite both current and initial states.
void slide::Cell_SPM::peekVoltage | ( | double | I | ) |
Peeks voltage for state for given I.
State related functions void validState() { slide::validState(st, s_ini); }
|
protected |
calculate the effect of SEI growth
< variables
< SEI side reaction current density of all models combined
< Loop for each model to use
< Use a switch to calculate the magnitude of the SEI growth according to this degradation model
< no SEI growth
< kinetics and diffusion according to Pinson & Bazant, Journal of the Electrochemical society 160 (2), 2013
< Arrhenius relation for the rate parameter at the cell temperature
< Add the effect of this model eta_sei = OCVneg + etaneg - OCVsei + rsei*I isei = nFk exp(-nF/RT alpha eta_sei) on charge, I < 0 and etan < 0. so higher charging current -> more negative term in exponential -> larger isei
< Kinetic model according to Ning & Popov, Journal of the Electrochemical Society 151 (10), 2004 #TODO -> In slidepack case1 paper and case2 paper are swapped.
< Arrhenius relation for the rate parameter at the cell temperature
< Arrhenius relation for the diffusion constant at the cell temperature
< derivation of the formula: start equations (use the symbols from Yang, Leng, Zhang, Ge, Wang, Journal of Power Sources 360, 2017 but with opposite sign for j j = nFk c exp(..) j/nF = - D/delta (c - c0) j/nF = D/delta (- j/(nFk exp(..)) + c0) j * ( 1/(nFk exp(..)) + delta/DnF ) = c0 j = c0 / ( 1/(nFk exp(..)) + delta/DnF )
< Add the effects of this model
< model from Christensen & Newmann, Journal of the Electrochemical Society 152 (4), 2005
< Arrhenius relation for the rate parameter at the cell temperature
< Arrhenius relation for the diffusion constant at the cell temperature
< Use equation [22] from the paper
< the parameter a_L_K is set to 0.134461 but this constant can be lumped into the rate- and diffusion constants
< Add the effects of this model
< model from the optimisation in the paper
< Arrhenius relation for the rate parameter at the cell temperature
< Arrhenius relation for the diffusion constant at the cell temperature
< Use equation [22] from the paper
< the parameter a_L_K is set to 0.134461 but this constant can be lumped into the rate- and diffusion constants
< note: model used in optimisation had kpt/knt or vice versa, here a fixed value
< Add the effects of this model
< unknown degradation model
< Make the output for the SEI side reaction current density
< Calculate how much we decrease the volume fraction due to SEI growth
< don't decrease volume fraction
< decrease volume fraction according to Ashwin, Chung, Wang, Journal of Power Sources 328, 2016
< molar flux on the negative particle
< note: they use J = volumetric current [A/m3] -> they multiply with 'an' but we already have density [A/m2]
< unknown degradation model
void slide::Cell_SPM::setC | ( | double | cp0, |
double | cn0 | ||
) |
void setCurrent(bool critical, bool check, double I); //!< set the cell's current to the specified value -> From old slide.
set the concentrations to the given (uniform) concentration
< #NOTHOTFUNCTION
< The second transformation is to the eigenspace: z = V * u with V the inverse of the matrix with the eigenvectors. As explained, we know that there is one eigenvector corresponding to a uniform concentration So we need to calculate only this one nonzero value for the (twice) transformed concentration The location of the uniform eigenvector (which has a 0 eigenvalue) is written in M->Input[3]
< (twice) transformed negative uniform concentration
< (twice) transformed positive uniform concentration
< Do the first transformation, to u(i) = radius(i) * concentration = x(i) * R * concentration(i)
< ------------------------------------------------------------------------— loop to calculate the row of V * u corresponding to the uniform concentration
< Make the arrays for the (twice) transformed concentration with all zero
< except the ind which will be set soon.
< Set the cell current to 0 to reflect the boundary condition for a fully uniform concentration
< #TODO we do not do this.
< the stress values stored in the class variables for stress are no longer valid because the state has changed
void slide::Cell_SPM::setCharacterisationParam | ( | double | Dp, |
double | Dn, | ||
double | kpi, | ||
double | kni, | ||
double | Rdc | ||
) |
sets the parameters related to the characterisation of the cell
void Cell_SPM::setRamping(double Istep, double tstep) { /*
< store the rate constants in the cell
< calculate the specific electrode resistance from the total DC resistance
< overwrite the cycling related parameters
|
overridevirtual |
< #TODO in future make it into states.
Implements slide::StorageUnit.
void slide::Cell_SPM::setGeometricParameters | ( | double | capnom, |
double | elec_surf, | ||
double | ep, | ||
double | en, | ||
double | thickp, | ||
double | thickn | ||
) |
void setRamping(double Istep, double tstep); //!< sets the ramping parameters
sets the geometric parameters related to the amount of active material
< If you change the volume fraction, you also have to change the effective surface area
< overwrite the original geometric parameters, use cell version
void slide::Cell_SPM::setInitialConcentration | ( | double | cmaxp, |
double | cmaxn, | ||
double | lifracp, | ||
double | lifracn | ||
) |
sets the initial concentration
< Store the maximum concentrations
< set the concentration
< std::cout << "Throw test: " << 33 << '
';
void slide::Cell_SPM::setOCVcurve | ( | const std::string & | namepos, |
const std::string & | nameneg | ||
) |
sets the OCV curve of the cell to the given value
< verbose = verbosei;
< Cell_Fit::Cell_Fit(const slide::Model_SPM &MM, int verbosei) #TODO important for ageing to set no ageing. : Cell(settings::path::Kokam::namepos, settings::path::Kokam::nameneg, settings::path::Kokam::nameentropicC, settings::path::Kokam::nameentropicCell) { /*
< Store the OCV curves
< the OCV curve of the anode, the first column gives the lithium fractions (increasing), the 2nd column gives the OCV vs li/li+
< std::cout << "Throw test: " << 32 << '
';
|
overridevirtual |
< get the voltage Does not throw anymore!
< Restore states here.
Also not used except test functions. |
Implements slide::Cell.
|
inline |
|
overridevirtual |
opposite of getStates, check the states are valid?
< Back-up values.
< Copy states.
< Remove first Nstates elements from span.
< Restore states here.
Implements slide::StorageUnit.
|
overridevirtual |
set the cell's temperature
< #TODO if we need to check if we are in limits. if T is in limits.
< the stress values stored in the class variables for stress are no longer valid because the state has changed
Implements slide::StorageUnit.
void slide::Cell_SPM::setTenv | ( | double | Tenv | ) |
void setStates(const State_SPM &si, double I); //!< set the cell's states to the states in the State object and the cell current to the given value
set the environmental temperature
void Cell_SPM::setVlimits(double VMAX, double VMIN) { /*
< check if the environmental temperature is in the allowed limits
< the temperature limits are defined in State.hpp
|
inlineoverridevirtual |
|
inlineoverridevirtualnoexcept |
returns the uniform battery temperature in [K]
Implements slide::StorageUnit.
double slide::Cell_SPM::thermal_getTotalHeat | ( | ) |
function for unit testing
double getR() noexcept //!< get the total cell DC resistance { /*
|
overridevirtual |
Calculate the thermal model for this cell and update its cell temperature. The heat exchange in [W] with element i is given by Qc = Kneighbours[i]*A*(Tneighbours[i] - getT()); We assume all temperatures were constant for the past period, such that the thermal energy in [J] is given by Ec = Qc * Therm_time This is added up with the internal heath generated, and the cell's temperature is returned
IN Nneigtbours the number of neighbouring cells / cooling systems etc. Tneighbours array with the temperature of the neighbouring elements [K] Kneighbours array with the heat transfer constants of the neighbouring elements, k or h Aneighb array with the surface area of the neigbouring cell the heat transfer happens over the smaller one of Aneighb[i] and the area of this cell tim the time since the last thermal model calculation [for verification]
throws 8 invalid time keeping 9 invalid module temperature
< #TODO thermal model implementation should be outside.
< indicate we have ran the thermal model
Reimplemented from slide::Cell.
|
protected |
< total heat generation since last time the temperature was updated
< cooling with the environment
< cooling with the environment [W m-3]
< Calculate the new temperature rho * cp * dT/dt = Qtot / V where Qtot is total power in W V is the cell's volume L * elec_surf so integrated over time this is rho * cp * (Tnew - Told) = Etot / V
< Check the new temperature is valid, and if so, set it
< Reset the cumulative thermal variables
|
protected |
cell to cell variations // #TODO why do we store variations?
total heat generation [J] since start of cell's life
< check the time since the last checkup is not too large. if the parent has not called the thermal model for a while, the equation becomes unstable cause E = time * kA dT, so even a small dT will cause a huge E, and therefore a very large temperature swint
< #TODO magic number.
< Check the new temperature is valid, and if so, set it
< Reset the cumulative thermal variables
|
overridevirtual |
setStates(std::move(states)); //!< store new states, checks if they are illegal (throws an error in that case)
void Cell_SPM::ETI(bool print, double dti, bool blockDegradation) { /*
< check the time step is positive
< Update the stress values stored in the attributes with the stress of the previous time step
< Dai's stress in the positive particle in the previous time step
< Dai's stress in the negative particle in the previous time step
< Laresgoiti's stress in the negative particle in the previous time step
< ********************************************* Resolve the diffusion model for every dt time step ****************************************************************************************
< Calculate the time derivatives
< forward Euler time integration: s(t+1) = s(t) + ds/dt * dt
< increase the cumulative variables of this cell
< Calculate the internal heat generation of the cell
< If this cell has a parent module, this parent will call the thermal model with the correct parameters which will include heat exchange with the cell's neighbours and cooling from the cooling system of the module.
< if there is no parent, this cell is a stand-alone cell. We assume convective cooling with an environment If there is no parent, assume we update T every nstep*dt. So update the temperature now
< else it is the responsibility of the parent to call the thermal model function with the correct settings
< so cooling Qc = Qch * SAV * dT / (rho*cp) = Qch * A * dT / (rho*cp)
< calculate the surface of this cell
< #TODO in slide-pack it does not check the temperature.
< Calculate the stress values stored in the attributes for the stress in this time step
< only a few degradation models actually need the stress according to Dai, so only calculate it if needed
< only a few degradation models need the stress according to Laresgoiti
< Calculate the time derivatives
Implements slide::StorageUnit.
|
noexcept |
updated the stored stress values for Dai's stress model
< Make variables to store the stress
< Get the stress // #TODO these variables are not used.
< indicate that the values in the class variables are updated
void slide::Cell_SPM::updateLaresgoitiStress | ( | bool | ) |
update the stored stress values for Laresgoiti's stress model
setters void setVlimits(double VMAX, double VMIN); //!< set the voltage limits of the cell
< Update the stored value
< indicate that the values in the class variables are updated
|
overridevirtual |
print is an optional argument
< If the stored value is the most up to date one, then simply return this value
< print if the (global) verbose-setting is above the threshold
< Get the surface concentrations
< check if the surface concentration is within the allowed range 0 < cp < Cmaxpos 0 < cn < Cmaxneg don't allow 0 or Cmax because in that case, i0 will be 0, and the overpotentials will have 1/0 = inf or nan
< print error message unless you want to suppress the output
< Surface concentration is out of bounds.
< Calculate the li-fraction (instead of the li-concentration)
< Calculate the electrode potentials
< in linear interpolation, throw an error if you are out of the allowed range
< entropic coefficient of the total cell voltage [V/K]
< anode potential [V]
< cathode potential [V]
< current density on the electrodes [I m-2]
< Calculate the cell voltage the cell OCV at the reference temperature is OCV_p - OCV_n this OCV is adapted to the actual cell temperature using the entropic coefficient dOCV * (T - Tref) then the overpotentials and the resistive voltage drop are added
< we now have the most up to date value stored
Implements slide::StorageUnit.
|
overridevirtual |
checks if a state array is valid
< print if the (global) verbose-setting is above the threshold
< Check if all fields are present & extract their values
< check whether SoC in the allowed range
< there is no range on the current
Implements slide::StorageUnit.
|
inlineoverridevirtual |
Only for cells to see individual states.
Reimplemented from slide::StorageUnit.
|
protected |
charge transfer constant for the plating reaction [-]
|
protected |
charge transfer coefficient of the SEI reaction [-]
|
protected |
Li- concentration in electrolyte [mol m-3] standard concentration of 1 molar.
|
protected |
bulk concentration of the electrolyte molecule participating in the SEI growth (e.g. EC) [mol m-3]
|
protected |
maximum lithium concentration in the anode [mol m-3] value for C
|
protected |
maximum lithium concentration in the cathode [mol m-3] value for NMC
|
protected |
thermal capacity of the battery #TODO = units missing.
Geometric parameters
|
protected |
structure with the fitting parameters of the different crack growth models
LAM parameters & constants
DEG_ID slide::Cell_SPM::deg_id |
structure with the identification of which degradation model(s) to use #TODO may be protected.
|
protected |
activation energy for the Arrhenius relation of Dn
Thermal model parameters
|
protected |
activation energy for the Arrhenius relation of Dp
|
protected |
other geometric parameters are part of State because they can change over the battery's lifetime
|
protected |
rate constant of main reaction at negative electrode at reference temperature
|
protected |
The diffusion constants at reference temperature are part of State because they can change over the battery's lifetime.
activation energy for the Arrhenius relation of kn
|
protected |
rate constant of main reaction at positive electrode at reference temperature
|
protected |
activation energy for the Arrhenius relation of kp
|
protected |
structure with the fitting parameters of the different LAM models
Li-plating parameters & constants //!< #TODO if we can make these static and speed gain.
|
protected |
OCV curves.
|
protected |
number of electrons involved in the main reaction [-] #TODO if really constant?
parameters of the main li-insertion reaction
|
protected |
number of electrons involved in the plating reaction [-]
|
protected |
number of electrons involved in the SEI reaction [-]
|
protected |
|
protected |
equilibrium potential of the NMC dissolution side reaction [V]
|
protected |
OCV of the plating reaction [V].
|
protected |
equilibrium potential of the SEI side reaction [V]
|
protected |
structure with the fitting parameters of the different plating models
Matrices for spatial discretisation of the solid diffusion model
|
protected |
convective heat transfer coefficient per volume [W K-1 m-3]
|
protected |
density of the battery
|
protected |
density of the plated lithium layer
|
protected |
partial molar volume of the SEI layer [m3 mol-1]
|
protected |
specific resistance times real surface area of the SEI film [Ohm m] ? #TODO if unit is correct. aging fit.
|
protected |
the battery current/initial state, grouping all parameter which change over the battery's lifetime (see State_SPM.hpp)
Battery model constants
|
protected |
structure with the fitting parameters of the different SEI growth models
surface crack parameters & constants
|
protected |
Stress parameters.
Constants and parameters for the SEI growth model //!< #TODO if we can make these static and speed gain.
|
protected |
< protected such that child classes can access the class variables
|
protected |
environment temperature [K]
|
protected |
reference temperature [K]
Qch: 90 gives very good cooling, as if there is a fan pointed at the cell. values of 30-60 are representative for a cell on a shelf without forced cooling 40 is representative for a cell on a shelf without forced cooling
|
protected |
total heat generation since the last update [J]
|
protected |
variable for unit testing, total heat generation since the beginning of this cell's life [J]
|
protected |
time since the last update of the thermal model
|
protected |
relative factor increasing the capacity of the cell
|
protected |
relative factor to speed up or slow down the rate of LAM
|
protected |
relative factor to speed up or slow down the rate of SEI growth
|
protected |
relative factor increasing the DC resistance
|
protected |
Functions.
|
protected |
partial molar volume of the main reaction, see Ashwin et al, 2016
|
protected |
partial molar volume of the SEI side reaction, see Ashwin et al., 2016