Fluxes#
Class Flux
#
- class multipy.flux.Flux(species_velocities)#
Supports computing and storing fluxes. This class assumes that the species velocities, \(\mathbf{u}_i\), are known.
Diffusive fluxes:
mass diffusive flux relative to a mass-averaged velocity, \(\mathbf{j}_i\)
mass diffusive flux relative to a molar-averaged velocity, \(\mathbf{j}_i^u\)
molar diffusive flux relative to a mass-averaged velocity, \(\mathbf{J}_i^v\)
molar diffusive flux relative to a molar-averaged velocity, \(\mathbf{J}_i\)
- Parameters
species_velocities – vector
numpy.ndarray
specifying the species velocities \(\mathbf{u}_i\) in \([m/s]\). It should be of size(n_species,n_observations)
.
Getters:
get_species_velocities
get_diffusive_molar_molar (is set to
None
at class init)get_diffusive_molar_mass (is set to
None
at class init)get_diffusive_mass_molar (is set to
None
at class init)get_diffusive_mass_mass (is set to
None
at class init)
Setters:
set_species_velocities
set_diffusive_molar_molar (is set to
None
at class init)set_diffusive_molar_mass (is set to
None
at class init)set_diffusive_mass_molar (is set to
None
at class init)set_diffusive_mass_mass (is set to
None
at class init)
Flux.plot_diffusive_flux
#
- multipy.flux.Flux.plot_diffusive_flux(self, species_names=None, colors=None, figsize=(10, 5), filename=None)#
Plots the computed diffusive fluxes.
Example:
- Parameters
species_names – (optional)
list
ofstr
specifying the species names.colors – (optional)
list
ofstr
specifying the plotting colors for each species. Example:colors=['#C7254E', '#BBBBBB', '#008CBA']
.figsize – (optional)
tuple
specifying the figure size.filename – (optional)
str
specifying the filename. If set toNone
, plot will not be saved to a file.
Flux.diffusive_molar_molar
#
- multipy.flux.Flux.diffusive_molar_molar(self, species_mole_fractions, species_molar_densities)#
Computes the molar diffusive flux relative to a molar-averaged velocity:
\[\mathbf{J}_i = c_i \mathbf{u}_i + c_i \mathbf{u}\]- Parameters
species_mole_fractions – scalar
numpy.ndarray
specifying the species mole fractions, \(X_i\), in \([-]\). It should be of size(n_species,n_observations)
.species_molar_densities – scalar
numpy.ndarray
specifying the molar densities of species, \(c_i\), in \([mole/m^3]\). It should be of size(n_species,n_observations)
.
- Returns
diffusive_flux - vector
numpy.ndarray
of molar diffusive fluxes relative to a molar-averaged velocity \(\mathbf{J}_i\) in \([mole/(m^2s)]\). It has size(n_species,n_observations)
.
Flux.diffusive_molar_mass
#
- multipy.flux.Flux.diffusive_molar_mass(self, species_mass_fractions, species_molar_densities)#
Computes the molar diffusive flux relative to a mass-averaged velocity:
\[\mathbf{J}_i^v = c_i \mathbf{u}_i + c_i \mathbf{v}\]- Parameters
species_mass_fractions – scalar
numpy.ndarray
specifying the species mass fractions, \(Y_i\), in \([-]\). It should be of size(n_species,n_observations)
.species_molar_densities – scalar
numpy.ndarray
specifying the species molar densities \(c_i\) in \([mole/m^3]\). It should be of size(n_species,n_observations)
.
- Returns
diffusive_flux - vector
numpy.ndarray
of molar diffusive fluxes relative to a mass-averaged velocity \(\mathbf{J}_i^v\) in \([mole/(m^2s)]\). It has size(n_species,n_observations)
.
Flux.diffusive_mass_molar
#
- multipy.flux.Flux.diffusive_mass_molar(self, species_mole_fractions, species_mass_densities)#
Computes the mass diffusive flux relative to a molar-averaged velocity:
\[\mathbf{j}_i^u = \rho_i \mathbf{u}_i + \rho_i \mathbf{u}\]- Parameters
species_mole_fractions – scalar
numpy.ndarray
specifying the species mole fractions \(X_i\) in \([-]\). It should be of size(n_species,n_observations)
.species_mass_densities – scalar
numpy.ndarray
specifying the species mass densities \(\mathbf{\rho}_i\) in \([kg/m^3]\). It should be of size(n_species,n_observations)
.
- Returns
diffusive_flux - vector
numpy.ndarray
of mass diffusive fluxes relative to a molar-averaged velocity \(\mathbf{j}_i^u\) in \([kg/(m^2s)]\). It has size(n_species,n_observations)
.
Flux.diffusive_mass_mass
#
- multipy.flux.Flux.diffusive_mass_mass(self, species_mass_fractions, species_mass_densities)#
Computes the mass diffusive flux relative to a mass-averaged velocity:
\[\mathbf{j}_i = \rho_i \mathbf{u}_i + \rho_i \mathbf{v}\]- Parameters
species_mass_fractions – scalar
numpy.ndarray
specifying the species mass fractions \(Y_i\) in \([-]\). It should be of size(n_species, n_observations)
.species_mass_densities – scalar
numpy.ndarray
specifying the species mass densities \(\mathbf{\rho}_i\) in \([kg/m^3]\). It should be of size(n_species, n_observations)
.
- Returns
diffusive_flux - vector
numpy.ndarray
of mass diffusive fluxes relative to a mass-averaged velocity \(\mathbf{j}_i\) in \([kg/(m^2s)]\). It has size(n_species, n_observations)
.