josie.euler package

Submodules

josie.euler.eos module

This module contains the different Equation of State (EOS) implementations

class josie.euler.eos.EOS

Bases: object

An Abstract Base Class representing an EOS for an Euler System

abstract p(rho, e)
Return type

Union[ndarray, float]

abstract rho(p, e)
Return type

Union[ndarray, float]

abstract rhoe(rho, p)
Return type

Union[ndarray, float]

abstract sound_velocity(rho, p)
Return type

Union[ndarray, float]

class josie.euler.eos.PerfectGas(gamma=1.4)

Bases: josie.euler.eos.EOS

This class embeds methods to compute states for the Euler problem using an EOS (Equation of State) for perfect gases

p = \rho \mathcal{R} T = \rho \left( \gamma - 1 \right)e

gamma

The adiabatic coefficient

p(rho, e)

This returns the pressure from density and internal energy

Parameters
  • rho (Union[ndarray, float]) – A ArrayAndScalar containing the values of the density on the mesh cells

  • e (Union[ndarray, float]) – A ArrayAndScalar containing the values of the internal energy on the mesh cells

Returns

A :class:`ArrayAndScalar containing the values of the pressure multiplied by the density

Return type

p

rho(p, e)

This returns the sound velocity from density and pressure

Parameters
  • p (Union[ndarray, float]) – A ArrayAndScalar containing the values of the pressure on the mesh cells

  • e (Union[ndarray, float]) – A ArrayAndScalar containing the values of the internal energy on the mesh cells

Returns

A ArrayAndScalar containing the values of the density on the mesh cells

Return type

rho

rhoe(rho, p)

This returns the internal energy multiplied by the density

Parameters
  • rho (Union[ndarray, float]) – A ArrayAndScalar containing the values of the density on the mesh cells

  • p (Union[ndarray, float]) – A ArrayAndScalar containing the values of the pressure on the mesh cells

Returns

A ArrayAndScalar containing the values of the internal energy multiplied by the density

Return type

rhoe

sound_velocity(rho, p)

This returns the sound velocity from density and pressure

Parameters
  • rho (Union[ndarray, float]) – A ArrayAndScalar containing the values of the density on the mesh cells

  • p (Union[ndarray, float]) – A ArrayAndScalar containing the values of the pressure on the mesh cells

Returns

A ArrayAndScalar containing the values of the sound velocity multiplied by the density

Return type

c

class josie.euler.eos.StiffenedGas(gamma=3, p0=1000000000.0)

Bases: josie.euler.eos.PerfectGas

p(rho, e)

This returns the pressure from density and internal energy

Parameters
  • rho (Union[ndarray, float]) – A ArrayAndScalar containing the values of the density on the mesh cells

  • e (Union[ndarray, float]) – A ArrayAndScalar containing the values of the internal energy on the mesh cells

Returns

A :class:`ArrayAndScalar containing the values of the pressure multiplied by the density

Return type

p

rho(p, e)

This returns the sound velocity from density and pressure

Parameters
  • p (Union[ndarray, float]) – A ArrayAndScalar containing the values of the pressure on the mesh cells

  • e (Union[ndarray, float]) – A ArrayAndScalar containing the values of the internal energy on the mesh cells

Returns

A ArrayAndScalar containing the values of the density on the mesh cells

Return type

rho

rhoe(rho, p)

This returns the internal energy multiplied by the density

Parameters
  • rho (Union[ndarray, float]) – A ArrayAndScalar containing the values of the density on the mesh cells

  • p (Union[ndarray, float]) – A ArrayAndScalar containing the values of the pressure on the mesh cells

Returns

A ArrayAndScalar containing the values of the internal energy multiplied by the density

Return type

rhoe

sound_velocity(rho, p)

This returns the sound velocity from density and pressure

Parameters
  • rho (Union[ndarray, float]) – A ArrayAndScalar containing the values of the density on the mesh cells

  • p (Union[ndarray, float]) – A ArrayAndScalar containing the values of the pressure on the mesh cells

Returns

A ArrayAndScalar containing the values of the sound velocity multiplied by the density

Return type

c

josie.euler.exact module

class josie.euler.exact.Exact(eos, Q_L, Q_R)

Bases: object

This class implements the exact solution of the Riemann problem as scheme.

See [Tor09] for a detailed view on compressible schemes and [GR96] for a deep study on hyperbolic systems and the Riemann problem. This is also based on [CG85] and [Kam15]

Parameters
  • eos – The EOS to use

  • Q_L (EulerState) – The left state

  • state (Q_R The right) –

eos

The EOS to use

Q_L

The left state

Q_R The right state
Q_star_L

The left star state

Q_star_R

The right star state

left_wave

The type of the left wave (it’s populated after calling solve()

Type

josie.euler.exact.WaveType

right_wave

The type of the right wave (it’s populated after calling solve()

Type

josie.euler.exact.WaveType

full_state(p_star, U_star)

Compute the full Euler state from the value of p^* and u^*

left_wave: josie.euler.exact.WaveType
rankine_hugoniot(rho, p, rho_k, p_k)

The Rankine-Hugoniot locus connecting the known state before the k shock wave to the “star” state after the shock wave

rarefaction(p, Q_k, wave)

Non linear function for the rarefaction

Parameters
  • p (float) – The after rarefaction fan pressure value

  • rho – An initial estimate for the after rarefaction density (unused)

  • Q_k (EulerState) – The state before the rarefaction fan

  • wave (Wave) – The wave to consider

Return type

float

Returns

  • The value of the non linear function f_k(p^*) which statisfies

  • u^* = U_k \pm f_k(p^*)

rarefaction_ode(p, y, wave)

Solve the rarefaction ODE for a vectorized set of states

Parameters
  • p – is the indipendent variable, i.e. the pressure

  • q – is the state variable. It’s an array of size [num_riemann x 2].

Returns

A tuple containing the derivatives of the state

Return type

derivatives

right_wave: josie.euler.exact.WaveType
sample(x, t, origin=0.5)

Sample the solution at given (x, t)

Returns

The state solution at the requested (x, t)

Return type

state

sample_rarefaction(U_c, V_k, wave)

Return the state within the rarefaction fan

Return type

EulerState

shock(p, Q_k, wave)

This function returns the after-shock velocity

Parameters
  • p (float) – The after shock pressure value

  • rho_guess – An initial estimate for the after shock density

  • Q_k (EulerState) – The state before the shock

  • wave (Wave) – The wave to consider

Returns

Return type

The value of the non linear function after the shock

solve()
class josie.euler.exact.RarefactionState(*args, **kwargs)

Bases: josie.state.State

fields: Type[josie.fields.RarefactionFields]
class josie.euler.exact.Wave(value)

Bases: enum.IntEnum

An enumeration.

LEFT = 1
RIGHT = -1
class josie.euler.exact.WaveType(value)

Bases: enum.Enum

An enumeration.

RAREFACTION = 2
SHOCK = 1

josie.euler.fields module

josie.euler.problem module

class josie.euler.problem.EulerProblem(eos, **kwargs)

Bases: josie.problem.Problem

A class representing an Euler system problem

eos

An instance of EOS, an equation of state for the fluid

F(cells)

This returns the tensor representing the flux for an Euler model

A general problem can be written in a compact way:

\pdv{\vb{q}}{t} + \div{\vb{F\qty(\vb{q})}} + \vb{B}\qty(\vb{q})
\cdot \gradient{\vb{q}} = \vb{s\qty(\vb{q})}

This function needs to return \vb{F}\qty(\vb{q})

Parameters

cells (CellSet) – A MeshCellSet that contains the cell data

Returns

An array of dimension Nx \times Ny \times 4 \times 2, i.e. an array that of each x cell and y cell stores the 4 \times 2 flux tensor

The flux tensor is:

\pdeConvective =
\begin{bmatrix}
    \rho u & \rho v \\
    \rho u^2 + p & \rho uv \\
    \rho vu & \rho v^ 2 + p \\
    (\rho E + p)U & (\rho E + p)V
\end{bmatrix}

Return type

F

josie.euler.solver module

class josie.euler.solver.EulerSolver(mesh, scheme)

Bases: josie.solver.Solver

A solver for the Euler system

t: float

josie.euler.state module

We create one big state that contains the actual conservative variables that are used in the flux together with the “auxiliary” variables that are instead needed, for example, to compute the speed of sound.

\eulerState

  • rho: density \rho

  • rhoU: component along x of the velocity \vb{u},
    multiplied by the density
  • rhoV: component along y of the velocity \vb{u},
    multiplied by the density
  • rhoE: total energy multiplied by the density \rho E

  • rhoe: internal energy multiplied by the density \rho e

  • U: component along x of the velocity u

  • V: component along y of the velocity v

  • p: pressure p

  • c: sound velocity c

class josie.euler.state.EulerConsState(*args, **kwargs)

Bases: josie.state.SubsetState

A State class representing the conservative state variables of the Euler system

fields

alias of josie.fields.ConsFields

full_state_fields

alias of josie.fields.EulerFields

class josie.euler.state.EulerState(*args, **kwargs)

Bases: josie.fluid.state.SingleFluidState

The class representing the state variables of the Euler system

\eulerState

cons_state

alias of josie.euler.state.EulerConsState

fields

alias of josie.fields.EulerFields

Module contents