Composition#

Class Composition#

class multipy.composition.Composition#

Supports computing and storing composition-related quantities.

Getters:

  • get_species_mole_fractions (is set to None at class init)

  • get_species_mass_fractions (is set to None at class init)

  • get_species_volume_fractions (is set to None at class init)

  • get_grad_species_mole_fractions (is set to None at class init)

  • get_grad_species_mass_fractions (is set to None at class init)

  • get_grad_species_volume_fractions (is set to None at class init)

  • get_species_molar_densities (is set to None at class init)

  • get_species_mass_densities (is set to None at class init)

  • get_mixture_molar_density (is set to None at class init)

  • get_mixture_molar_volume (is set to None at class init)

  • get_mixture_molar_mass (is set to None at class init)

  • get_mixture_mass_density (is set to None at class init)

Setters:

  • set_species_mole_fractions setter for get_species_mole_fractions

  • set_species_mass_fractions setter for get_species_mass_fractions

  • set_species_volume_fractions setter for get_species_volume_fractions

  • set_grad_species_mole_fractions setter for get_grad_species_mole_fractions

  • set_grad_species_mass_fractions setter for get_grad_species_mass_fractions

  • set_grad_species_volume_fractions setter for get_grad_species_volume_fractions


Composition.plot_species_mole_fractions#

multipy.composition.Composition.plot_species_mole_fractions(self, species_names=None, custom_coordinates=None, colors=None, figsize=(10, 5), filename=None)#

Plots the computed species mole fractions, \(\mathbf{X}_i\).

Example:

../_images/species-mole-fractions.svg
Parameters
  • species_names – (optional) list of str specifying the species names.

  • custom_coordinates – (optional) dict specifying the custom coordinates on the \(x\) axis. It should be of the format { label : array } where label is a string that will be plotted as an \(x\) axis label and array is a vector defining the custom grid. If not specified, a generic uniform grid between 0 and 1 will be used.

  • colors – (optional) list of str 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 to None, plot will not be saved to a file.

Returns

  • plot_handle - matplotlib.pyplot plot handle.

Composition.plot_species_mass_fractions#

multipy.composition.Composition.plot_species_mass_fractions(self, species_names=None, custom_coordinates=None, colors=None, figsize=(10, 5), filename=None)#

Plots the computed species mass fractions, \(\mathbf{Y}_i\).

Example:

../_images/species-mass-fractions.svg
Parameters
  • species_names – (optional) list of str specifying the species names.

  • custom_coordinates – (optional) dict specifying the custom coordinates on the \(x\) axis. It should be of the format { label : array } where label is a string that will be plotted as an \(x\) axis label and array is a vector defining the custom grid. If not specified, a generic uniform grid between 0 and 1 will be used.

  • colors – (optional) list of str 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 to None, plot will not be saved to a file.

Returns

  • plot_handle - matplotlib.pyplot plot handle.

Composition.plot_species_volume_fractions#

multipy.composition.Composition.plot_species_volume_fractions(self, species_names=None, custom_coordinates=None, colors=None, figsize=(10, 5), filename=None)#

Plots the computed species volume fractions, \(\mathbf{V}_i\).

Example:

../_images/species-volume-fractions.svg
Parameters
  • species_names – (optional) list of str specifying the species names.

  • custom_coordinates – (optional) dict specifying the custom coordinates on the \(x\) axis. It should be of the format { label : array } where label is a string that will be plotted as an \(x\) axis label and array is a vector defining the custom grid. If not specified, a generic uniform grid between 0 and 1 will be used.

  • colors – (optional) list of str 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 to None, plot will not be saved to a file.

Returns

  • plot_handle - matplotlib.pyplot plot handle.

Composition.plot_grad_species_mole_fractions#

multipy.composition.Composition.plot_grad_species_mole_fractions(self, species_names=None, custom_coordinates=None, colors=None, figsize=(10, 5), filename=None)#

Plots the computed species mole fraction gradients, \(\nabla \mathbf{X}_i\).

Example:

../_images/species-mole-fraction-gradients.svg
Parameters
  • species_names – (optional) list of str specifying the species names.

  • custom_coordinates – (optional) dict specifying the custom coordinates on the \(x\) axis. It should be of the format { label : array } where label is a string that will be plotted as an \(x\) axis label and array is a vector defining the custom grid. If not specified, a generic uniform grid between 0 and 1 will be used.

  • colors – (optional) list of str 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 to None, plot will not be saved to a file.

Returns

  • plot_handle - matplotlib.pyplot plot handle.

Composition.plot_grad_species_mass_fractions#

multipy.composition.Composition.plot_grad_species_mass_fractions(self, species_names=None, custom_coordinates=None, colors=None, figsize=(10, 5), filename=None)#

Plots the computed species mass fraction gradients, \(\nabla \mathbf{Y}_i\).

Example:

../_images/species-mass-fraction-gradients.svg
Parameters
  • species_names – (optional) list of str specifying the species names.

  • custom_coordinates – (optional) dict specifying the custom coordinates on the \(x\) axis. It should be of the format { label : array } where label is a string that will be plotted as an \(x\) axis label and array is a vector defining the custom grid. If not specified, a generic uniform grid between 0 and 1 will be used.

  • colors – (optional) list of str 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 to None, plot will not be saved to a file.

Returns

  • plot_handle - matplotlib.pyplot plot handle.

Composition.plot_grad_species_volume_fractions#

multipy.composition.Composition.plot_grad_species_volume_fractions(self, species_names=None, custom_coordinates=None, colors=None, figsize=(10, 5), filename=None)#

Plots the computed species volume fraction gradients, \(\nabla \mathbf{V}_i\).

Example:

../_images/species-volume-fraction-gradients.svg
Parameters
  • species_names – (optional) list of str specifying the species names.

  • custom_coordinates – (optional) dict specifying the custom coordinates on the \(x\) axis. It should be of the format { label : array } where label is a string that will be plotted as an \(x\) axis label and array is a vector defining the custom grid. If not specified, a generic uniform grid between 0 and 1 will be used.

  • colors – (optional) list of str 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 to None, plot will not be saved to a file.

Returns

  • plot_handle - matplotlib.pyplot plot handle.


Composition.species_mole_fractions#

multipy.composition.Composition.species_mole_fractions(self, species_molar_densities, T, p)#

Computes the species mole fractions:

\[X_i = \frac{c_i}{c}\]

under the assumption of the ideal gas law, where the mixture molar density, \(c\), is constant at constant temperature, \(T\), and pressure, \(p\):

\[c = p/RT\]
Parameters
  • species_molar_densities – scalar numpy.ndarray specifying the species molar densities, \(\mathbf{c}_i\), in \([mole/m^3]\). It should be of size (n_species,n_observations). Note that n_species can be equal to 1, if you are computing the mole fraction for only one of the species in the mixture.

  • Tint or float specifying the temperature, \(T\), in \([K]\).

  • pint or float specifying the pressure, \(p\), in \([Pa]\).

Returns

  • species_mole_fractions - scalar numpy.ndarray of species mole fractions, \(\mathbf{X}_i\), in \([-]\). It has size (n_species,n_observations).

Composition.species_mass_fractions#

multipy.composition.Composition.species_mass_fractions(self, species_mass_densities, mixture_mass_density)#

Computes the species mass fractions:

\[Y_i = \frac{\rho_i}{\rho}\]
Parameters
  • species_mass_densities – scalar numpy.ndarray specifying the species mass densities, \(\pmb{\rho}_i\), in \([kg/m^3]\). It should be of size (n_species,n_observations). Note that n_species can be equal to 1, if you are computing the mass fraction for only one of the species in the mixture.

  • mixture_mass_density – scalar numpy.ndarray specifying the mixture mass density, \(\pmb{\rho}\), in \([kg/m^3]\). It should be of size (1,n_observations).

Returns

  • species_mass_fractions - scalar numpy.ndarray of species mass fractions, \(\mathbf{Y}_i\), in \([-]\). It has size (n_species,n_observations).

Composition.species_volume_fractions#

multipy.composition.Composition.species_volume_fractions(self, species_volume, mixture_volume)#

Computes the species mass fractions:

\[V_i = \frac{v_i}{V}\]
Parameters
  • species_volume – scalar numpy.ndarray specifying the species volumes, \(v_i\), in \([m^3]\). It should be of size (n_species,n_observations). Note that n_species can be equal to 1, if you are computing the volume fraction for only one of the species in the mixture.

  • mixture_volume – scalar numpy.ndarray specifying the mixture volume, \(V\), in \([m^3]\). It should be of size (1,n_observations).

Returns

  • species_volume_fractions - scalar numpy.ndarray of species volume fractions, \(\mathbf{V}_i\), in \([-]\). It has size (n_species,n_observations).

Composition.grad_species_mole_fractions#

multipy.composition.Composition.grad_species_mole_fractions(self, species_mole_fractions, delta, edge_order=1)#

Computes species mole fraction gradients, \(\nabla \mathbf{X}_i\), numerically from the species mole fractions, \(X_i\):

\[\nabla X_i = \partial_d X_i\]

assuming a uniform grid spacing in a spatial dimension \(d\).

Note

numpy.gradient is used to compute gradients. Gradients are computed using second order accurate central differences in the interior points and either first or second order accurate (forward or backward) differences at the boundaries.

Parameters
  • species_mole_fractions – scalar numpy.ndarray specifying the species mole fractions, \(\mathbf{X}_i\), in \([-]\). It should be of size (n_species,n_observations).

  • deltaint or float specifying the grid spacing in the dimension \(d\).

  • edge_orderint specifying whether first or second order accurate differences are computed at the boundaries. It should be 1 or 2.

Returns

  • species_mole_fractions_gradients - vector numpy.ndarray of species mole fractions gradients, \(\nabla \mathbf{X}_i\), in \([-]\). It has size (n_species,n_observations).

Composition.grad_species_mass_fractions#

multipy.composition.Composition.grad_species_mass_fractions(self, species_mass_fractions, delta, edge_order=1)#

Computes species mass fraction gradients, \(\nabla \mathbf{Y}_i\), numerically from the species mass fractions, \(Y_i\):

\[\nabla Y_i = \partial_d Y_i\]

assuming a uniform grid spacing in a spatial dimension \(d\).

Note

numpy.gradient is used to compute gradients. Gradients are computed using second order accurate central differences in the interior points and either first or second order accurate (forward or backward) differences at the boundaries.

Parameters
  • species_mass_fractions – scalar numpy.ndarray specifying the species mass fractions, \(\mathbf{Y}_i\), in \([-]\). It should be of size (n_species,n_observations).

  • deltaint or float specifying the grid spacing in the dimension \(d\).

  • edge_orderint specifying whether first or second order accurate differences are computed at the boundaries. It should be 1 or 2.

Returns

  • species_mass_fractions_gradients - vector numpy.ndarray of species mass fractions gradients, \(\nabla \mathbf{Y}_i\), in \([-]\). It has size (n_species,n_observations).

Composition.grad_species_volume_fractions#

multipy.composition.Composition.grad_species_volume_fractions(self, species_volume_fractions, delta, edge_order=1)#

Computes species volume fraction gradients, \(\nabla \mathbf{V}_i\), numerically from the species volume fractions, \(V_i\):

\[\nabla V_i = \partial_d V_i\]

assuming a uniform grid spacing in a spatial dimension \(d\).

Note

numpy.gradient is used to compute gradients. Gradients are computed using second order accurate central differences in the interior points and either first or second order accurate (forward or backward) differences at the boundaries.

Parameters
  • species_volume_fractions – scalar numpy.ndarray specifying the species volume fractions, \(\mathbf{V}_i\), in \([-]\). It should be of size (n_species,n_observations).

  • deltaint or float specifying the grid spacing in the dimension \(d\).

  • edge_orderint specifying whether first or second order accurate differences are computed at the boundaries. It should be 1 or 2.

Returns

  • species_volume_fractions_gradients - vector numpy.ndarray of species volume fractions gradients, \(\nabla \mathbf{V}_i\), in \([-]\). It has size (n_species,n_observations).

Composition.species_molar_densities#

multipy.composition.Composition.species_molar_densities(self, species_mole_fractions, T, p)#

Computes the species molar densities:

\[c_i = c X_i\]

under the assumption of the ideal gas law, where the mixture molar density, \(c\), is constant at constant temperature, \(T\), and pressure, \(p\):

\[c = p/RT\]
Parameters
  • species_mole_fractions – scalar numpy.ndarray specifying the species mole fractions \(X_i\) in \([-]\). It should be of size (n_observations,n_species). Note that n_species can be equal to 1, if you are computing the molar density for only one of the species in the mixture.

  • Tint or float specifying the temperature \(T\) in \([K]\).

  • pint or float specifying the pressure \(p\) in \([Pa]\).

Returns

  • species_molar_densities - scalar numpy.ndarray of species molar densities \(c_i\) in \([mole/m^3]\). It has size (n_observations,n_species).

Composition.species_mass_densities#

multipy.composition.Composition.species_mass_densities(self, species_molar_densities, species_molar_masses)#

Computes the species mass densities:

\[\mathbf{\rho}_i = c_i M_i\]
Parameters
  • species_molar_densities – scalar numpy.ndarray specifying the species molar densities, \(\mathbf{c}_i\), in \([mole/m^3]\). It should be of size (n_species,n_observations). Note that n_species can be equal to 1, if you are computing the mass density for only one of the species in the mixture.

  • species_molar_masses – scalar numpy.ndarray specifying the species molar masses, \(\mathbf{M}_i\), in \([kg/mole]\). It should be of size (n_species,1).

Returns

  • species_mass_densities - scalar numpy.ndarray of species mass densities, \(\pmb{\mathbf{\rho}}_i\), in \([kg/m^3]\). It has size (n_species,n_observations).


Composition.mixture_molar_density#

multipy.composition.Composition.mixture_molar_density(self, T, p)#

Computes the mixture molar density:

\[c = p/RT\]

under the assumption of the ideal gas law, where the mixture molar density, \(c\), is constant at constant temperature, \(T\), and pressure, \(p\).

Parameters
  • Tint or float specifying the temperature, \(T\), in \([K]\).

  • pint or float specifying the pressure, \(p\), in \([Pa]\).

Returns

  • mixture_molar_density - mixture molar density, \(c\), in \([mole/m^3]\).

Composition.mixture_molar_volume#

multipy.composition.Composition.mixture_molar_volume(self, T, p)#

Computes the mixture molar volume, \(\bar{V}\), from:

\[\bar{V} = \frac{1}{c}\]

under the assumption of the ideal gas law, where the mixture molar density, \(c\), is constant at constant temperature, \(T\), and pressure, \(p\):

\[c = p/RT\]
Parameters
  • T – (optional) int or float specifying the temperature, \(T\), in \([K]\).

  • p – (optional) int or float specifying the pressure, \(p\), in \([Pa]\).

Returns

  • mixture_molar_volume - mixture molar volume, \(\bar{V}\), in \([m^3/mole]\).

Composition.mixture_molar_mass#

multipy.composition.Composition.mixture_molar_mass(self, species_fractions, basis, species_molar_masses)#

Computes the mixture molar mass at each observation, either from the species mole fractions (when basis='molar'):

\[\mathbf{M} = \sum_{i=1}^{n} X_i M_i\]

or from species mass fractions (when basis='mass'):

\[\mathbf{M} = \Big( \sum_{i=1}^{n} Y_i M_i \Big)^{-1}\]
Parameters
  • species_fractions – scalar numpy.ndarray specifying all species mole (or mass) fractions, \(X_i\) (or \(Y_i\)), in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

  • basisstr specifying whether the molar or mass equation should be used. Can be 'molar' or 'mass'.

  • species_molar_masses – scalar numpy.ndarray specifying the species molar masses, \(\mathbf{M}_i\), in \([kg/mole]\). It should be of size (n_species,1).

Returns

  • mixture_molar_mass - scalar numpy.ndarray of mixture molar masses, \(\pmb{M}\), in \([kg/mole]\). It has size (1,n_observations).

Composition.mixture_mass_density#

multipy.composition.Composition.mixture_mass_density(self, species_mass_densities)#

Computes the mixture mass density:

\[\rho = \sum_{i=1}^{n} \rho_i\]
Parameters

species_mass_densities – scalar numpy.ndarray specifying all species mass densities, \(\pmb{\rho}_i\), in \([kg/m^3]\). It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • mixture_mass_density - scalar numpy.ndarray of mixture mass density, \(\pmb{\rho}\), in \([kg/m^3]\). It has size (1,n_observations).