josie.bn package

Submodules

josie.bn.closure module

class josie.bn.closure.Classical

Bases: josie.bn.closure.Closure

This is the classical choice for p_I and \vb{u}_I described in [BN86]

p_I = p_2 \\
\vb{u}_I = \vb{u}_1

pI(state_array)
Parameters

Q – A Nx \times Ny \times 9 array containing the values for all the state variables

Returns

A Nx \times Ny \times 1 array containing the value of the p_I

Return type

pI

uI(state_array)
Parameters

Q – A Nx \times Ny \times 9 array containing the values for all the state variables

Returns

A Nx \times Ny \times 2 array that contains the components of the velocity \vb{u}_I.

Return type

uI

class josie.bn.closure.Closure

Bases: object

A class representing the closure relation for p_I and \vb{u}_I. Use them as mixin with the Equation of State in order to provide full closure for the system

abstract pI(state_array)
Parameters

Q – A Nx \times Ny \times 9 array containing the values for all the state variables

Returns

A Nx \times Ny \times 1 array containing the value of the p_I

Return type

pI

abstract uI(state_array)
Parameters

Q – A Nx \times Ny \times 9 array containing the values for all the state variables

Returns

A Nx \times Ny \times 2 array that contains the components of the velocity \vb{u}_I.

Return type

uI

josie.bn.eos module

class josie.bn.eos.TwoPhaseEOS(phase1, phase2)

Bases: josie.twofluid.state.PhasePair

An Abstract Base Class representing en EOS for a twophase system. In particular two euler.eos.EOS instances for each phase need to be provided.

You can access the EOS for a specified phase using the __getitem__()

josie.bn.problem module

class josie.bn.problem.TwoPhaseProblem(eos, closure)

Bases: josie.problem.Problem

A class representing a two-phase system problem governed by the equations first described in [BN86]

B(cells)

This returns the tensor that pre-multiplies the non-conservative term of the problem.

A general problem can be written in a compact way:

\pdeFull

This method needs in general to return \pdeNonConservativeMultiplier that for this case would be

\pdeNonConservativeMultiplier_r =
\begin{bmatrix}
{u_r}_I & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 \\
0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 \\
-p_I & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 \\
0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 \\
-p_I {u_r}_I & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 \\
0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 \\
p_I & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 \\
0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 \\
p_I {u_r}_I & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 & 0 \\
\end{bmatrix} \qquad r = 1 \dotso N_\text{dim}

But since most of the \pdeNonConservativeMultiplier is zero, since we just have the terms that pre-multiply \gradient{\alpha} we just return B_{p1r} that is:

\tilde{\vb{B}}\qty(\pdeState) =
\begin{bmatrix}
u_I & v_I \\
0 & 0 \\
-p_I & 0 \\
0 & -p_I \\
-p_I u_I & -p_I v_I \\
0 & 0 \\
p_I & 0 \\
0 & p_I \\
p_I u_I & p_I v_I \\
\end{bmatrix}

Parameters
  • cells (Union[CellSet, MeshCellSet]) – A MeshCellSet that contains the cell data

  • eos – An implementation of the equation of state. In particular it needs to implement the Closure trait in order to be able to return pI and uI (in addition to the EOS)

F(cells)

This returns the tensor representing the flux for a two-fluid model as described originally by [BN86]

Parameters

values – A np.ndarray that has dimension Nx \times Ny \times
19 containing the values for all the values variables in all the mesh points

Returns

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

The flux tensor is:

\pdeConvective =
\begin{bmatrix}
    0 & 0 \\
    \alpha_1 \rho u_1 & \alpha_1 \rho v_1 \\
    \alpha_1(\rho_1 u_1^2 + p_1) & \alpha_1 \rho_1 u_1 v_1 \\
    \alpha_1 \rho_1 v_1 u_1 & \alpha_1(\rho v_1^ 2 + p_1) \\
    \alpha_1(\rho_1 E_1 + p_1)u_1 &
        \alpha_1 (\rho_1 E + p)v_1 \\
    \alpha_2 \rho u_2 & \alpha_2 \rho v_2 \\
    \alpha_2(\rho_2 u_2^2 + p_2) & \alpha_2 \rho_2 u_2 v_2 \\
    \alpha_2 \rho_2 v_2 u_2 & \alpha_2(\rho v_2^ 2 + p_2) \\
    \alpha_2(\rho_2 E_2 + p_2)u_2 & \alpha_2 (\rho_2 E + p)v_2
\end{bmatrix}

Return type

F

josie.bn.schemes module

class josie.bn.schemes.BaerScheme(eos, closure)

Bases: josie.scheme.scheme.Scheme

A base class for a twophase scheme

post_step(cells)

During the step we update the conservative values. After the step we update the non-conservative variables. This method updates the values of the non-conservative (auxiliary) variables using the EOS

problem: josie.bn.problem.TwoPhaseProblem
class josie.bn.schemes.Rusanov(eos, closure)

Bases: josie.scheme.convective.ConvectiveScheme, josie.bn.schemes.BaerScheme

CFL(cells, CFL_value)

This method returns the optimal dt value that fulfills the CFL condition for the concrete the given scheme

Parameters
  • cells (MeshCellSet) – A MeshCellSet containing the cell data at the current time step

  • CFL_value – The value of the CFL coefficient to impose

Returns

The Optimal dt fulfilling the CFL condition for the given CFL number

Return type

dt

F(cells, neighs)

This schemes implements the Rusanov scheme for a TwoPhaseProblem. It applies the Rusanov scheme indipendently for each phase (with the \sigma correctly calculated among all the two phases state)

Parameters
  • cells (MeshCellSet) – A MeshCellSet containing the state of the mesh cells

  • neighs (NeighboursCellSet) – A NeighboursCellSet containing data of neighbour cells corresponding to the values

Returns

The value of the numerical convective flux multiplied by the surface value \numConvective

Return type

F

problem: josie.bn.problem.TwoPhaseProblem
class josie.bn.schemes.Upwind(eos, closure)

Bases: josie.scheme.nonconservative.NonConservativeScheme, josie.bn.schemes.BaerScheme

An optimized upwind scheme that reduces the size of the \pdeNonConservativeMultiplier knowing that for [BN86] the only state variable appearing in the non conservative term is \alpha. It concentratres the numerical flux computation into G().

Check also B.

G(cells, neighs)

This is the non-conservative flux implementation of the scheme. See [Tor09] for a great overview on numerical methods for hyperbolic problems.

A general problem can be written in a compact way:

\pdeFull

The non-conservative term is discretized as follows:

\numNonConservativeFull

A concrete instance of this class needs to implement the discretization of the numerical flux on one face of a cell. It needs to implement the term \numNonConservative

Parameters
  • values – The values of the state fields in each cell

  • neighs (NeighboursCellSet) – A NeighboursCellSet containing data of neighbour cells corresponding to the values

Returns

The value of the numerical nonconservative flux multiplied by the surface value \numNonConservative

Return type

G

problem: josie.bn.problem.TwoPhaseProblem

josie.bn.solver module

class josie.bn.solver.BaerSolver(mesh, scheme)

Bases: josie.solver.Solver

A solver for the TwoPhase system

t: float

josie.bn.state module

class josie.bn.state.BaerConsState(*args, **kwargs)

Bases: josie.state.SubsetState

State array for conservative part of the state

fields

alias of josie.fields.BaerConsFields

full_state_fields

alias of josie.fields.BaerFields

class josie.bn.state.BaerPhaseState(*args, **kwargs)

Bases: josie.twofluid.state.PhaseState

State array for one single phase

fields

alias of josie.fields.BaerPhaseFields

full_state_fields

alias of josie.fields.BaerFields

class josie.bn.state.Q(*args, **kwargs)

Bases: josie.twofluid.state.TwoFluidState

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.

The state of system described in [BN86] is actually two Euler states together with the state associated to the volume fraction \alpha

cons_state

alias of josie.bn.state.BaerConsState

fields

alias of josie.fields.BaerFields

get_conservative()

Returns the conservative part of the state

Return type

BaerConsState

get_phase(phase)

Returns the part of the state associated to a specified phase as an instance of PhaseQ

Warning

This does not return the first variable of the state, i.e. \alpha

Parameters

phase (Phases) – A Phases instance identifying the requested phase partition of the state

Returns

A Q instance corresponding to the partition of the system associated to the requested phase

Return type

state

phase_state

alias of josie.bn.state.BaerPhaseState

Module contents