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