Name | Description |
---|---|
PhaseChange | Examples of phase change |
Reactions | Examples of phase change |
AirColumn | Gas in a vertical array of subregions, affected by gravity |
BinaryDiffusion | Pressure-driven flow of H2, dragging H2O |
Echo | Two regions of gas with an initial pressure difference |
EchoCentral | Two regions of gas with an initial pressure difference, with central difference scheme |
ElectricalConduction | Migration of e- through C+ |
InternalFlow | Internal, laminar flow of liquid water |
Subregion | Single subregion, with H2 by default |
Subregions | Horizontal array of subregions with an initial pressure difference (H2 included by default) |
ThermalConduction | Thermal conduction (through solid) |
ThermalConductionConvection | Thermal conduction through a solid alongside conduction and convection through a gas |
This is a model of a vertical column of 10 × 10 × 10 m regions with N2 gas. The upper boundary is held at 1 bar and the lower boundary has zero velocity. The initial pressure difference is zero, but a gas enters the upper boundary and a pressure difference develops due to gravity. There are oscillations due to the inertance and compression of the gas. After about 1.5 s, the pressure difference settles to Ly ay m ρ as expected.
A temperature gradient is created due to the thermodynamics of the expanding and contracting gases. It takes much longer (over a year!) for the temperatures to equalize due to the size of the system and the low thermal conductivity of the gas. With a stiff solver, the model should simulate at this time scale as well.
The damping factor (k) can be used to scale the continuity (ζ) of the gas in the regions. The oscillations are dampened considerably at k = 100. However, with high values of the factor, the boundary pressures are decoupled from the pressures in the region because the nonequilibrium force is considerable.
Assumptions:
Type | Name | Default | Description |
---|---|---|---|
Integer | n_y | 3 | Number of discrete subregions along the y axis, besides the 2 boundary subregions |
Pressure | Deltap_IC | 0 | Initial pressure difference [M/(L.T2)] |
NumberAbsolute | k | 5 | Damping factor [1] |
model AirColumn "Gas in a vertical array of subregions, affected by gravity" import FCSys.Utilities.round; extends Modelica.Icons.Example; parameter Integer n_y=3 "Number of discrete subregions along the y axis, besides the 2 boundary subregions"; parameter Q.Pressure Deltap_IC=0 "Initial pressure difference"; output Q.Pressure Deltap=subregion2.gas.N2.p - subregion1.gas.N2.p "Measured pressure difference"; output Q.Pressure Deltap_ex=-(n_y + 1)*10*U.m*environment.a[Axis.y]* Characteristics.N2.Gas.m*subregions[round(n_y/2)].gas.N2.rho "Expected pressure difference"; parameter Q.NumberAbsolute k=5 "Damping factor"; inner Conditions.Environment environment(p=U.bar); FCSys.Subregions.Subregion subregion1( L={10,10,10}*U.m, inclTransX=false, inclTransY=true, inclTransZ=false, gas(inclN2=true, N2( zeta=k*FCSys.Characteristics.H2.Gas.zeta(), p_IC=environment.p - Deltap_IC/2, upstreamY=false))); FCSys.Subregions.Subregion subregions[n_y]( each L={10,10,10}*U.m, each inclTransZ=false, gas(each inclN2=true, N2( each zeta=k*FCSys.Characteristics.H2.Gas.zeta(), p_IC={environment.p - Deltap_IC/2 - i*Deltap_IC/(n_y + 1) for i in 1: n_y}, each upstreamY=false)), each inclTransX=false, each inclTransY=true) if n_y > 0; FCSys.Subregions.Subregion subregion2( L={10,10,10}*U.m, inclTransZ=false, inclTransX=false, inclTransY=true, gas(inclN2=true, N2( zeta=k*FCSys.Characteristics.H2.Gas.zeta(), p_IC=environment.p + Deltap_IC/2, upstreamY=false))); Conditions.ByConnector.BoundaryBus.Single.Sink BC(gas(inclN2=true, N2( materialSet(y=environment.p)))); equation connect(subregion1.yPositive, subregions[1].yNegative); for i in 1:n_y - 1 loop connect(subregions[1:n_y - 1].yPositive, subregions[2:n_y].yNegative) "Not shown in the diagram"; end for; if n_y == 0 then connect(subregion1.yPositive, subregion2.yNegative) "Not shown in the diagram"; end if; connect(subregions[n_y].yPositive, subregion2.yNegative); connect(BC.boundary, subregion2.yPositive); end AirColumn;
Type | Name | Default | Description |
---|---|---|---|
Species | |||
Boolean | 'inclC+' | true | Carbon plus (C+) |
Boolean | 'inclSO3-' | false | Nafion sulfonate (C19HF37O5S-, abbreviated as SO3-) |
Boolean | 'incle-' | false | Electrons (e-) |
Boolean | 'inclH+' | false | Protons (H+) |
Boolean | inclH2 | true | Hydrogen (H2) |
Boolean | inclH2O | true | Water vapor (H2O) |
Boolean | inclN2 | false | Nitrogen (N2) |
Boolean | inclO2 | false | Oxygen (O2) |
model BinaryDiffusion "Pressure-driven flow of H2, dragging H2O" extends Subregion( inclH2O=true, 'inclC+'=true, subregion(common(k_Phi={1e-5,1e-5,1e-5}),gas( common(k_Phi={1e4,1e4,1e4}), H2(T(stateSelect=StateSelect.always), phi(each stateSelect=StateSelect.always)), H2O(initEnergy=Init.none, boundaries(Ndot(each stateSelect=StateSelect.always))))), environment(RH=0.6)); Conditions.ByConnector.BoundaryBus.Single.Source source(gas( inclH2=true, inclH2O=true, H2( thermalSet(y=environment.T), redeclare function materialSpec = Conditions.ByConnector.Boundary.Single.Material.pressure, redeclare Modelica.Blocks.Sources.Ramp materialSet( height=-U.Pa, offset=environment.p_dry, duration=10)), H2O( redeclare function materialSpec = Conditions.ByConnector.Boundary.Single.Material.pressure, materialSet(y=environment.p_H2O), thermalSet(y=environment.T)))); Conditions.ByConnector.BoundaryBus.Single.Source sink(gas( inclH2=true, inclH2O=true, H2(materialSet(y=-source.gas.H2.boundary.Ndot), thermalSet(y=environment.T)), H2O( redeclare function materialSpec = Conditions.ByConnector.Boundary.Single.Material.pressure, materialSet(y=environment.p_H2O), thermalSet(y=environment.T)))); equation connect(source.boundary, subregion.xNegative); connect(subregion.xPositive, sink.boundary); end BinaryDiffusion;
Type | Name | Default | Description |
---|---|---|---|
NumberAbsolute | k | 1 | Damping factor [1] |
Integer | n_x | 0 | Number of discrete subregions along the x axis, besides the 2 side subregions |
Pressure | Deltap_IC | -100*U.Pa | Initial pressure difference [M/(L.T2)] |
Species | |||
Boolean | 'inclC+' | false | Carbon plus (C+) |
Boolean | 'inclSO3-' | false | Nafion sulfonate (C19HF37O5S) |
Boolean | 'incle-' | false | Electrons (e-) |
Boolean | 'inclH+' | false | Protons (H+) |
Boolean | inclH2 | true | Hydrogen (H2) |
Boolean | inclH2O | false | Water vapor (H2O) |
Boolean | inclN2 | false | Nitrogen (N2) |
Boolean | inclO2 | false | Oxygen (O2) |
model Echo "Two regions of gas with an initial pressure difference" parameter Q.NumberAbsolute k=1 "Damping factor"; extends Subregions(subregion1(gas(H2(zeta=k*Characteristics.H2.Gas.zeta()))), subregion2(gas(H2(zeta=k*FCSys.Characteristics.H2.Gas.zeta())))); end Echo;
Type | Name | Default | Description |
---|---|---|---|
NumberAbsolute | k | 1 | Damping factor [1] |
Integer | n_x | 0 | Number of discrete subregions along the x axis, besides the 2 side subregions |
Pressure | Deltap_IC | -100*U.Pa | Initial pressure difference [M/(L.T2)] |
Species | |||
Boolean | 'inclC+' | false | Carbon plus (C+) |
Boolean | 'inclSO3-' | false | Nafion sulfonate (C19HF37O5S) |
Boolean | 'incle-' | false | Electrons (e-) |
Boolean | 'inclH+' | false | Protons (H+) |
Boolean | inclH2 | true | Hydrogen (H2) |
Boolean | inclH2O | false | Water vapor (H2O) |
Boolean | inclN2 | false | Nitrogen (N2) |
Boolean | inclO2 | false | Oxygen (O2) |
model EchoCentral "Two regions of gas with an initial pressure difference, with central difference scheme" extends Echo( subregion1(gas( H2(upstreamX=false), H2O(upstreamX=false), N2(upstreamX=false), O2(upstreamX=false))), subregions(gas( H2(each upstreamX=false), H2O(each upstreamX=false), N2(each upstreamX=false), O2(each upstreamX=false))), subregion2(gas( H2(upstreamX=false), H2O(upstreamX=false), N2(upstreamX=false), O2(upstreamX=false)))); end EchoCentral;
This is an example of Ohm's law. The 1 cm × 1 mm × 1 mm subregion contains carbon (C+) and electrons. An electrical current of 50 mA (5 A/cm2) is delivered into the negative boundary and exits the positive boundary. Due to the finite mobility of the electrons a force is required to support the current; this maps directly to electrical potential. The example shows that the measured resistance is R = L/(A ρ μ) as expected, where ρ is electronic density and μ is electronic mobility.
The measured rate of heat generation (subregion.graphite.'e-'.Edot_DT
)
is equal to P = (zI)2 R as expected, where
zI = 50 mA is the electrical current. This heat is conducted through the carbon
to the boundaries, which are held at 25 °C. The measured steady state temperature
is T = T0 + θ L P/(4 A) as expected, where
T0 = 25 °C is the boundary temperature and
θ is the thermal resistance. The factor of one fourth
is due to the boundary conditions; the conduction length is half of the total length
and the heat is rejected to both sides. There is no thermal convection or
radiation—only conduction to the sides.
Type | Name | Default | Description |
---|---|---|---|
Species | |||
Boolean | 'inclC+' | true | Carbon plus (C+) |
Boolean | 'inclSO3-' | false | Nafion sulfonate (C19HF37O5S-, abbreviated as SO3-) |
Boolean | 'incle-' | true | Electrons (e-) |
Boolean | 'inclH+' | false | Protons (H+) |
Boolean | inclH2 | false | Hydrogen (H2) |
Boolean | inclH2O | false | Water vapor (H2O) |
Boolean | inclN2 | false | Nitrogen (N2) |
Boolean | inclO2 | false | Oxygen (O2) |
model ElectricalConduction "Migration of e- through C+" output Q.Potential w=subregion.graphite.'e-'.Deltag[1] "Electrical potential"; output Q.Current zI=-subregion.graphite.'e-'.I[1] "Electrical current"; output Q.ResistanceElectrical R=w/zI "Measured electrical resistance"; output Q.ResistanceElectrical R_ex=subregion.L[Axis.x]/(subregion.graphite. 'e-'.sigma*subregion.A[Axis.x]) "Expected electrical resistance"; output Q.Power P=subregion.graphite.'e-'.Edot_AT "Measured rate of heat generation"; output Q.Power P_ex=zI^2*R "Expected rate of heat generation"; output Q.TemperatureAbsolute T=subregion.graphite.'C+'.T "Measured temperature"; output Q.TemperatureAbsolute T_ex=environment.T + subregion.graphite.'C+'.theta *U.cm*P/(4*subregion.A[Axis.x]) "Expected temperature"; extends Examples.Subregion( 'inclC+'=true, 'incle-'=true, inclH2=false, subregion(L={U.cm,U.mm,U.mm}, graphite('C+'(epsilon=1),'e-'(sigma=1e2*U.S/U.m)))); Conditions.ByConnector.BoundaryBus.Single.Source BC1(graphite( 'incle-'='incle-', 'inclC+'=true, 'C+'(set(y=environment.T)), 'e-'( materialSet(y=-0.05*U.A), redeclare function thermalSpec = Conditions.ByConnector.Boundary.Single.Thermal.heatRate, thermalSet(y=0)))); Conditions.ByConnector.BoundaryBus.Single.Sink BC2(graphite( 'inclC+'=true, 'incle-'='incle-', redeclare Conditions.ByConnector.ThermalDiffusive.Single.Temperature 'C+' (set(y=environment.T)))); equation connect(BC1.boundary, subregion.xNegative); connect(subregion.xPositive, BC2.boundary); end ElectricalConduction;
A small-signal variation is added to the time-average flow rate in order to demonstrate the effects of inertance.
Note that the temperature increases due to viscous dissipation, but the increase is limited by thermal convection. The walls are adiabatic.
Extends from Examples.Subregion (Single subregion, with H2 by default).
Type | Name | Default | Description |
---|---|---|---|
VolumeRate | Vdot_large | -1.5*U.cc/U.s | Prescribed large signal volumetric flow rate [L3/T] |
Species | |||
Boolean | 'inclC+' | false | Carbon plus (C+) |
Boolean | 'inclSO3-' | false | Nafion sulfonate (C19HF37O5S-, abbreviated as SO3-) |
Boolean | 'incle-' | false | Electrons (e-) |
Boolean | 'inclH+' | false | Protons (H+) |
Boolean | inclH2 | false | Hydrogen (H2) |
Boolean | inclH2O | false | Water vapor (H2O) |
Boolean | inclN2 | false | Nitrogen (N2) |
Boolean | inclO2 | false | Oxygen (O2) |
model InternalFlow "Internal, laminar flow of liquid water" import FCSys.Utilities.Delta; final parameter Q.Area A=subregion.A[Axis.x] "Cross-sectional area"; // Conditions parameter Q.VolumeRate Vdot_large=-1.5*U.cc/U.s "Prescribed large signal volumetric flow rate"; // Measurements output Q.Pressure Deltap=Delta(subregion.liquid.H2O.boundaries[1, :].p) "Measured pressure difference"; output Q.Length D=2*A/(subregion.L[Axis.y] + subregion.L[Axis.z]); output Q.Number Re=subregion.liquid.H2O.phi[Axis.x]*D*subregion.liquid.H2O.eta *subregion.liquid.H2O.Data.m*subregion.liquid.H2O.rho if environment.analysis "Reynolds number"; output Q.Pressure Deltap_Poiseuille=-32*subregion.L[Axis.x]*subregion.liquid.H2O.phi[ Axis.x]/(D^2*subregion.liquid.H2O.eta) "Pressure difference according to Poiseuille's law"; Real x=Deltap_Poiseuille/Deltap; output Q.Power Qdot_gen_Poiseuille=-Deltap_Poiseuille*Vdot "Rate of heat generation according to Poiseuille's law"; output Q.VolumeRate Vdot=BC1.liquid.H2O.materialSet.y "Total volumetric flow rate"; extends Examples.Subregion(inclH2=false, subregion( L={U.m,U.mm,U.mm}, inclTransY=true, inclTransZ=true, liquid(inclH2O=true, H2O(initMaterial=Init.none)))); Conditions.ByConnector.BoundaryBus.Single.Source BC1(liquid(inclH2O=true, H2O( redeclare Modelica.Blocks.Sources.Sine materialSet( amplitude=0.2*Vdot_large, offset=Vdot_large, freqHz=1), redeclare function materialSpec = Conditions.ByConnector.Boundary.Single.Material.volumeRate (redeclare package Data = FCSys.Characteristics.H2O.Liquid), thermalSet(y=environment.T)))); Conditions.ByConnector.BoundaryBus.Single.Sink BC2(liquid(inclH2O=true, H2O( materialSet(y=U.atm)))); Conditions.ByConnector.BoundaryBus.Single.Source BC3(liquid(inclH2O=true, H2O( redeclare function thermalSpec = Conditions.ByConnector.Boundary.Single.Thermal.heatRate, thermalSet(y=0)))); Conditions.ByConnector.BoundaryBus.Single.Source BC4(liquid(inclH2O=true, H2O( redeclare function thermalSpec = Conditions.ByConnector.Boundary.Single.Thermal.heatRate, thermalSet(y=0)))); Conditions.ByConnector.BoundaryBus.Single.Source BC5(liquid(inclH2O=true, H2O( redeclare function thermalSpec = Conditions.ByConnector.Boundary.Single.Thermal.heatRate, thermalSet(y=0)))); Conditions.ByConnector.BoundaryBus.Single.Source BC6(liquid(inclH2O=true, H2O( redeclare function thermalSpec = Conditions.ByConnector.Boundary.Single.Thermal.heatRate, thermalSet(y=0)))); equation connect(BC1.boundary, subregion.xNegative); connect(BC2.boundary, subregion.xPositive); connect(subregion.yNegative, BC3.boundary); connect(BC4.boundary, subregion.yPositive); connect(subregion.zPositive, BC6.boundary); connect(subregion.zNegative, BC5.boundary); end InternalFlow;
This model is boring. It just establishes a single subregion with H2 by default. There are no boundary conditions other than those implied by the open connectors (no diffusion current, no forces, no thermal conduction). Other examples in this package extend from this one.
Extends from Modelica.Icons.Example (Icon for runnable examples).
Type | Name | Default | Description |
---|---|---|---|
Species | |||
Boolean | 'inclC+' | false | Carbon plus (C+) |
Boolean | 'inclSO3-' | false | Nafion sulfonate (C19HF37O5S-, abbreviated as SO3-) |
Boolean | 'incle-' | false | Electrons (e-) |
Boolean | 'inclH+' | false | Protons (H+) |
Boolean | inclH2 | true | Hydrogen (H2) |
Boolean | inclH2O | false | Water vapor (H2O) |
Boolean | inclN2 | false | Nitrogen (N2) |
Boolean | inclO2 | false | Oxygen (O2) |
model Subregion "Single subregion, with H2 by default" extends Modelica.Icons.Example; parameter Boolean 'inclC+'=false "Carbon plus (C+)"; parameter Boolean 'inclSO3-'=false "Nafion sulfonate (C19HF37O5S-, abbreviated as SO3-)"; parameter Boolean 'incle-'=false "Electrons (e-)"; parameter Boolean 'inclH+'=false "Protons (H+)"; parameter Boolean inclH2=true "Hydrogen (H2)"; parameter Boolean inclH2O=false "Water vapor (H2O)"; parameter Boolean inclN2=false "Nitrogen (N2)"; parameter Boolean inclO2=false "Oxygen (O2)"; inner Conditions.Environment environment(RH=0); FCSys.Subregions.Subregion subregion( L={1,1,1}*U.cm, inclTransY=false, inclTransZ=false, graphite(final 'inclC+'='inclC+', final 'incle-'='incle-'), ionomer(final 'inclSO3-'='inclSO3-', final 'inclH+'='inclH+'), gas( final inclH2=inclH2, final inclH2O=inclH2O, final inclN2=inclN2, final inclO2=inclO2)); end Subregion;
Type | Name | Default | Description |
---|---|---|---|
Integer | n_x | 0 | Number of discrete subregions along the x axis, besides the 2 side subregions |
Pressure | Deltap_IC | -100*U.Pa | Initial pressure difference [M/(L.T2)] |
Species | |||
Boolean | 'inclC+' | false | Carbon plus (C+) |
Boolean | 'inclSO3-' | false | Nafion sulfonate (C19HF37O5S) |
Boolean | 'incle-' | false | Electrons (e-) |
Boolean | 'inclH+' | false | Protons (H+) |
Boolean | inclH2 | true | Hydrogen (H2) |
Boolean | inclH2O | false | Water vapor (H2O) |
Boolean | inclN2 | false | Nitrogen (N2) |
Boolean | inclO2 | false | Oxygen (O2) |
model Subregions "Horizontal array of subregions with an initial pressure difference (H2 included by default)" extends Modelica.Icons.Example; parameter Integer n_x=0 "Number of discrete subregions along the x axis, besides the 2 side subregions"; parameter Q.Pressure Deltap_IC=-100*U.Pa "Initial pressure difference"; parameter Boolean 'inclC+'=false "Carbon plus (C+)"; parameter Boolean 'inclSO3-'=false "Nafion sulfonate (C19HF37O5S)"; parameter Boolean 'incle-'=false "Electrons (e-)"; parameter Boolean 'inclH+'=false "Protons (H+)"; parameter Boolean inclH2=true "Hydrogen (H2)"; parameter Boolean inclH2O=false "Water vapor (H2O)"; parameter Boolean inclN2=false "Nitrogen (N2)"; parameter Boolean inclO2=false "Oxygen (O2)"; inner Conditions.Environment environment; FCSys.Subregions.Subregion subregion1( L={1,1,1}*U.cm, inclTransY=false, inclTransZ=false, gas( final inclH2=inclH2, final inclH2O=inclH2O, final inclN2=inclN2, final inclO2=inclO2, H2(p_IC=environment.p - Deltap_IC/2, phi(each stateSelect=StateSelect.always, each fixed=true)), H2O(p_IC=environment.p - Deltap_IC/2, phi(each stateSelect=StateSelect.always, each fixed=true)), N2(p_IC=environment.p - Deltap_IC/2, phi(each stateSelect=StateSelect.always, each fixed=true)), O2(p_IC=environment.p - Deltap_IC/2, phi(each stateSelect=StateSelect.always, each fixed=true))), graphite(final 'inclC+'='inclC+', final 'incle-'='incle-'), ionomer(final 'inclSO3-'='inclSO3-', final 'inclH+'='inclH+'), liquid(H2O(epsilon_IC=0.25))); FCSys.Subregions.Subregion subregions[n_x]( each L={1,1,1}*U.cm, each inclTransY=false, each inclTransZ=false, gas( each final inclH2=inclH2, each final inclH2O=inclH2O, each final inclN2=inclN2, each final inclO2=inclO2, H2(p_IC={environment.p - Deltap_IC/2 - i*Deltap_IC/(n_x + 1) for i in 1: n_x}, each phi(each stateSelect=StateSelect.always, each fixed=true)), H2O(p_IC={environment.p - Deltap_IC/2 - i*Deltap_IC/(n_x + 1) for i in 1: n_x}, each phi(each stateSelect=StateSelect.always, each fixed=true)), N2(p_IC={environment.p - Deltap_IC/2 - i*Deltap_IC/(n_x + 1) for i in 1: n_x}, each phi(each stateSelect=StateSelect.always, each fixed=true)), O2(p_IC={environment.p - Deltap_IC/2 - i*Deltap_IC/(n_x + 1) for i in 1: n_x}, each phi(each stateSelect=StateSelect.always, each fixed=true))), graphite(each final 'inclC+'='inclC+', each final 'incle-'='incle-'), ionomer(each final 'inclSO3-'='inclSO3-', each final 'inclH+'='inclH+'), liquid(H2O(each epsilon_IC=0.25))) if n_x > 0; FCSys.Subregions.Subregion subregion2( L={1,1,1}*U.cm, inclTransY=false, inclTransZ=false, gas( final inclH2=inclH2, final inclH2O=inclH2O, final inclN2=inclN2, final inclO2=inclO2, H2(p_IC=environment.p + Deltap_IC/2), H2O(p_IC=environment.p + Deltap_IC/2), N2(p_IC=environment.p + Deltap_IC/2), O2(p_IC=environment.p + Deltap_IC/2)), graphite(final 'inclC+'='inclC+', final 'incle-'='incle-'), ionomer(final 'inclSO3-'='inclSO3-', final 'inclH+'='inclH+'), liquid(H2O(each epsilon_IC=0.25))); equation connect(subregion1.xPositive, subregions[1].xNegative); for i in 1:n_x - 1 loop connect(subregions[i].xPositive, subregions[i + 1].xNegative) "Not shown in the diagram"; end for; if n_x == 0 then connect(subregion1.xPositive, subregion2.xNegative) "Not shown in the diagram"; end if; connect(subregions[n_x].xPositive, subregion2.xNegative); end Subregions;
Type | Name | Default | Description |
---|---|---|---|
Integer | n_x | 6 | Number of discrete subregions along the x axis, besides the 2 side subregions |
Pressure | Deltap_IC | -100*U.Pa | Initial pressure difference [M/(L.T2)] |
Species | |||
Boolean | 'inclC+' | true | Carbon plus (C+) |
Boolean | 'inclSO3-' | false | Nafion sulfonate (C19HF37O5S) |
Boolean | 'incle-' | false | Electrons (e-) |
Boolean | 'inclH+' | false | Protons (H+) |
Boolean | inclH2 | false | Hydrogen (H2) |
Boolean | inclH2O | false | Water vapor (H2O) |
Boolean | inclN2 | false | Nitrogen (N2) |
Boolean | inclO2 | false | Oxygen (O2) |
model ThermalConduction "Thermal conduction (through solid)" extends Examples.Subregions( n_x=6, 'inclC+'=true, inclH2=false, subregion1(graphite('C+'(T_IC=environment.T + 30*U.K, epsilon=1))), subregions(each graphite('C+'(epsilon=1))), subregion2(graphite('C+'(epsilon=1)))); end ThermalConduction;
Type | Name | Default | Description |
---|---|---|---|
Integer | n_x | 6 | Number of discrete subregions along the x axis, besides the 2 side subregions |
Pressure | Deltap_IC | -100*U.Pa | Initial pressure difference [M/(L.T2)] |
Species | |||
Boolean | 'inclC+' | true | Carbon plus (C+) |
Boolean | 'inclSO3-' | false | Nafion sulfonate (C19HF37O5S) |
Boolean | 'incle-' | false | Electrons (e-) |
Boolean | 'inclH+' | false | Protons (H+) |
Boolean | inclH2 | false | Hydrogen (H2) |
Boolean | inclH2O | false | Water vapor (H2O) |
Boolean | inclN2 | true | Nitrogen (N2) |
Boolean | inclO2 | false | Oxygen (O2) |
model ThermalConductionConvection "Thermal conduction through a solid alongside conduction and convection through a gas" extends ThermalConduction( inclN2=true, subregion1(gas(N2( T_IC=subregion1.graphite.'C+'.T_IC, phi(stateSelect=StateSelect.always, displayUnit="mm/s"), phi_boundaries(each displayUnit="mm/s"))), graphite('C+'(epsilon=0.5))), subregions( common(each k_Phi={10,10,10}), gas(N2(each phi(each stateSelect=StateSelect.always, displayUnit="mm/s"), each phi_boundaries(each displayUnit="mm/s"))), each graphite('C+'(epsilon=0.5))), subregion2(gas(N2(phi(displayUnit="mm/s"), phi_boundaries(each displayUnit="mm/s"))), graphite('C+'(epsilon=0.5))), environment(psi_O2_dry=0, RH=0)); end ThermalConductionConvection;