logoimage#

https://img.shields.io/badge/GitHub-multipy-blue.svg?style=flat http://img.shields.io/badge/license-MIT-blue.svg?style=flat https://readthedocs.org/projects/multipy-lib/badge/?version=latest https://mybinder.org/badge_logo.svg

Motivation#

Study hard what interests you the most in the most undisciplined, irreverent and original manner possible.

- Richard P. Feynman

This library is intended to support your learning of multicomponent mass transfer. The goal was to create a set of functions that are the essential building blocks with which you can play to get more intuition and understanding of quantities involved in multicomponent flows. With a fair amount of linear algebra and matrix operations, we paid special attention to document the equations along with the sizes of matrices. With these tools you can set-up your own problems such as the Stefan tube or the two-bulb diffusion without a whole lot of coding. We wish you a lot of joy in studying multicomponent mass transfer!


Credits#

Many materials helped in creating this library and we list them below. These are also great references to support your learning and broaden your knowledge!

  • [1] J. C. Sutherland - Multicomponent mass transfer course, CHEN-6603, The University of Utah, 2012

  • [2] R. Taylor, R. Krishna - Multicomponent mass transfer, Wiley, 1993

  • [3] R. B. Bird, W. E. Stewart, E. N. Lightfoot, D. J. Klingenberg - Introductory transport phenomena, Wiley, 2015

  • [4] H. M. Shey - Div, Grad, Curl, and All that: An Informal Text on Vector Calculus, W.W. Norton & Company, 2005

  • [5] A. Ern, V. Giovangigli - Multicomponent transport algorithms, Springer Science & Business Media, 1994

  • [6] V. Giovangigli - Multicomponent flow modeling, Birkhäuser Boston, 1999

  • [7] H. Brenner - Kinematics of volume transport, Physica A: Statistical Mechanics and its Applications, 349(1-2), pp.11-59

  • [8] H. Brenner - Navier-Stokes revisited, Physica A: Statistical Mechanics and its Applications, 349(1-2), pp.60-132

  • [9] H. T. Cullinan, Jr. - Analysis of the flux equations of multicomponent diffusion, Industrial & Engineering Chemistry Fundamentals, 4(2) (1965) 133-139

  • [10] J. C. Sutherland, C. A. Kennedy - Improved boundary conditions for viscous, reacting, compressible flows, Journal of Computational Physics 191 (2003) 502-524

  • [11] M. A. Hansen, J. C. Sutherland - On the consistency of state vectors and Jacobian matrices, Combustion and Flame 193 (2018) 257-271

  • [12] J. C. Sutherland, P. J. Smith, J. H. Chen - Quantification of differential diffusion in nonpremixed systems, Combustion Theory and Modelling 9(2) (2005) 365-383

  • [13] G. Sanderson - The diffusion equation


Documentation & User Guide#

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).

Velocities#

Class Velocity#

class multipy.Velocity(species_velocities=None)#

Supports computing and storing velocities.

Species convective velocities:

  • species velocities, \(\mathbf{u}_i\)

Mixture-averaged convective velocities:

  • mass-averaged velocity, \(\mathbf{v}\)

  • molar-averaged velocity, \(\mathbf{u}\)

  • volume-averaged velocity, \(\mathbf{u}^V\)

  • arbitrarily-averaged velocity, \(\mathbf{u}^a\)

Parameters

species_velocities – (optional) vector numpy.ndarray specifying the species velocities \(\mathbf{u}_i\) in \([m/s]\). It should be of size (n_observations,n_species) where n_species is at least 2.

Getters:

  • get_species_velocities

  • get_molar_averaged (is set to None at class init)

  • get_mass_averaged (is set to None at class init)

  • get_volume_averaged (is set to None at class init)

  • get_arbitrarily_averaged (is set to None at class init)

Setters:

  • set_species_velocities setter for get_species_velocities

  • set_molar_averaged setter for get_molar_averaged

  • set_mass_averaged setter for get_mass_averaged

  • set_volume_averaged setter for get_volume_averaged

  • set_arbitrarily_averaged setter for get_arbitrarily_averaged

Velocity.plot_species_velocities#
multipy.velocity.Velocity.plot_species_velocities(self, species_names=None, colors=None, figsize=(10, 5), filename=None)#

Plots the species velocities.

Example:

_images/stefan-tube-species-velocities.svg
Parameters
  • species_names – (optional) list of str specifying the species names.

  • 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.

Velocity.plot_averaged_velocities#
multipy.velocity.Velocity.plot_averaged_velocities(self, colors=None, figsize=(10, 5), filename=None)#

Plots the averaged velocities.

Example:

_images/stefan-tube-averaged-velocities.svg
Parameters
  • colors – (optional) list of str specifying the plotting colors for each averaged velocity. Example: colors=['#222222', '#858585'].

  • 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.

Velocity.species_velocities#
multipy.velocity.Velocity.species_velocities(self, total_flux, species_fractions, basis='molar', mixture_molar_density=None, mixture_mass_density=None)#

Computes the species velocities, \(\mathbf{u}_i\).

If basis is set to 'molar', species velocities are computed using the total molar fluxes, \(\mathbf{N}_i\), species mole fractions, \(X_i\) and the mixture molar density, \(c\):

\[\mathbf{u}_i = \frac{\mathbf{N}_i}{c X_i}\]

If basis is set to 'mass', species velocities are computed using the total mass fluxes, \(\mathbf{n}_i\), species mass fractions, \(Y_i\) and the mixture mass density, \(\rho\):

\[\mathbf{u}_i = \frac{\mathbf{n}_i}{\rho Y_i}\]
Parameters
  • total_flux – vector numpy.ndarray of total molar fluxes \(\mathbf{N}_i\) in \([mole/(m^2s)]\) or total mass fluxes \(\mathbf{n}_i\) in \([kg/(m^2s)]\). It should be of size (n_species, n_observations).

  • species_fractions – scalar numpy.ndarray specifying the species mole fractions \(X_i\) in \([-]\) or species mass fractions \(Y_i\) in \([-]\). It should be of size (n_species, n_observations).

  • basis – (optional) str specifying whether the molar or mass total flux equation should be used. Can be 'molar' or 'mass'.

  • mixture_molar_density – (optional) mixture molar density \(c\) in \([mole/m^3]\). Has to be specified if basis is set to molar.

  • mixture_mass_density – (optional) mixture mass density \(\rho\) in \([kg/m^3]\). Has to be specified if basis is set to mass.

Returns

  • species_velocities - vector numpy.ndarray of species velocities \(\mathbf{u}_i\) in \([m/s]\). It has size (n_species, n_observations).

Velocity.molar_averaged#
multipy.velocity.Velocity.molar_averaged(self, species_mole_fractions)#

Computes the molar-averaged velocity:

\[\mathbf{u} = \sum_{i=1}^{n} X_i \mathbf{u}_i\]

where \(n\) is the number of species.

Parameters

species_mole_fractions – scalar numpy.ndarray specifying all species mole fractions \(X_i\) in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • molar_averaged_velocity - vector numpy.ndarray of molar-averaged velocity \(\mathbf{u}\) in \([m/s]\). It has size (1,n_observations).

Velocity.mass_averaged#
multipy.velocity.Velocity.mass_averaged(self, species_mass_fractions)#

Computes the mass-averaged velocity:

\[\mathbf{v} = \sum_{i=1}^{n} Y_i \mathbf{u}_i\]

where \(n\) is the number of species.

Parameters

species_mass_fractions – scalar numpy.ndarray specifying all species mass fractions \(Y_i\) in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • mass_averaged_velocity - vector numpy.ndarray of mass-averaged velocity in \([m/s]\). It has size (1,n_observations).

Velocity.volume_averaged#
multipy.velocity.Velocity.volume_averaged(self, species_volume_fractions)#

Computes the volume-averaged velocity:

\[\mathbf{u}^V = \sum_{i=1}^{n} V_i \mathbf{u}_i\]

where \(n\) is the number of species.

Parameters

species_volume_fractions – scalar numpy.ndarray specifying all species volume fractions \(V_i\) in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • volume_averaged_velocity - vector numpy.ndarray of volume-averaged velocity \(\mathbf{u}^V\) in \([m/s]\). It has size (1,n_observations).

Velocity.arbitrarily_averaged#
multipy.velocity.Velocity.arbitrarily_averaged(self, arbitrary_weighting_factors)#

Computes the arbitrarily-averaged velocity:

\[\mathbf{u}^a = \sum_{i=1}^{n} a_i \mathbf{u}_i\]

where \(n\) is the number of species and \(a_i\) are the arbitrary weighting factors, such that \(\sum_{i=1}^{n} a_i = 1\).

Parameters

arbitrary_weighting_factors – scalar numpy.ndarray specifying arbitrary weighting factors, \(a_i\) in \([-]\), for all species. It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • arbitrarily_averaged_velocity - vector numpy.ndarray of arbitrarily-averaged velocity \(\mathbf{u}^a\) in \([m/s]\). It has size (1,n_observations).

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:

_images/stefan-tube-diffusive-flux-molar-diff-molar-avg.svg
Parameters
  • species_names – (optional) list of str specifying the species names.

  • 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.

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).

Diffusion#

Class Diffusion#

class multipy.diffusion.Diffusion(binary_diffusion_coefficients, species_names=None)#

Supports computing and storing diffusion-related quantities.

Parameters
  • binary_diffusion_coefficients – scalar numpy.ndarray specifying the binary diffusion coefficients, \(\pmb{\mathcal{D}}\), in \([m^2/s]\) for all species. It should be a symmetric matrix of size (n_species,n_species).

  • species_names – (optional) list of str specifying the names for all species. It should match the number and ordering of species as per the binary_diffusion_coefficients parameter. If not specified, species will be tagged with consecutive integers, i.e. '1', '2' and so on.

Getters:

  • get_binary_diffusion_coefficients (is set at class init)

  • get_n_species read only - (is set to 0 at class init)

  • get_fickian_diffusion_coefficients_molar_molar (is set to None at class init)

  • get_fickian_diffusion_coefficients_mass_mass (is set to None at class init)

  • get_fickian_diffusion_coefficients_molar_volume (is set to None at class init)

Setters:

  • set_binary_diffusion_coefficients setter for get_binary_diffusion_coefficients

  • set_fickian_diffusion_coefficients_molar_molar setter for get_fickian_diffusion_coefficients_molar_molar

  • set_fickian_diffusion_coefficients_mass_mass setter for get_fickian_diffusion_coefficients_mass_mass

  • set_fickian_diffusion_coefficients_molar_volume setter for get_fickian_diffusion_coefficients_molar_volume

Diffusion.print_binary_diffusion_coefficients#
multipy.diffusion.Diffusion.print_binary_diffusion_coefficients(self, table_format='pandas', float_format='.8f')#

Prints the binary diffusion coefficients matrix, \(\pmb{\mathcal{D}}\).

If table_format is set to 'pandas', a table will be printed in the pandas.DataFrame format:

_images/print_binary_diffusion_coefficients-pandas.png

If table_format is set to 'raw', a table will be printed in the raw text format:

|            | Acetone    | Methanol   | Air        |
| Acetone    | 0.0        | 8.48e-06   | 1.372e-05  |
| Methanol   | 8.48e-06   | 0.0        | 1.991e-05  |
| Air        | 1.372e-05  | 1.991e-05  | 0.0        |
Parameters
  • table_format – (optional) str specifying the printing format. It can be 'pandas' or 'raw'.

  • float_format – (optional) str specifying the float formatting when printing the table.

Diffusion.fickian_diffusion_coefficients_molar_molar#
multipy.diffusion.Diffusion.fickian_diffusion_coefficients_molar_molar(self, species_mole_fractions)#

Computes the molar Fickian diffusion coefficients expressed in a molar-averaged velocity reference frame, \(\mathbf{D}\), from:

\[\mathbf{D} = \mathbf{B}^{-1}\]

where the elements of \(\mathbf{B}\) are given by:

\[ \begin{align}\begin{aligned}B_{i, i} = \frac{X_i}{\mathcal{D}_{i, n}} + \sum_{j \neq i}^{n} \frac{X_j}{\mathcal{D}_{i,j}}\\B_{i, j} = - X_i \Big( \frac{1}{\mathcal{D}_{i, j}} - \frac{1}{\mathcal{D}_{i, n}} \Big)\end{aligned}\end{align} \]

and \(n\) is the number of species.

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).

Returns

  • fickian_diffusion_coefficients - scalar numpy.ndarray of molar Fickian diffusion coefficients expressed in a molar-averaged velocity reference frame, \(\mathbf{D}\), in \([m^2/s]\). It has size (n_species-1,n_species-1,n_observations).

Diffusion.fickian_diffusion_coefficients_mass_mass#
multipy.diffusion.Diffusion.fickian_diffusion_coefficients_mass_mass(self, species_mole_fractions, species_mass_fractions)#

Computes the mass Fickian diffusion coefficients expressed in a mass-averaged velocity reference frame, \(\mathbf{D}^o\), from:

\[\mathbf{D}^o =\]
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).

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

Returns

  • fickian_diffusion_coefficients - scalar numpy.ndarray of mass Fickian diffusion expressed in a mass-averaged velocity reference frame, \(\mathbf{D}^o\), in \([m^2/s]\). It has size (n_species-1,n_species-1,n_observations).

Diffusion.fickian_diffusion_coefficients_molar_volume#
multipy.diffusion.Diffusion.fickian_diffusion_coefficients_molar_volume(self, species_volume_fractions)#

Computes the molar Fickian diffusion coefficients expressed in a volume-averaged velocity reference frame, \(\mathbf{D}^V\), from:

\[\mathbf{D}^V =\]
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).

Returns

  • fickian_diffusion_coefficients - scalar numpy.ndarray of molar Fickian diffusion expressed in a volume-averaged velocity reference frame, \(\mathbf{D}^V\), in \([m^2/s]\). It has size (n_species-1,n_species-1,n_observations).

Diffusion.effective_diffusivity#
multipy.diffusion.Diffusion.effective_diffusivity(self, case=None)#

Computes the effective diffusivity using one of the selected assumptions:

Diffusion.diffusive_flux_molar_molar#
multipy.diffusion.Diffusion.diffusive_flux_molar_molar(self, grad_species_mole_fractions, mixture_molar_density)#

Computes the molar diffusive flux relative to a molar-averaged velocity, \(\mathbf{J}_i\), from the generalized Fick’s law:

\[\mathbf{J}_i = - c \mathbf{D} \nabla \mathbf{X}_i\]

where \(\mathbf{D}\) are the molar Fickian diffusion coefficients expressed in a molar-averaged velocity reference frame.

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

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

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).

Diffusion.diffusive_flux_mass_mass#
multipy.diffusion.Diffusion.diffusive_flux_mass_mass(self, grad_species_mass_fractions, mixture_mass_density)#

Computes the mass diffusive flux relative to a mass-averaged velocity, \(\mathbf{j}_i\), from the generalized Fick’s law:

\[\mathbf{j}_i = - \rho \mathbf{D}^o \nabla \mathbf{Y}_i\]

where \(\mathbf{D}^o\) are the mass Fickian diffusion coefficients expressed in a molar-averaged velocity reference frame.

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

  • 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

  • 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).

Diffusion.diffusive_flux_molar_volume#
multipy.diffusion.Diffusion.diffusive_flux_molar_volume(self, grad_species_molar_densities)#

Computes the molar diffusive flux relative to a volume-averaged velocity, \(\mathbf{J}_i^V\), from the generalized Fick’s law:

\[\mathbf{J}_i^V = - \mathbf{D}^V \nabla \mathbf{c}_i\]

where \(\mathbf{D}^V\) are the molar Fickian diffusion coefficients expressed in a volume-averaged velocity reference frame.

Parameters

grad_species_molar_densities – vector numpy.ndarray specifying the species molar densities gradients, \(\nabla \mathbf{c}_i\), in \([-]\). It should be of size (n_species,n_observations).

Returns

  • diffusive_flux - vector numpy.ndarray of molar diffusive fluxes relative to a volume-averaged velocity, \(\mathbf{J}_i^V\), in \([mole/(m^2s)]\). It has size (n_species,n_observations).

Transformations#

Class Transform#

class multipy.transform.Transform#

Supports performing transformations of multicomponent quantities to other bases or reference frames.

Transform.species_fractions_mole_to_mass#
multipy.transform.Transform.species_fractions_mole_to_mass(self, species_mole_fractions, species_molar_masses)#

Computes the species mass fractions, \(\mathbf{Y}_i\), from the species mole fractions, \(\mathbf{X}_i\), using the relation:

\[Y_i = \frac{M_i}{M} X_i\]
Parameters
  • species_mole_fractions – scalar numpy.ndarray specifying all species mole fractions, \(X_i\), in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

  • 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) where n_species is at least 2.

Returns

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

Transform.species_fractions_mass_to_mole#
multipy.transform.Transform.species_fractions_mass_to_mole(self, species_mass_fractions, species_molar_masses)#

Computes the species mole fractions, \(\mathbf{X}_i\), from the species mass fractions, \(\mathbf{Y}_i\), using the relation:

\[X_i = \frac{M}{M_i} Y_i\]
Parameters
  • species_mass_fractions – scalar numpy.ndarray specifying all species mass fractions, \(Y_i\), in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

  • 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) where n_species is at least 2.

Returns

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

Transform.species_gradients_mole_to_mass#
multipy.transform.Transform.species_gradients_mole_to_mass(self, species_mass_fractions, species_molar_masses)#

Computes an invertible, \(n-1\) dimensional transformation matrix, \(\mathbf{J}^{XY}\), that allows to transform from the species mole fraction gradients, \(\nabla \mathbf{X}_i\), to the species mass fraction gradients, \(\nabla \mathbf{Y}_i\), according to:

\[\nabla \mathbf{Y}_i = \mathbf{J}^{XY} \nabla \mathbf{X}_i\]

where:

\[J_{i,j}^{XY} = \frac{M_i}{M} \Bigg( \delta_{i,j} + \frac{Y_i}{M_i} (M_n - M_j) \Bigg)\]

Note

\(\mathbf{J}^{XY} = (\mathbf{J}^{YX})^{-1}\).

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

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

Returns

  • transformation_matrix - scalar numpy.ndarray transformation matrix, \(\mathbf{J}^{XY}\), in \([-]\). It has size (n_species-1,n_species-1,n_observations).

Transform.species_gradients_mass_to_mole#
multipy.transform.Transform.species_gradients_mass_to_mole(self, species_mass_fractions, species_molar_masses)#

Computes an invertible, \(n-1\) dimensional transformation matrix, \(\mathbf{J}^{YX}\), that allows to transform from the species mass fraction gradients, \(\nabla \mathbf{Y}_i\), to the species mole fraction gradients, \(\nabla \mathbf{X}_i\), according to:

\[\nabla \mathbf{X}_i = \mathbf{J}^{YX} \nabla \mathbf{Y}_i\]

where:

\[J_{i,j}^{YX} = \frac{M}{M_i} \Bigg( \delta_{i,j} + M Y_i \Big( \frac{1}{M_n} - \frac{1}{M_j} \Big) \Bigg)\]

Note

\(\mathbf{J}^{YX} = (\mathbf{J}^{XY})^{-1}\).

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

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

Returns

  • transformation_matrix - scalar numpy.ndarray transformation matrix, \(\mathbf{J}^{YX}\), in \([-]\). It has size (n_species-1,n_species-1,n_observations).

Transform.diffusive_flux_mass_mass_to_mass_molar#
multipy.transform.Transform.diffusive_flux_mass_mass_to_mass_molar(self, species_mole_fractions, species_mass_fractions)#

Computes an invertible, \(n-1\) dimensional transformation matrix, \(\mathbf{B}^{uo}\), that allows to transform from the mass diffusive flux relative to a mass-averaged velocity, \(\mathbf{j}_i\), to the mass diffusive flux relative to a molar-averaged velocity, \(\mathbf{j}_i^u\), according to:

\[\mathbf{j}_i^u = \mathbf{B}^{uo} \mathbf{j}_i\]

where:

\[B_{i,j}^{uo} = \delta_{i,j} - Y_i \Big( \frac{X_j}{Y_j} - \frac{X_n}{Y_n} \Big)\]

Note

\(\mathbf{B}^{uo} = (\mathbf{B}^{ou})^{-1}\).

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

  • species_mass_fractions – scalar numpy.ndarray specifying all species mass fractions, \(\mathbf{Y}_i\), in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • transformation_matrix - scalar numpy.ndarray transformation matrix, \(\mathbf{B}^{uo}\), in \([-]\). It has size (n_species-1,n_species-1,n_observations).

Transform.diffusive_flux_mass_molar_to_mass_mass#
multipy.transform.Transform.diffusive_flux_mass_molar_to_mass_mass(self, species_mole_fractions, species_mass_fractions)#

Computes an invertible, \(n-1\) dimensional transformation matrix, \(\mathbf{B}^{ou}\), that allows to transform from the mass diffusive flux relative to a molar-averaged velocity, \(\mathbf{j}_i^u\), to the mass diffusive flux relative to a mass-averaged velocity, \(\mathbf{j}_i\), according to:

\[\mathbf{j}_i = \mathbf{B}^{ou} \mathbf{j}_i^u\]

where:

\[B_{i,j}^{ou} = \delta_{i,j} - Y_i \Big( 1 - \frac{Y_n X_j}{X_n Y_j} \Big)\]

Note

\(\mathbf{B}^{ou} = (\mathbf{B}^{uo})^{-1}\).

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

  • species_mass_fractions – scalar numpy.ndarray specifying all species mass fractions, \(\mathbf{Y}_i\), in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • transformation_matrix - scalar numpy.ndarray transformation matrix, \(\mathbf{B}^{ou}\), in \([-]\). It has size (n_species-1,n_species-1,n_observations).

Transform.diffusive_flux_molar_molar_to_molar_volume#
multipy.transform.Transform.diffusive_flux_molar_molar_to_molar_volume(self, T, p, species_mole_fractions, species_partial_molar_volumes)#

Computes an invertible, \(n-1\) dimensional transformation matrix, \(\mathbf{B}^{Vu}\), that allows to transform from the molar diffusive flux relative to a molar-averaged velocity, \(\mathbf{J}_i\), to the molar diffusive flux relative to a volume-averaged velocity, \(\mathbf{J}_i^V\), according to:

\[\mathbf{J}_i^V = \mathbf{B}^{Vu} \mathbf{J}_i\]

where:

\[B_{i,j}^{Vu} = \delta_{i,j} - X_i (\bar{V}_j - \bar{V}_n) / \bar{V}\]

Note

\(\mathbf{B}^{Vu} = (\mathbf{B}^{uV})^{-1}\).

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

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

  • species_mole_fractions – scalar numpy.ndarray specifying all species mole fractions, \(\mathbf{X}_i\), in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

  • species_partial_molar_volumes – scalar numpy.ndarray specifying all species partial molar volumes, \(\bar{\mathbf{V}}_i\), in \([m^3/mole]\). It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • transformation_matrix - scalar numpy.ndarray transformation matrix, \(\mathbf{B}^{Vu}\), in \([-]\). It has size (n_species-1,n_species-1,n_observations).

Transform.diffusive_flux_molar_volume_to_molar_molar#
multipy.transform.Transform.diffusive_flux_molar_volume_to_molar_molar(self, species_mole_fractions, species_partial_molar_volumes)#

Computes an invertible, \(n-1\) dimensional transformation matrix, \(\mathbf{B}^{uV}\), that allows to transform from the molar diffusive flux relative to a volume-averaged velocity, \(\mathbf{J}_i^V\), to the molar diffusive flux relative to a molar-averaged velocity, \(\mathbf{J}_i\), according to:

\[\mathbf{J}_i = \mathbf{B}^{uV} \mathbf{J}_i^V\]

where:

\[B_{i,j}^{uV} = \delta_{i,j} - X_i (1 - \bar{V}_j / \bar{V}_n)\]

Note

\(\mathbf{B}^{uV} = (\mathbf{B}^{Vu})^{-1}\).

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

  • species_partial_molar_volumes – scalar numpy.ndarray specifying all species partial molar volumes, \(\bar{\mathbf{V}}_i\), in \([m^3/mole]\). It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • transformation_matrix - scalar numpy.ndarray transformation matrix, \(\mathbf{B}^{uV}\), in \([-]\). It has size (n_species-1,n_species-1,n_observations).

Transform.fickian_diffusion_coefficients_molar_molar_to_mass_mass#
multipy.transform.Transform.fickian_diffusion_coefficients_molar_molar_to_mass_mass(self, species_mole_fractions, species_mass_fractions)#

Computes an invertible, \(n-1\) dimensional transformation matrix, \((\mathbf{B}^{uo})^{-1} \mathrm{diag}(\mathbf{Y}_i) (\mathrm{diag}(\mathbf{X}_i))^{-1}\), that allows to transform from the molar Fickian diffusion coefficients in the molar-averaged velocity reference frame, \(\mathbf{D}\), to the mass Fickian diffusion coefficients in the volume-averaged velocity reference frame, \(\mathbf{D}^o\), according to:

\[\mathbf{D}^o = (\mathbf{B}^{uo})^{-1} \mathrm{diag}(\mathbf{Y}_i) (\mathrm{diag}(\mathbf{X}_i))^{-1} \mathbf{D} \mathrm{diag}(\mathbf{X}_i) (\mathrm{diag}(\mathbf{Y}_i))^{-1} \mathbf{B}^{uo}\]

where:

\[B_{i,j}^{uo} = \delta_{i,j} - Y_i \Big( \frac{X_j}{Y_j} - \frac{X_n}{Y_n} \Big)\]

\(\mathrm{diag}(\mathbf{X}_i)\) and \(\mathrm{diag}(\mathbf{Y}_i)\) are diagonal matrices whose non-zero entries are the mole or mass fractions respectively of \(n-1\) species.

Note

\((\mathbf{B}^{uo})^{-1} \mathrm{diag}(\mathbf{Y}_i) (\mathrm{diag}(\mathbf{X}_i))^{-1} = \Big( \mathrm{diag}(\mathbf{X}_i) (\mathrm{diag}(\mathbf{Y}_i))^{-1} \mathbf{B}^{uo} \Big)^{-1}\).

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

  • species_mass_fractions – scalar numpy.ndarray specifying all species mass fractions, \(\mathbf{Y}_i\), in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • transformation_matrix - scalar numpy.ndarray transformation matrix, \((\mathbf{B}^{uo})^{-1} \mathrm{diag}(\mathbf{Y}_i) (\mathrm{diag}(\mathbf{X}_i))^{-1}\), in \([-]\). It has size (n_species-1,n_species-1,n_observations).

Transform.fickian_diffusion_coefficients_molar_molar_to_molar_volume#
multipy.transform.Transform.fickian_diffusion_coefficients_molar_molar_to_molar_volume(self, T, p, species_mole_fractions, species_partial_molar_volumes)#

Computes an invertible, \(n-1\) dimensional transformation matrix, \(\mathbf{B}^{Vu}\), that allows to transform the molar Fickian diffusion coefficients from the molar-averaged velocity reference frame, \(\mathbf{D}\), to the volume-averaged velocity reference frame, \(\mathbf{D}^V\), according to:

\[\mathbf{D}^V = \mathbf{B}^{Vu} \mathbf{D} (\mathbf{B}^{Vu})^{-1}\]

where:

\[B_{i,j}^{Vu} = \delta_{i,j} - X_i (\bar{V}_j - \bar{V}_n) / \bar{V}\]

Note

\(\mathbf{B}^{Vu} = (\mathbf{B}^{uV})^{-1}\).

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

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

  • species_mole_fractions – scalar numpy.ndarray specifying all species mole fractions, \(\mathbf{X}_i\), in \([-]\). It should be of size (n_species,n_observations) where n_species is at least 2.

  • species_partial_molar_volumes – scalar numpy.ndarray specifying all species partial molar volumes, \(\bar{\mathbf{V}}_i\), in \([m^3/mole]\). It should be of size (n_species,n_observations) where n_species is at least 2.

Returns

  • transformation_matrix - scalar numpy.ndarray transformation matrix, \(\mathbf{B}^{Vu}\), in \([-]\). It has size (n_species-1,n_species-1,n_observations).

Basic checks#

Class Check#

class multipy.check.Check#

Supports performing basic checks of the computed quantities.

Check.sum_of_species_fractions#
multipy.check.Check.sum_of_species_fractions(self, species_fractions, tolerance=1e-12, verbose=False)#

Checks if all species mole/mass/volume fractions sum to 1.0 for every observation within a specified tolerance.

For mole fractions:

\[\sum_{i=1}^{n} X_i = 1.0\]

For mass fractions:

\[\sum_{i=1}^{n} Y_i = 1.0\]

For volume fractions:

\[\sum_{i=1}^{n} V_i = 1.0\]

where \(n\) is the number of species.

Parameters
  • species_fractions – scalar numpy.ndarray specifying all species mole/mass/volume fractions in \([-]\). It should be of size (n_species, n_observations) where n_species is at least 2.

  • tolerance – (optional) float specifying the tolerance. It should be larger than 0.0 and smaller than 1.0.

  • verbose – (optional) bool for printing verbose information.

Returns

  • idx - indices of observations where species mole/mass/volume fractions do not sum to 1.0 within a specified tolerance.

Check.range_of_species_fractions#
multipy.check.Check.range_of_species_fractions(self, species_fractions, tolerance=1e-12, verbose=False)#

Checks if all species mole/mass/volume fraction values are bounded between 0 and 1.

For mole fractions:

\[X_i \in \langle 0, 1 \rangle\]

For mass fractions:

\[Y_i \in \langle 0, 1 \rangle\]

For volume fractions:

\[V_i \in \langle 0, 1 \rangle\]
Parameters
  • species_fractions – scalar numpy.ndarray specifying all species mole/mass/volume fractions in \([-]\). It should be of size (n_observations,n_species) where n_species is at least 2.

  • verbose – (optional) bool for printing verbose information.

Returns

  • idx_below_zero - indices of observations where species mole/mass/volume fractions are less than 0.0 within a specified tolerance.

  • idx_above_one - indices of observations where species mole/mass/volume fractions are larger than 1.0 within a specified tolerance.

Check.sum_of_species_gradients#
multipy.check.Check.sum_of_species_gradients(self, species_gradients, tolerance=1e-12, verbose=False)#

Checks if all species mole/mass/volume fraction gradients sum to 0.0 for every observation within a specified tolerance.

For mole fractions:

\[\sum_{i=1}^{n} \nabla X_i = 0.0\]

For mass fractions:

\[\sum_{i=1}^{n} \nabla Y_i = 0.0\]

For volume fractions:

\[\sum_{i=1}^{n} \nabla V_i = 0.0\]

where \(n\) is the number of species.

Parameters
  • species_gradients – scalar numpy.ndarray specifying all species mole/mass/volume fraction gradients in \([-]\). It should be of size (n_species, n_observations) where n_species is at least 2.

  • tolerance – (optional) float specifying the tolerance. It should be larger than 0.0 and smaller than 1.0.

  • verbose – (optional) bool for printing verbose information.

Returns

  • idx - indices of observations where species mole/mass/volume fraction gradients do not sum to 0.0 within a specified tolerance.

Check.sum_of_species_production_rates#
multipy.check.Check.sum_of_species_production_rates(self, species_production_rates, tolerance=1e-12, verbose=False)#

Checks if all species production rates sum to 0.0 for every observation within a specified tolerance:

For net molar production rates:

\[\sum_{i=1}^{n} s_i = 0.0\]

For net mass production rates:

\[\sum_{i=1}^{n} \omega_i = 0.0\]

where \(n\) is the number of species.

Parameters
  • species_production_rates – scalar numpy.ndarray specifying all species production rates, \(s_i\) in \(mole/(m^3s)\) or \(\omega_i\) in \(kg/(m^3s)\). It should be of size (n_species,n_observations) where n_species is at least 2.

  • tolerance – (optional) float specifying the tolerance. It should be larger than 0.0 and smaller than 1.0.

  • verbose – (optional) bool for printing verbose information.

Returns

  • idx - indices of observations where species source terms do not sum to 0.0 within a specified tolerance.

Solution templates#

Class Templates#

class multipy.templates.Templates#

Contains solution templates and building blocks for common multicomponent mass trasfer problems.

Templates.stefan_diffusion#
multipy.templates.Templates.stefan_diffusion(self, alpha, beta)#

Computes a matrix \(\pmb{\Phi}\) and a vector \(\pmb{\phi}\) such that:

\[ \begin{align}\begin{aligned}\Phi_{i, i} = \frac{\alpha_i}{\beta_{i, n}} + \sum_{j \neq i}^{n} \frac{\alpha_j}{\beta_{i,j}}\\\Phi_{i, j} = - \alpha_i \Big( \frac{1}{\beta_{i, j}} - \frac{1}{\beta_{i, n}} \Big)\end{aligned}\end{align} \]

and

\[\phi_i = - \frac{\alpha_i}{\beta_{i, n}}\]

where \(n\) is the number of species and \(\alpha_i\) and \(\beta_{i,j}\) are free user-specified coefficients.

This template can be used in Stefan diffusion -type problems.

Parameters
  • alpha – scalar numpy.ndarray specifying the cofficients \(\alpha_{i}\). It should be of size (n_species,1).

  • beta – scalar numpy.ndarray specifying the coefficients \(\beta_{i,j}\). It should be of size (n_species,n_species).

Returns

  • phi - scalar numpy.ndarray \(\pmb{\phi}\). It has size (n_species-1,1).

  • Phi - scalar numpy.ndarray \(\pmb{\Phi}\). It has size (n_species-1,n_species-1).

Notation#

This part of the documentation presents notations for all multipy quantities.


Indexing#

  • \(i\) denotes that a quantity is related to an \(i^{th}\) component of the mixture.

  • \(n\) represents the number of components (species) of the mixture.


Bases and reference frames#

Multicomponent quantities can be written in various bases (e.g. mass or molar) and in reference to different average velocity of the mixture (e.g. mass-averaged or molar-averged). To keep track of a basis and a reference frame related to mixture velocity, we adopted the following notation in naming functions or parameters:

_images/notation-reference-frames.svg

As an example, Flux.diffusive_mass_molar should be interpreted as a mass diffusive flux relative to a molar-averaged velocity.


Shape of a general multicomponent quantity matrix#

We assume that a matrix \(\pmb{\phi}_i \in \mathcal{R}^{n \times N}\) describes any multicomponent quantity, where \(n\) is the number of components (species) in the mixture and \(N\) is the number of observations of that quantity (for instance spatial positions). Mixture components are thus stored in rows and their observations are stored in columns:

_images/notation-matrix-shape.svg

Vector calculus primer#

Gradient#

For a scalar function \(f\), the gradient of \(f\) is denoted:

\[\nabla f\]

In three dimensions we can compute this as:

\[\begin{split}\begin{gather} f \begin{bmatrix} \frac{\partial}{\partial x} \\ \frac{\partial}{\partial y} \\ \frac{\partial}{\partial z} \\ \end{bmatrix} = \begin{bmatrix} \frac{\partial f}{\partial x} \\ \frac{\partial f}{\partial y} \\ \frac{\partial f}{\partial z} \\ \end{bmatrix} \end{gather}\end{split}\]

This results in a vector.

For a vector function \(\mathbf{f}\), the gradient of \(\mathbf{f}\) is denoted:

\[\nabla \mathbf{f}\]

In three dimensions we can compute this as:

\[\begin{split}\begin{gather} \begin{bmatrix} \frac{\partial}{\partial x} \\ \frac{\partial}{\partial y} \\ \frac{\partial}{\partial z} \\ \end{bmatrix} \cdot \begin{bmatrix} f_x & f_y & f_z \\ \end{bmatrix} = \begin{bmatrix} \frac{\partial f_x}{\partial x} & \frac{\partial f_y}{\partial x} & \frac{\partial f_z}{\partial x}\\ \frac{\partial f_x}{\partial y} & \frac{\partial f_y}{\partial y} & \frac{\partial f_z}{\partial y}\\ \frac{\partial f_x}{\partial z} & \frac{\partial f_y}{\partial z} & \frac{\partial f_z}{\partial z}\\ \end{bmatrix} \end{gather}\end{split}\]

This results in a tensor.

Divergence#

For a vector function \(\mathbf{f}\), the divergence of \(\mathbf{f}\) is denoted:

\[\nabla \cdot \mathbf{f}\]

In three dimensions we can compute this as:

\[\begin{split}\begin{gather} \begin{bmatrix} \frac{\partial}{\partial x} & \frac{\partial}{\partial y} & \frac{\partial}{\partial z} \end{bmatrix} \cdot \begin{bmatrix} f_x \\ f_y \\ f_z \\ \end{bmatrix} = \frac{\partial f_x}{\partial x} + \frac{\partial f_y}{\partial y} + \frac{\partial f_z}{\partial z} \end{gather}\end{split}\]

This results in a scalar.

Outer product#

The outer product between matrices A and B can be computed using numpy as:

numpy.outer(A,B)
Tensor contraction#

For two tensors, \(\mathbf{A}\) and \(\mathbf{B}\), tensor contraction (scalar product) is denoted:

\[\mathbf{A} : \mathbf{B}\]

The tensor contraction (scalar product) between matrices A and B can be computed using numpy as:

numpy.tensordot(A,B,axes=2)

which achieves the same thing as:

numpy.sum(numpy.multiply(A,B))
Various forms of the divergence theorem#
  • For a scalar field \(\phi\): \(\int_{S(t)} \phi \mathbf{a} dS = \int_{V(t)} \nabla \phi dV\)

  • For a vector field \(\mathbf{q}\): \(\int_{S(t)} \mathbf{q} \cdot \mathbf{a} dS = \int_{V(t)} \nabla \cdot \mathbf{q} dV\)

  • For a tensor field \(\pmb{\tau}\): \(\int_{S(t)} \pmb{\tau} \cdot \mathbf{a} dS = \int_{V(t)} \nabla \cdot \pmb{\tau} dV\)

Quantities#

This part of the documentation presents notations for all multipy quantities.


Composition#

Species mole fractions#

Notation

Unit

Code

Type

Shape

\(\mathbf{X}_i\)

\([-]\)

Composition.species_mole_fractions

scalar numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{X}_i = \begin{bmatrix} \cdots & X_1 & \cdots \\ \cdots & X_2 & \cdots \\ & \vdots & \\ \cdots & X_n & \cdots \\ \end{bmatrix}\end{split}\]
Species mole fractions gradients#

Notation

Unit

Code

Type

Shape

\(\nabla \mathbf{X}_i\)

\([-]\)

Composition.grad_species_mole_fractions

vector numpy.ndarray

(n_species,n_observations)

\[\begin{split}\nabla \mathbf{X}_i = \begin{bmatrix} \cdots & \nabla X_1 & \cdots \\ \cdots & \nabla X_2 & \cdots \\ & \vdots & \\ \cdots & \nabla X_n & \cdots \\ \end{bmatrix}\end{split}\]
Species mass fractions#

Notation

Unit

Code

Type

Shape

\(\mathbf{Y}_i\)

\([-]\)

Composition.species_mass_fractions

scalar numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{Y}_i = \begin{bmatrix} \cdots & Y_1 & \cdots \\ \cdots & Y_2 & \cdots \\ & \vdots & \\ \cdots & Y_n & \cdots \\ \end{bmatrix}\end{split}\]
Species mass fractions gradients#

Notation

Unit

Code

Type

Shape

\(\nabla \mathbf{Y}_i\)

\([-]\)

Composition.grad_species_mass_fractions

vector numpy.ndarray

(n_species,n_observations)

\[\begin{split}\nabla \mathbf{Y}_i = \begin{bmatrix} \cdots & \nabla Y_1 & \cdots \\ \cdots & \nabla Y_2 & \cdots \\ & \vdots & \\ \cdots & \nabla Y_n & \cdots \\ \end{bmatrix}\end{split}\]
Species volume fractions#

Notation

Unit

Code

Type

Shape

\(\mathbf{V}_i\)

\([-]\)

Composition.species_volume_fractions

scalar numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{V}_i = \begin{bmatrix} \cdots & V_1 & \cdots \\ \cdots & V_2 & \cdots \\ & \vdots & \\ \cdots & V_n & \cdots \\ \end{bmatrix}\end{split}\]
Species volume fractions gradients#

Notation

Unit

Code

Type

Shape

\(\nabla \mathbf{V}_i\)

\([-]\)

Composition.grad_species_volume_fractions

vector numpy.ndarray

(n_species,n_observations)

\[\begin{split}\nabla \mathbf{V}_i = \begin{bmatrix} \cdots & \nabla V_1 & \cdots \\ \cdots & \nabla V_2 & \cdots \\ & \vdots & \\ \cdots & \nabla V_n & \cdots \\ \end{bmatrix}\end{split}\]
Species molar densities#

Notation

Unit

Code

Type

Shape

\(\mathbf{c}_i\)

\([mole/m^3]\)

Composition.species_molar_densities

scalar numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{c}_i = \begin{bmatrix} \cdots & c_1 & \cdots \\ \cdots & c_2 & \cdots \\ & \vdots & \\ \cdots & c_n & \cdots \\ \end{bmatrix}\end{split}\]
Species mass densities#

Notation

Unit

Code

Type

Shape

\(\pmb{\rho}_i\)

\([kg/m^3]\)

Composition.species_mass_densities

scalar numpy.ndarray

(n_species,n_observations)

\[\begin{split}\pmb{\rho}_i = \begin{bmatrix} \cdots & \rho_1 & \cdots \\ \cdots & \rho_2 & \cdots \\ & \vdots & \\ \cdots & \rho_n & \cdots \\ \end{bmatrix}\end{split}\]
Species molar masses#

Notation

Unit

Code

Type

Shape

\(\pmb{M}_i\)

\([kg/mole]\)

Composition.species_molar_masses

scalar numpy.ndarray

(n_species,1)

\[\begin{split}\pmb{M}_i = \begin{bmatrix} M_1 \\ M_2 \\ \vdots \\ M_n \\ \end{bmatrix}\end{split}\]
Species molar production rates#

Notation

Unit

Code

Type

Shape

\(\mathbf{s}_i\)

\([mole/(m^3s)]\)

Composition.get_species_molar_production_rates

scalar numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{s}_i = \begin{bmatrix} \cdots & s_1 & \cdots \\ \cdots & s_2 & \cdots \\ & \vdots & \\ \cdots & s_n & \cdots \\ \end{bmatrix}\end{split}\]
Species mass production rates#

Notation

Unit

Code

Type

Shape

\(\pmb{\omega}_i\)

\([kg/(m^3s)]\)

Composition.get_species_mass_production_rates

scalar numpy.ndarray

(n_species,n_observations)

\[\begin{split}\pmb{\omega}_i = \begin{bmatrix} \cdots & \omega_1 & \cdots \\ \cdots & \omega_2 & \cdots \\ & \vdots & \\ \cdots & \omega_n & \cdots \\ \end{bmatrix}\end{split}\]
Mixture molar density#

Notation

Unit

Code

Type

Shape

\(c\)

\([mole/m^3]\)

Composition.mixture_molar_density

scalar float

N/A

Mixture molar volume#

Notation

Unit

Code

Type

Shape

\(\bar{V}\)

\([m^3/mole]\)

Composition.mixture_molar_volume

scalar float

N/A

Mixture mass density#

Notation

Unit

Code

Type

Shape

\(\pmb{\rho}\)

\([kg/m^3]\)

Composition.mixture_mass_density

scalar numpy.ndarray

(1,n_observations)

\[\begin{split}\pmb{\rho} = \begin{bmatrix} \cdots & \rho & \cdots \\ \end{bmatrix}\end{split}\]
Mixture molar mass#

Notation

Unit

Code

Type

Shape

\(\pmb{M}\)

\([kg/mole]\)

Composition.mixture_molar_mass

scalar numpy.ndarray

(1,n_observations)

\[\begin{split}\pmb{M} = \begin{bmatrix} \cdots & M & \cdots \\ \end{bmatrix}\end{split}\]

Velocity#

Species velocities#

Notation

Unit

Code

Type

Shape

\(\mathbf{u}_i\)

\([m/s]\)

Velocity.species_velocities

vector numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{u}_i = \begin{bmatrix} \cdots & \mathbf{u}_1 & \cdots \\ \cdots & \mathbf{u}_2 & \cdots \\ & \vdots & \\ \cdots & \mathbf{u}_n & \cdots \\ \end{bmatrix}\end{split}\]
Molar-averaged mixture velocity#

Notation

Unit

Code

Type

Shape

\(\mathbf{u}\)

\([m/s]\)

Velocity.molar_averaged

vector numpy.ndarray

(1,n_observations)

\[\begin{split}\mathbf{u} = \begin{bmatrix} \cdots & \mathbf{u} & \cdots \\ \end{bmatrix}\end{split}\]
Mass-averaged mixture velocity#

Notation

Unit

Code

Type

Shape

\(\mathbf{v}\)

\([m/s]\)

Velocity.mass_averaged

vector numpy.ndarray

(1,n_observations)

\[\begin{split}\mathbf{v} = \begin{bmatrix} \cdots & \mathbf{v} & \cdots \\ \end{bmatrix}\end{split}\]
Volume-averaged mixture velocity#

Notation

Unit

Code

Type

Shape

\(\mathbf{u}^V\)

\([m/s]\)

Velocity.volume_averaged

vector numpy.ndarray

(1,n_observations)

\[\begin{split}\mathbf{u}^V = \begin{bmatrix} \cdots & \mathbf{u}^V & \cdots \\ \end{bmatrix}\end{split}\]
Arbitrarily-averaged mixture velocity#

Notation

Unit

Code

Type

Shape

\(\mathbf{u}^a\)

\([m/s]\)

Velocity.arbitrarily_averaged

vector numpy.ndarray

(1,n_observations)

\[\begin{split}\mathbf{u}^a = \begin{bmatrix} \cdots & \mathbf{u}^a & \cdots \\ \end{bmatrix}\end{split}\]

Flux#

Total molar flux#

Notation

Unit

Code

Type

Shape

\(\mathbf{N}_i\)

\([mole/(m^2 s)]\)

Flux.total_molar_flux

vector numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{N}_i = \begin{bmatrix} \cdots & \mathbf{N}_1 & \cdots \\ \cdots & \mathbf{N}_2 & \cdots \\ & \vdots & \\ \cdots & \mathbf{N}_n & \cdots \\ \end{bmatrix}\end{split}\]
Total mass flux#

Notation

Unit

Code

Type

Shape

\(\mathbf{n}_i\)

\([kg/(m^2 s)]\)

Flux.total_mass_flux

vector numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{n}_i = \begin{bmatrix} \cdots & \mathbf{n}_1 & \cdots \\ \cdots & \mathbf{n}_2 & \cdots \\ & \vdots & \\ \cdots & \mathbf{n}_n & \cdots \\ \end{bmatrix}\end{split}\]
Molar diffusive flux relative to a molar-averaged velocity#

Notation

Unit

Code

Type

Shape

\(\mathbf{J}_i\)

\([mole/(m^2 s)]\)

Flux.diffusive_molar_molar

vector numpy.ndarray

(n_species,n_observations)

\(\mathbf{J}_i\)

\([mole/(m^2 s)]\)

Diffusion.diffusive_flux_molar_molar

vector numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{J}_i = \begin{bmatrix} \cdots & \mathbf{J}_1 & \cdots \\ \cdots & \mathbf{J}_2 & \cdots \\ & \vdots & \\ \cdots & \mathbf{J}_n & \cdots \\ \end{bmatrix}\end{split}\]
Molar diffusive flux relative to a mass-averaged velocity#

Notation

Unit

Code

Type

Shape

\(\mathbf{J}_i^v\)

\([mole/(m^2 s)]\)

Flux.diffusive_molar_mass

vector numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{J}_i^v = \begin{bmatrix} \cdots & \mathbf{J}_1^v & \cdots \\ \cdots & \mathbf{J}_2^v & \cdots \\ & \vdots & \\ \cdots & \mathbf{J}_n^v & \cdots \\ \end{bmatrix}\end{split}\]
Mass diffusive flux relative to a molar-averaged velocity#

Notation

Unit

Code

Type

Shape

\(\mathbf{j}_i^u\)

\([kg/(m^2 s)]\)

Flux.diffusive_mass_molar

vector numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{j}_i^u = \begin{bmatrix} \cdots & \mathbf{j}_1^u & \cdots \\ \cdots & \mathbf{j}_2^u & \cdots \\ & \vdots & \\ \cdots & \mathbf{j}_{n}^u & \cdots \\ \end{bmatrix}\end{split}\]
Mass diffusive flux relative to a mass-averaged velocity#

Notation

Unit

Code

Type

Shape

\(\mathbf{j}_i\)

\([kg/(m^2 s)]\)

Flux.diffusive_mass_mass

vector numpy.ndarray

(n_species,n_observations)

\(\mathbf{j}_i\)

\([kg/(m^2 s)]\)

Diffusion.diffusive_flux_mass_mass

vector numpy.ndarray

(n_species,n_observations)

\[\begin{split}\mathbf{j}_i = \begin{bmatrix} \cdots & \mathbf{j}_1 & \cdots \\ \cdots & \mathbf{j}_2 & \cdots \\ & \vdots & \\ \cdots & \mathbf{j}_{n} & \cdots \\ \end{bmatrix}\end{split}\]

Diffusion#

Binary diffusion coefficients#

Notation

Unit

Code

Type

Shape

\(\pmb{\mathcal{D}}\)

\([m^2/s]\)

Diffusion.get_binary_diffusion_coefficients

scalar numpy.ndarray

(n_species,n_species)

\[\begin{split}\pmb{\mathcal{D}} = \begin{bmatrix} - & \mathcal{D}_{1,2} & \dots & \mathcal{D}_{1,n} \\ \mathcal{D}_{2,1} & - & \dots & \mathcal{D}_{2,n} \\ \vdots & \vdots & \ddots & \vdots \\ \mathcal{D}_{n,1} & \mathcal{D}_{n,2} & \dots & - \\ \end{bmatrix}\end{split}\]

where \(\mathcal{D}_{i,j} = \mathcal{D}_{j,i} \,\,\, \forall_{i \neq j}\).

Molar Fickian diffusion coefficients#

Notation

Unit

Code

Type

Shape

\(\mathbf{D}\)

\([m^2/s]\)

Diffusion.fickian_diffusion_coefficients_molar_molar

scalar numpy.ndarray

(n_species-1,n_species-1,n_observations)

\[\begin{split}\mathbf{D} = \begin{bmatrix} D_{1,1} & D_{1,2} & \dots & D_{1,n-1} \\ D_{2,1} & D_{2,2} & \dots & D_{2,n-1} \\ \vdots & \vdots & \ddots & \vdots \\ D_{n-1,1} & D_{n-1,2} & \dots & D_{n-1,n-1} \\ \end{bmatrix}\end{split}\]

where, in general, \(D_{i,j} \neq D_{j,i}\) and \(D_{i,j} \neq 0 \,\,\, \forall_{i, j}\).

Governing equations#


Multicomponent mixture#

Let \(\mathbf{v}\) be the mass-averaged velocity of the mixture, defined as:

\[\mathbf{v} := \sum_{i = 1}^{n} Y_i \mathbf{u}_i\]

where:

  • \(Y_i\) is the mass fraction of species \(i\)

  • \(\mathbf{u}_i\) is the velocity of species \(i\)

  • \(n\) is the number of species (components) in the mixture

In an analogous way, we can define the molar-averaged velocity of the mixture, \(\mathbf{u}\), as:

\[\mathbf{u} := \sum_{i = 1}^{n} X_i \mathbf{u}_i\]

where:

  • \(X_i\) is the mole fraction of species \(i\)

  • \(\mathbf{u}_i\) is the velocity of species \(i\)

  • \(n\) is the number of species (components) in the mixture

At a given point in space and time, transport of physical quantities in a multicomponent mixture can be described by the following set of governing equations written in the conservative (strong) form:

Continuity equation#
\[\frac{\partial \rho}{\partial t} = - \nabla \cdot \rho \mathbf{v}\]

where:

  • \(\rho\) is the mixture mass density

  • \(\mathbf{v}\) is the mass-averaged velocity of the mixture

Species mass conservation equation#
\[\frac{\partial \rho Y_i}{\partial t} = - \nabla \cdot \rho Y_i \mathbf{v} - \nabla \cdot \mathbf{j}_i + \omega_i\]

where:

  • \(\rho\) is the mixture mass density

  • \(Y_i\) is the mass fraction of species \(i\)

  • \(\mathbf{v}\) is the mass-averaged velocity of the mixture

  • \(\mathbf{j}_i\) is the mass diffusive flux of species \(i\) relative to a mass-averaged velocity

  • \(\omega_i\) is the net mass production rate of species \(i\)

Species moles conservation equation#
\[\frac{\partial c X_i}{\partial t} = - \nabla \cdot c X_i \mathbf{u} - \nabla \cdot \mathbf{J}_i + s_i\]

where:

  • \(c\) is the mixture molar density

  • \(X_i\) is the mole fraction of species \(i\)

  • \(\mathbf{u}\) is the molar-averaged velocity of the mixture

  • \(\mathbf{J}_i\) is the molar diffusive flux of species \(i\) relative to a molar-averaged velocity

  • \(s_i\) is the net molar production rate of species \(i\)

Momentum equation#
\[\frac{\partial \rho \mathbf{v}}{\partial t} = - \nabla \cdot \rho \mathbf{v} \mathbf{v} - \nabla \cdot \pmb{\tau} - \nabla \cdot p \mathbf{I} + \rho \sum_{i=1}^{n} Y_i \mathbf{f}_i\]

where:

  • \(\rho\) is the mixture mass density

  • \(\mathbf{v}\) is the mass-averaged velocity of the mixture

  • \(\pmb{\tau}\) is the viscous momentum flux tensor

  • \(p\) is the pressure

  • \(\mathbf{I}\) is the identity matrix

  • \(Y_i\) is the mass fraction of species \(i\)

  • \(\mathbf{f}_i\) is the net acceleration from body forces acting on species \(i\)

  • \(n\) is the number of species (components) in the mixture

Total internal energy equation#
\[\frac{\partial \rho e_0}{\partial t} = - \nabla \cdot \rho e_0 \mathbf{v} - \nabla \cdot \mathbf{q} - \nabla \cdot \pmb{\tau} \cdot \mathbf{v} - \nabla \cdot p \mathbf{v} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{n}_i\]

where:

  • \(\rho\) is the mixture mass density

  • \(e_0\) is the mixture total internal energy

  • \(\mathbf{v}\) is the mass-averaged velocity of the mixture

  • \(\mathbf{q}\) is the heat flux

  • \(\pmb{\tau}\) is the viscous momentum flux tensor

  • \(p\) is the pressure

  • \(\mathbf{f}_i\) is the net acceleration from body forces acting on species \(i\)

  • \(\mathbf{n}_i\) is the total mass flux of species \(i\)

  • \(n\) is the number of species (components) in the mixture

Internal energy equation#
\[\frac{\partial \rho e}{\partial t} = - \nabla \cdot \rho e \mathbf{v} - \nabla \cdot \mathbf{q} - \pmb{\tau} : \nabla \mathbf{v} - p \nabla \cdot \mathbf{v} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{j}_i\]

where:

  • \(\rho\) is the mixture mass density

  • \(e\) is the mixture internal energy

  • \(\mathbf{v}\) is the mass-averaged velocity of the mixture

  • \(\mathbf{q}\) is the heat flux

  • \(\pmb{\tau}\) is the viscous momentum flux tensor

  • \(p\) is the pressure

  • \(\mathbf{f}_i\) is the net acceleration from body forces acting on species \(i\)

  • \(\mathbf{j}_i\) is the mass diffusive flux of species \(i\) relative to a mass-averaged velocity

  • \(n\) is the number of species (components) in the mixture

Enthalpy equation#
\[\frac{\partial \rho h}{\partial t} = - \nabla \cdot \rho h \mathbf{v} - \nabla \cdot \mathbf{q} - \pmb{\tau} : \nabla \mathbf{v} + \frac{Dp}{Dt} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{j}_i\]

where:

  • \(\rho\) is the mixture mass density

  • \(h\) is the mixture enthalpy

  • \(\mathbf{v}\) is the mass-averaged velocity of the mixture

  • \(\mathbf{q}\) is the heat flux

  • \(\pmb{\tau}\) is the viscous momentum flux tensor

  • \(p\) is the pressure

  • \(\mathbf{f}_i\) is the net acceleration from body forces acting on species \(i\)

  • \(\mathbf{j}_i\) is the mass diffusive flux of species \(i\) relative to a mass-averaged velocity

  • \(n\) is the number of species (components) in the mixture

Temperature equation#
\[\rho c_p \frac{DT}{D t} = - \nabla \cdot \mathbf{q} + \alpha T \frac{Dp}{Dt} - \pmb{\tau} : \nabla \mathbf{v} + \sum_{i=1}^{n} \big( h_i (\nabla \cdot \mathbf{j}_i - \omega_i) + \mathbf{f}_i \cdot \mathbf{j}_i \big)\]

where:

  • \(\rho\) is the mixture mass density

  • \(c_p\) is the mixture isobaric specific heat capacity

  • \(T\) is the mixture temperature

  • \(\mathbf{q}\) is the heat flux

  • \(\alpha\) is the coefficient of thermal expansion

  • \(p\) is the pressure

  • \(\pmb{\tau}\) is the viscous momentum flux tensor

  • \(\mathbf{v}\) is the mass-averaged velocity of the mixture

  • \(h_i\) is the enthalpy of species \(i\)

  • \(\mathbf{j}_i\) is the mass diffusive flux of species \(i\) relative to a mass-averaged velocity

  • \(\omega_i\) is the net mass production rate of species \(i\)

  • \(\mathbf{f}_i\) is the net acceleration from body forces acting on species \(i\)

  • \(n\) is the number of species (components) in the mixture

Entropy equation#
\[\frac{\partial \rho s}{\partial t} = - \nabla \cdot \rho s \mathbf{v} - \nabla \Big( \frac{1}{T} \big( \mathbf{q} - \sum_{i=1}^{n} \tilde{\mu}_i \mathbf{j}_i \big) \Big) + \mathbf{q} \cdot \nabla \Big( \frac{1}{T} \Big) - \sum_{i=1}^{n} \mathbf{j}_i \cdot \nabla \Big( \frac{\tilde{\mu}_i}{T} \Big) - \frac{1}{T} \pmb{\tau} : \nabla \mathbf{v} + \frac{1}{T} \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{j}_i - \frac{1}{T} \sum_{i=1}^{n} \tilde{\mu}_i \omega_i\]

where:

  • \(\rho\) is the mixture mass density

  • \(s\) is the mixture entropy

  • \(\mathbf{v}\) is the mass-averaged velocity of the mixture

  • \(T\) is the mixture temperature

  • \(\mathbf{q}\) is the heat flux

  • \(\tilde{\mu}_i\) is the chemical potential of species \(i\)

  • \(\mathbf{j}_i\) is the mass diffusive flux of species \(i\) relative to a mass-averaged velocity

  • \(\pmb{\tau}\) is the viscous momentum flux tensor

  • \(n\) is the number of species (components) in the mixture

Derivations#

This section documents few less commonly presented derivations.


Governing equations#

Summing over all \(n\) species mass conservation equations yields the continuity equation#

If we sum over all \(n\) species mass conservation equations, we should arrive at the continuity equation for the entire mixture. Below, we show that for the Eulerian differential and integral forms of the species mass conservation equation.

The Eulerian differential form of the species mass conservation equation#

We begin with:

\[\frac{\partial \rho Y_i}{\partial t} = - \nabla \cdot \mathbf{n}_i + \omega_i\]

Summing over all \(n\) equations:

\[\sum_{i=1}^{n} \frac{\partial \rho Y_i}{\partial t} = - \sum_{i=1}^{n} \nabla \cdot \mathbf{n}_i + \sum_{i=1}^{n} \omega_i\]

Since \(\sum_{i=1}^{n} \omega_i = 0\):

\[\sum_{i=1}^{n} \frac{\partial \rho Y_i}{\partial t} = - \sum_{i=1}^{n} \nabla \cdot \mathbf{n}_i\]

Which we can also write as:

\[\frac{\partial}{\partial t} \sum_{i=1}^{n} \rho Y_i = - \nabla \cdot \sum_{i=1}^{n} \mathbf{n}_i\]

The density \(\rho\) can now be taken out of the summation. Since \(\sum_{i=1}^{n} Y_i = 1\) and \(\sum_{i=1}^{n} \mathbf{n}_i = \mathbf{n}_t\):

\[\frac{\partial \rho}{\partial t} = - \nabla \cdot \mathbf{n}_t\]

which is the continuity equation in the differential form that we were after. In addition, we may express the total fluxes as: \(\mathbf{n}_t = \rho \mathbf{v}\) which gives:

\[\frac{\partial \rho}{\partial t} = - \nabla \cdot \rho \mathbf{v}\]
The Eulerian integral form of the species mass conservation equation#

We begin with:

\[\int_{V(t)} \omega_i dV = \int_{V(t)} \frac{\partial \rho Y_i}{\partial t} dV + \int_{S(t)} \mathbf{n}_i \cdot \mathbf{a} dS\]

Summing over \(n\) equations:

\[\sum_{i=1}^{n} \int_{V(t)} \omega_i dV = \sum_{i=1}^{n} \int_{V(t)} \frac{\partial \rho Y_i}{\partial t} dV + \sum_{i=1}^{n} \int_{S(t)} \mathbf{n}_i \cdot \mathbf{a} dS\]

which we can write as:

\[\int_{V(t)} \sum_{i=1}^{n} \omega_i dV = \int_{V(t)} \sum_{i=1}^{n} \frac{\partial \rho Y_i}{\partial t} dV + \int_{S(t)} \sum_{i=1}^{n} \mathbf{n}_i \cdot \mathbf{a} dS\]

We can thus use the same reasoning for terms \(\sum_{i=1}^{n} \omega_i\) and \(\sum_{i=1}^{n} \frac{\partial \rho Y_i}{\partial t}\) as previously to get:

\[0 = \int_{V(t)} \frac{\partial \rho}{\partial t} dV + \int_{S(t)} \sum_{i=1}^{n} \mathbf{n}_i \cdot \mathbf{a} dS\]

Finally, we need to think whether \(\sum_{i=1}^{n} (\mathbf{n}_i \cdot \mathbf{a}) = \mathbf{n}_t \cdot \mathbf{a}\)? Once we obtain a discretization of surface \(S(t)\) into multiple \(dS\), each \(dS\) has got a fixed normal vector, \(\mathbf{a}\), associated with it. Since the same Eulerian boundary encapsulates the total flux of species \(i\) in the above equation and the total flux of the mixture, we can conclude that \(\mathbf{a}\) is independent of \(i\) for a given \(S(t)\) and pull it out of the summation. This gives us:

\[\sum_{i=1}^{n} (\mathbf{n}_i \cdot \mathbf{a}) = \mathbf{a} \cdot \sum_{i=1}^{n} \mathbf{n}_i = \mathbf{a} \cdot \mathbf{n}_t\]

which brings us to the continuity equation in the integral form:

\[\int_{V(t)} \frac{\partial \rho}{\partial t} dV = - \int_{S(t)} \mathbf{n}_t \cdot \mathbf{a} dS\]

We may at this point substitute \(\mathbf{n}_t = \rho \mathbf{v}\) and use the divergence theorem to get:

\[\int_{V(t)} \frac{\partial \rho}{\partial t} dV = - \int_{V(t)} \nabla \cdot \rho \mathbf{v} dS\]

Total energy equation from the momentum equation#

We begin with the total internal energy equation:

\[\frac{\partial \rho e_0}{\partial t} = - \nabla \cdot \rho e_0 \mathbf{v} - \nabla \cdot \mathbf{q} - \nabla \cdot \pmb{\tau} \cdot \mathbf{v} - \nabla \cdot p \mathbf{v} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{n}_i\]

and with the momentum equation:

\[\frac{\partial \rho \mathbf{v}}{\partial t} = - \nabla \cdot \rho \mathbf{v} \mathbf{v} - \nabla \cdot \pmb{\tau} - \nabla p + \rho \sum_{i=1}^{n} Y_i \mathbf{f}_i\]

We multiply the momentum equation by \(\mathbf{v}\) to obtain the kinetic energy equation:

\[\boxed{\mathbf{v} \cdot \frac{\partial \rho \mathbf{v}}{\partial t} = - \mathbf{v} \cdot \nabla \cdot \rho \mathbf{v} \mathbf{v} - \mathbf{v} \cdot \nabla \cdot \pmb{\tau} - \mathbf{v} \cdot \nabla p + \mathbf{v} \rho \sum_{i=1}^{n} Y_i \mathbf{f}_i}\]

We note that the total internal energy is equal to the internal energy and the kinetic energy:

\[e_0 = e + \frac{1}{2} \mathbf{v} \cdot \mathbf{v}\]

We first substitute the expression for \(e_0\) in the total internal energy equation:

\[\boxed{\frac{\partial \rho (e + \frac{1}{2} \mathbf{v} \cdot \mathbf{v})}{\partial t} = - \nabla \cdot \rho (e + \frac{1}{2} \mathbf{v} \cdot \mathbf{v}) \mathbf{v} - \nabla \cdot \mathbf{q} - \nabla \cdot \pmb{\tau} \cdot \mathbf{v} - \nabla \cdot p \mathbf{v} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{n}_i}\]

We then subtract the kinetic energy equation from the total internal energy equation, since \(e = e_0 - \frac{1}{2} \mathbf{v} \cdot \mathbf{v}\).

LHS#

LHS after subtraction is:

\[\frac{\partial \rho (e + \frac{1}{2} \mathbf{v} \cdot \mathbf{v})}{\partial t} - \mathbf{v} \cdot \frac{\partial \rho \mathbf{v}}{\partial t}\]

Using the chain rule on all terms:

\[\frac{\partial \rho e}{\partial t} + \frac{1}{2} \rho \mathbf{v} \frac{\partial \mathbf{v}}{\partial t} + \frac{1}{2} \rho \mathbf{v} \frac{\partial \mathbf{v}}{\partial t} + \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \frac{\partial \rho}{\partial t} - \rho \mathbf{v} \cdot \frac{\partial \mathbf{v}}{\partial t} - \mathbf{v} \cdot \mathbf{v} \frac{\partial \rho}{\partial t}\]

This is equal to:

\[\frac{\partial \rho e}{\partial t} - \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \frac{\partial \rho}{\partial t}\]
RHS#

We will split RHS after subtraction into few groups of terms. The first group are terms involving density:

  • \(- \nabla \cdot \rho (e + \frac{1}{2} \mathbf{v} \cdot \mathbf{v}) \mathbf{v} + \mathbf{v} \cdot \nabla \cdot \rho \mathbf{v} \mathbf{v} = - \nabla \cdot \rho e \mathbf{v} - \nabla \cdot \frac{1}{2} \rho \mathbf{v} \cdot \mathbf{v} \mathbf{v} + \mathbf{v} \cdot \nabla \cdot \rho \mathbf{v} \mathbf{v}\)

The second group is the heat flux term:

  • \(- \nabla \cdot \mathbf{q}\)

The third group are terms involving the viscous molecular flux tensor:

  • \(- \nabla \cdot \pmb{\tau} \cdot \mathbf{v} + \mathbf{v} \cdot \nabla \cdot \pmb{\tau} = - \pmb{\tau} : \nabla \mathbf{v} - \mathbf{v} \cdot \nabla \cdot \pmb{\tau} + \mathbf{v} \cdot \nabla \cdot \pmb{\tau} = - \pmb{\tau} : \nabla \mathbf{v}\)

The fourth group are terms involving pressure:

  • \(- \nabla \cdot p \mathbf{v} + \mathbf{v} \cdot \nabla p = - p \nabla \cdot \mathbf{v} - \mathbf{v} \cdot \nabla p + \mathbf{v} \cdot \nabla p = - p \nabla \cdot \mathbf{v}\)

The fifth group are the terms involving body forces:

  • \(\sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{n}_i - \mathbf{v} \rho \sum_{i=1}^{n} Y_i \mathbf{f}_i = \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{n}_i - \sum_{i=1}^{n} \rho Y_i \mathbf{v} \cdot \mathbf{f}_i = \sum_{i=1}^{n} \underbrace{(\mathbf{n}_i - \rho_i \mathbf{v})}_\text{diffusive flux of $i$} \cdot \mathbf{f}_i = \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{j}_i\)

LHS with RHS#

We now put LHS and RHS together to get:

\[\frac{\partial \rho e}{\partial t} \underbrace{- \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \frac{\partial \rho}{\partial t}} = - \nabla \cdot \rho e \mathbf{v} \underbrace{- \nabla \cdot \frac{1}{2} \rho \mathbf{v} \cdot \mathbf{v} \mathbf{v}} + \underbrace{\mathbf{v} \cdot \nabla \cdot \rho \mathbf{v} \mathbf{v}} - \nabla \cdot \mathbf{q} - \pmb{\tau} : \nabla \mathbf{v} - p \nabla \cdot \mathbf{v} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{j}_i\]

We tackle the underbraced terms below:

\[- \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \frac{\partial \rho}{\partial t} = - \nabla \cdot \frac{1}{2} \rho \mathbf{v} \cdot \mathbf{v} \mathbf{v} + \mathbf{v} \cdot \nabla \cdot \rho \mathbf{v} \mathbf{v}\]

Applying the chain rule, we get:

\[- \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \frac{\partial \rho}{\partial t} = - \frac{1}{2} \rho \mathbf{v} \cdot \nabla \cdot \mathbf{v} \mathbf{v} - \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \cdot \nabla \cdot \rho \mathbf{v} + \mathbf{v} \cdot \mathbf{v} \cdot \nabla \cdot \rho \mathbf{v} + \rho \mathbf{v} \cdot \mathbf{v} \cdot \nabla \cdot \mathbf{v}\]

Applying the chain rule one more time on the \(- \frac{1}{2} \rho \mathbf{v} \cdot \nabla \cdot \mathbf{v} \mathbf{v}\) term we get \(- \rho \mathbf{v} \cdot \mathbf{v} \cdot \nabla \cdot \mathbf{v}\). This term cancels out with the \(\rho \mathbf{v} \cdot \mathbf{v} \cdot \nabla \cdot \mathbf{v}\) term. We are thus left with:

\[- \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \frac{\partial \rho}{\partial t} = - \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \cdot \nabla \cdot \rho \mathbf{v} + \mathbf{v} \cdot \mathbf{v} \cdot \nabla \cdot \rho \mathbf{v} = \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \cdot \nabla \cdot \rho \mathbf{v}\]

Collecting all terms on one side and factoring out \(\mathbf{v} \cdot \mathbf{v}\) we get:

\[\frac{1}{2} \mathbf{v} \cdot \mathbf{v} \frac{\partial \rho}{\partial t} + \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \cdot \nabla \cdot \rho \mathbf{v} = \frac{1}{2} \mathbf{v} \cdot \mathbf{v} \underbrace{\Big( \frac{\partial \rho}{\partial t} + \nabla \cdot \rho \mathbf{v} \Big)}_\text{equal to 0 from continuity} = 0\]

Finally, the internal energy equation is:

\[\boxed{\frac{\partial \rho e}{\partial t} = - \nabla \cdot \rho e \mathbf{v} - \nabla \cdot \mathbf{q} - \pmb{\tau} : \nabla \mathbf{v} - p \nabla \cdot \mathbf{v} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{j}_i}\]

Enthalpy equation from the internal energy equation#

Now that we’ve obtained the internal energy equation, we can use relation:

\[h = e + \frac{p}{\rho}\]

to obtain the enthalpy equation. Multiplying the above by \(\rho\):

\[\rho h = \rho e + p\]

and differentiating with respect to time:

\[\frac{\partial \rho h}{\partial t} = \frac{\partial \rho e}{\partial t} + \frac{\partial p}{\partial t}\]

If we insert the two above relationships into the internal energy equation, we get:

\[\frac{\partial \rho h}{\partial t} \underbrace{- \frac{\partial p}{\partial t}} = \underbrace{- \nabla \cdot (\rho h - p) \mathbf{v}} - \nabla \cdot \mathbf{q} - \pmb{\tau} : \nabla \mathbf{v} \underbrace{- p \nabla \cdot \mathbf{v}} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{j}_i\]

With some re-arrangements of the underbraced terms, this is equivalent to:

\[\frac{\partial \rho h}{\partial t} = \underbrace{- \nabla \cdot \rho h \mathbf{v} + \nabla \cdot p \mathbf{v} - p \nabla \cdot \mathbf{v} + \frac{\partial p}{\partial t}} - \nabla \cdot \mathbf{q} - \pmb{\tau} : \nabla \mathbf{v} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{j}_i\]

Applying the chain rule on the \(\nabla \cdot p \mathbf{v}\) term:

\[\frac{\partial \rho h}{\partial t} = - \nabla \cdot \rho h \mathbf{v} + p \nabla \cdot \mathbf{v} + \mathbf{v} \cdot \nabla p - p \nabla \cdot \mathbf{v} + \frac{\partial p}{\partial t} - \nabla \cdot \mathbf{q} - \pmb{\tau} : \nabla \mathbf{v} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{j}_i\]

We can cancel out the \(p \nabla \cdot \mathbf{v}\) term and write \(\frac{\partial p}{\partial t} + \mathbf{v} \cdot \nabla p = \frac{D p}{D t}\).

The enthalpy equation thus becomes:

\[\boxed{\frac{\partial \rho h}{\partial t} = - \nabla \cdot \rho h \mathbf{v} + \frac{D p}{D t} - \nabla \cdot \mathbf{q} - \pmb{\tau} : \nabla \mathbf{v} + \sum_{i=1}^{n} \mathbf{f}_i \cdot \mathbf{j}_i}\]

Non-reacting Stefan tube: species mole fractions#

In this tutorial, we demonstrate how multipy functionalities can be used to compute species mole fractions numerically and analytically in a one-dimensional, non-reacting Stefan tube.

Note

This tutorial has been computed using a Jupyter notebook that can be accessed here.


We start by importing the necessary modules:

import multipy
import numpy as np
from scipy.sparse.linalg import expm
from scipy.integrate import odeint
from scipy.optimize import fmin
from plotly import graph_objects as go

%run -i styles.py
%matplotlib inline

We instantiate objects of the Composition class (that will be necessary for computing composition-related quantities), of the Templates class (that contains useful building blocks that we will use here) and of the Check class (that will be useful for performing some basic checks on the computed quantities):

composition = multipy.Composition()
templates = multipy.Templates()
check = multipy.Check()
filename_prefix = 'non-reacting-stefan-tube'

Problem set-up#

In this tutorial, we will consider the Stefan tube problem - a one-dimensional flow of species inside a long tube. At the bottom of the tube there is a liquid mixture of acetone and methanol, which both evaporate and move upwards the tube. The tube is open to the atmosphere at the top. Air can enter the tube and is thus the third component of the mixture, moving downwards the tube. We will consider a pseudo steady-state condition in which the liquid is fully saturated with air. We also assume that no chemical reactions take place.

The Stefan tube problem is schematically presented below:

_images/stefan-tube-setup.svg

We assume the following species ordering:

  1. Acetone

  1. Methanol

  1. Air

We thus have three species:

n_species = 3

We can create a list of species names:

species_names = ['Acetone', 'Methanol', 'Air']

We set the temperature \(T\) in \([K]\):

T = 328.5

Pressure \(p\) in \([Pa]\):

p = 101325

Mixture molar density \(c\) in \([mole/m^3]\):

mixture_molar_density = composition.mixture_molar_density(T, p)

Species molar masses \(M_i\) in \([kg/mole]\):

species_molar_masses = np.array([[58.08/1000], [32.04/1000], [28.9628/1000]])

Binary diffusivities \(\mathcal{D}\) in \([m^2/s]\):

D12 = 8.48 / 1000**2
D21 = D12
D13 = 13.72 / 1000**2
D31 = D13
D23 = 19.91 / 1000**2
D32 = D23

We can now create a symmetric matrix of size (3, 3) that contains the binary diffusivities:

D_binary = np.array([[0, D12, D13],
                     [D21, 0, D23],
                     [D31, D32, 0]])

If we print the D_binary matrix we should see:

print(D_binary)

Next, we are going to discretize the spatial domain in the \(z\) direction. The length of the domain \(l\) in \([m]\) is:

length = 0.238

We are going to non-dimensionalize the domain and instead of representing it in the \(z\) coordinates, we will represent it in the normalized \(\eta\) coordinates which change from 0 to 1. We are free to select the number of points in the spatial domain (observations) that will be generated.

Discretization of the domain in the \(\eta\) \([-]\) space:

n_points = 1000
eta_coordinates = np.linspace(0,1,n_points)
delta_eta = eta_coordinates[1] - eta_coordinates[0]

We also define the initial condition, \(\mathbf{X}_i(\eta=0)\), which are the known mole fractions of acetone and methanol at the liquid surface, \(\eta=0\):

initial_condition = np.array([0.319, 0.528])

and define the boundary condition, \(\mathbf{X}_i(\eta=1)\), which are the known mole fractions at the tube outlet, \(\eta=1\):

boundary_condition = np.array([0,0,1])

Pause and ponder#

Before we dive into solving the Stefan tube problem, let’s pause for a while and think a little more about what is happening inside the tube. At different positions \(z\), we have different proportions of air to acetone and methanol. This means that the mass density of the mixture, \(\rho\), does not have to be constant throughout the tube. In fact, we know it’s not. Since air has much lower mass density than acetone and methanol, close to the tube outlet the mixture density will be lower than closer to the liquid surface. The mixture density is thus a function of position in the tube: \(\rho = \rho(z)\). This can also be shown by looking at the ideal gas law:

\[\rho(z) = \frac{p M(z)}{R T}\]

since for different compositions at different positions in the tube we will have different values for the molar mass of the mixture, \(M = M(z)\).

The same is not true about the molar density of the mixture, \(c\). Since the mixture is an ideal gas at standard temperature and pressure, we can show that \(c\) is constant:

\[c = \frac{n}{V} = \frac{p}{RT}\]

as long as the temperature and pressure are constant. This might seem weird at first, but yet another way to think about this is that there is the same amount of stuff (moles) throughout the tube, it’s just different stuff at different positions. This is conceptually presented in the figure below:

_images/stefan-tube-moles.svg

Finally, let’s take a closer look at the governing equations that can help us solve the Stefan tube problem.

We can use the species moles conservation equation:

\[\frac{d (c \mathbf{X}_i)}{dt} = - \nabla \cdot \mathbf{N}_i + s_i\]

where \(\mathbf{N}_i\) is the total molar flux of species \(i\) and \(s_i\) is the net molar production rate of species \(i\).

Since we’ve assumed a (pseudo) steady-state, the time derivative seen on the left hand side evaluates to 0. We also assume no reactions, so \(s_i = 0\). The species mole conservation equation thus simplifies to:

\[- \nabla \cdot \mathbf{N}_i = 0\]

We consider only one spatial direction, so the divergence operator can be written for the three species as:

\[\begin{split}\begin{equation} \begin{cases} - \frac{d N_1}{dz} = 0 \\ - \frac{d N_2}{dz} = 0 \\ - \frac{d N_3}{dz} = 0 \end{cases} \end{equation}\end{split}\]

which must mean that all three total molar fluxes are constants: \(N_1 = C_1\), \(N_1 = C_2\), \(N_1 = C_3\), where in general all constants \(C_1\), \(C_2\) and \(C_3\) can be different.

Compute species mole fractions, \(\mathbf{X}_i\)#

We start with computing the mole fraction profiles of all three species inside the Stefan tube.

The system of linear ODEs that we need to solve is:

\[\frac{d \mathbf{X}_i}{d \eta} = \mathbf{\Phi} \mathbf{X}_i + \pmb{\phi}\]

where:

\[\Phi_{ii} = \frac{\mathbf{N}_i}{c_t \mathcal{D}_{i,n}/l} + \sum_{k \neq i}^{n} \frac{\mathbf{N}_k}{c_t \mathcal{D}_{i,k}/l}\]
\[\Phi_{ij} = \mathbf{N}_i \Big( \frac{1}{c_t \mathcal{D}_{i,n}/l} - \frac{1}{c_t \mathcal{D}_{i,j}/l} \Big)\]
\[\phi_i = - \frac{\mathbf{N}_i}{c_t \mathcal{D}_{i,n}/l}\]

We will use the available function Templates.stefan_diffusion that will compute entries in the matrix \(\mathbf{\Phi}\) and in the vector \(\pmb{\phi}\) according to the above formulas. While this function computes these using a generic vector \(\pmb{\alpha}\) and a generic matrix \(\pmb{\beta}\), in our case:

  • \(\pmb{\alpha} = \mathbf{N}_i\)

  • \(\pmb{\beta} = \frac{c}{l}\pmb{\mathcal{D}}\)

Even though an analytic solution exists, let’s begin by solving this system of equations numerically.

Solve numerically using scipy.integrate.odeint#

We will use the Python function scipy.integrate.odeint to compute the mole fractions numerically.

Note that odeint can be used for initial value first order systems of ODEs. The independent variable allowed for odeint can be monotonically increasing or decreasing - in our case this is going to be the normalized spatial coordinates \(\eta\). The initial value corresponds to the first element of the independent variable and in our case this will be the mole fractions at the liquid surface \(x(\eta=0)\).

For our Stefan tube problem, the system of ODEs in the matrix form can be re-written as a system of two ODEs:

\[\begin{split}\begin{equation} \begin{cases} \frac{dX_1}{d \eta} = \Phi_{11} X_1 + \Phi_{12} X_2 + \phi_1 \\ \frac{dX_2}{d \eta} = \Phi_{12} X_1 + \Phi_{22} X_2 + \phi_2 \end{cases} \end{equation}\end{split}\]

subject to initial condition:

\[\begin{split}\begin{equation} \begin{cases} X_1(\eta=0) = 0.319 \\ X_2(\eta=0) = 0.528 \end{cases} \end{equation}\end{split}\]

Note that we solve only two ODEs, since we can get the third mole fraction from the relation:

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

Define the function that will compute the right hand side (RHS) of the system of ODEs:

def RHS_ODE(x, time_vector, N, D_binary, mixture_molar_density, length):

    x1, x2 = x

    (phi, Phi) = templates.stefan_diffusion(N[:,None], mixture_molar_density / length * D_binary)

    gradient_x_acetone = Phi[0,0] * x1 + Phi[0,1] * x2 + phi[0]
    gradient_x_methanol = Phi[1,0] * x1 + Phi[1,1] * x2 + phi[1]

    dxdeta = np.array([gradient_x_acetone, gradient_x_methanol])

    return dxdeta.ravel()

Note that in this problem, we don’t know the total molar fluxes, \(\mathbf{N}_i\). One way to go about this is to take an initial guess and see if after solving the system of ODEs we satisfy the boundary condition. The second way is to use an optimization algorithm that will search for the correct values of the total molar fluxes. We will present the latter approach next. For now, we define the total molar flux vector using the optimized solution that we will soon obtain.

Define the total molar flux \(\mathbf{N}_i\) in \([mole/(m^2s)]\):

N1 = 0.0018175355801046177
N2 = 0.0031885688320410326
N3 = 0

N = np.array([N1, N2, N3])

Solve the two ODEs numerically:

numerical_solution = odeint(RHS_ODE, initial_condition, eta_coordinates, args=(N, D_binary, mixture_molar_density, length))

Get the mole fractions computed numerically:

mole_fractions_numerical = np.vstack((numerical_solution[:,0], numerical_solution[:,1], 1 - numerical_solution[:,0] - numerical_solution[:,1]))

composition.set_species_mole_fractions = mole_fractions_numerical

We check if the values at the boundary match the expected mole fractions:

\[\begin{split}\begin{equation} \begin{cases} X_1(\eta=1) = 0 \\ X_2(\eta=1) = 0 \\ X_3(\eta=1) = 1 \end{cases} \end{equation}\end{split}\]
print('Mole fraction of acetone at the tube outlet:\t' + str(round(mole_fractions_numerical[0,-1], 5)))
print('Mole fraction of methanol at the tube outlet:\t' + str(round(mole_fractions_numerical[1,-1], 5)))
print('Mole fraction of air at the tube outlet:\t' + str(round(mole_fractions_numerical[2,-1],5)))

We also check if the mole fractions of all three species sum to 1.0 and if the mole fractions are bounded between 0 and 1 at any given position inside the tube:

idx = check.sum_of_species_fractions(mole_fractions_numerical, tolerance=0.000001, verbose=True)
idx = check.range_of_species_fractions(mole_fractions_numerical, tolerance=0.000001, verbose=True)

Note that the only reason why we obtained a correct solution is that we already knew the correct total molar fluxes.

Plot the species mole fractions:

composition.plot_species_mole_fractions(species_names=species_names,
                                        colors=colors,
                                        figsize=(6,3),
                                        filename='../images/stefan-tube-mole-fractions-numerical.svg');

Here, we will extend the numerical solution by including an optimization algorithm that will compute the correct total molar fluxes for us. We will use the Python function scipy.optimize.fmin to search for a minimum of an error function that we will define ourselves.

We start by the turning numerical computation of the mole fractions into a function that returns an error. The error will be measured as the difference between the expected boundary condition and the numerically obtained values for mole fractions at the tube outlet:

def error_function(fluxes):

    (N1, N2) = fluxes

    N = np.array([N1, N2, 0])

    numerical_solution = odeint(RHS_ODE, initial_condition, eta_coordinates, args=(N, D_binary, mixture_molar_density, length))

    mole_fractions_numerical = np.vstack((numerical_solution[:,0], numerical_solution[:,1], 1 - numerical_solution[:,0] - numerical_solution[:,1]))

    # Calculate an error between the current mole fractions at the boundary and the boundary condition:
    error = np.linalg.norm(boundary_condition - mole_fractions_numerical[:,-1])

    return error

We can first visualize the error surface as computed from the error_function for a range of values for the first two total molar fluxes, \(N_1\) and \(N_2\):

n_grid = 100
grid_range = np.linspace(0.001,0.004,n_grid)
N1_vector, N2_vector = np.meshgrid(grid_range, grid_range)
N1_vector = N1_vector.ravel()
N2_vector = N2_vector.ravel()
error_vector = np.zeros(len(N1_vector))
for k in range(0, len(N1_vector)):
        error = error_function((N1_vector[k], N2_vector[k]))
        error_vector[k] = error
error_matrix_numerical = np.log(np.reshape(error_vector, (n_grid, n_grid)))

fig = go.Figure(data=[go.Surface(x=grid_range, y=grid_range, z=error_matrix_numerical, colorscale='inferno', showscale=False)])
fig.update_layout(autosize=False,
                width=1000, height=600,
                margin=dict(l=65, r=50, b=65, t=90),
                scene = dict(
                xaxis_title='N1',
                yaxis_title='N2',
                zaxis_title='log(Error)'))
fig.show()
_images/stefan-tube-error-3d-plot.svg

From the plot above we could already take a reasonable guess for the values of \(N_1\) and \(N_2\) where the error function attains the minimum value.

Finally, we compute the optimized solution numerically:

minimum = fmin(func=error_function, x0=(0.001, 0.001), xtol=10**-8, ftol=10**-8)
print('Total molar flux of acetone should be:\t' + str(minimum[0]))
print('Total molar flux of methanol should be:\t' + str(minimum[1]))

Note that the optimized values for \(N_1\) and \(N_2\) are the ones that we used at the start of this numerical solution.

Solve analytically#

The analytic solution is:

\[\mathbf{X} = \exp[\mathbf{\Phi} \cdot \eta] \mathbf{X}_0 + (\exp[\mathbf{\Phi} \cdot \eta] - \mathbf{I}) \mathbf{\Phi}^{-1} \phi\]

where \(\exp[\mathbf{\Phi} \cdot \eta]\) is the exponential of the matrix \(\mathbf{\Phi} \cdot \eta\) defined as:

\[\exp[\mathbf{\Phi} \cdot \eta] = \mathbf{I} + \mathbf{\Phi} \cdot \eta + \frac{1}{2!} (\mathbf{\Phi} \cdot \eta)^2 + \dots = \sum_{k=0}^{\infty} \frac{1}{k!} (\mathbf{\Phi} \cdot \eta)^k\]

which is a generalization of the Maclaurin series for matrices. See reference [2] (Appendix A-B) for the full derivation of this form of the solution.

We are going to use the Python function scipy.sparse.linalg.expm to compute the matrix exponentials seen in the analytic solution.

(phi, Phi) = templates.stefan_diffusion(N[:,None], mixture_molar_density / length * D_binary)

species_mole_fractions = np.zeros((3,n_points))

for i, eta in enumerate(eta_coordinates):

    species_mole_fractions[0:2,i:i+1] = np.dot(expm(Phi*eta), initial_condition[:,None]) + np.dot((expm(Phi*eta) - np.identity(2)), np.dot(np.linalg.inv(Phi), phi))
    species_mole_fractions[2,i] = 1 - np.sum(species_mole_fractions[0:2,i])

Similarly as for the numerical solution, we check if the values at the boundary match the expected mole fractions:

\[\begin{split}X_1(\eta=1) = 0 \\ X_2(\eta=1) = 0 \\ X_3(\eta=1) = 1\end{split}\]
print('Mole fraction of acetone at the tube outlet:\t' + str(round(species_mole_fractions[0,-1], 5)))
print('Mole fraction of methanol at the tube outlet:\t' + str(round(species_mole_fractions[1,-1], 5)))
print('Mole fraction of air at the tube outlet:\t' + str(round(species_mole_fractions[2,-1],5)))

and check if the mole fractions of all three species sum to 1.0 at any given position inside the tube:

idx = check.sum_of_species_fractions(species_mole_fractions, tolerance=0.000001, verbose=True)
idx = check.range_of_species_fractions(species_mole_fractions, tolerance=0.000001, verbose=True)

Set the computed mole fractions as the Composition class atribute:

composition.set_species_mole_fractions = species_mole_fractions

Note, that also for the analytic solution we could have used an optimization algorithm to compute the total molar fluxes, if we didn’t know them already.

Plot the mole fractions:

composition.plot_species_mole_fractions(species_names=species_names,
                                        colors=colors,
                                        figsize=(6,3),
                                        filename='../images/' + filename_prefix + '_species-mole-fractions.svg');
_images/non-reacting-stefan-tube_species-mole-fractions.svg

Save the computed quantities#

Finally, we save the computed quantities for use in other tutorials:

np.savetxt('csv/' + filename_prefix + '_species-mole-fractions.csv', (species_mole_fractions), delimiter=',', fmt='%.16e')

References#

  • [1] J. C. Sutherland - Multicomponent mass transfer course, CHEN-6603, The University of Utah

  • [2] R. Taylor, R. Krishna - Multicomponent mass transfer, Wiley, 1993