tilupy.analytic_sol
Classes
Abstract base class representing simulation results for flow depth and velocity. |
|
Dam-break solution on a dry domain using shallow water theory. |
|
Dam-break solution on a wet domain using shallow water theory. |
|
Dam-break solution on a wet domain using shallow water theory. |
|
Dam-break solution on an inclined dry domain with friction using shallow water theory. |
|
Dam-break solution on a dry domain with friction using shallow water theory. |
|
Dam-break solution on a dry domain with friction using shallow water theory. |
|
Abstract base class representing shape results of a simulated flow. |
|
Shape solution on an inclined dry domain without friction. |
|
Class computing front position of a simulated flow. |
Module Contents
Classes
- class tilupy.analytic_sol.Depth_result(theta: float = None)
Bases:
abc.ABCAbstract base class representing simulation results for flow depth and velocity.
This class defines a common interface for analytical solution that compute flow height h(x,t) and flow velocity u(x,t).
- Parameters:
theta (float, optional) – Angle of the surface, in radian, by default 0.
- _g = 9.81
Gravitational constant.
- Type:
float
- _theta
Angle of the surface, in radian.
- Type:
float
- _x
Spatial coordinates.
- Type:
float or np.ndarray
- _t
Time instant.
- Type:
float or np.ndarray
- _h
Flow height depending on space at a moment.
- Type:
np.ndarray
- _u
Flow velocity depending on space at a moment.
- Type:
np.ndarray
- abstractmethod compute_h(x: float | numpy.ndarray, t: float | numpy.ndarray) None
Virtual function that compute the flow height
_hat given space and time.- Parameters:
x (float or np.ndarray) – Spatial coordinates.
t (float or np.ndarray) – Time instant.
- abstractmethod compute_u(x: float | numpy.ndarray, t: float | numpy.ndarray) None
Virtual function that compute the flow velocity
_uat given space and time.- Parameters:
x (float or np.ndarray) – Spatial coordinates.
t (float or np.ndarray) – Time instant.
- property h
Accessor of h(x,t) solution.
- Returns:
Attribute
_h. If None, no solution computed.- Return type:
numpy.ndarray
- property u
Accessor of u(x,t) solution.
- Returns:
Attribute
_u. If None, no solution computed.- Return type:
numpy.ndarray
- property x
Accessor of the spatial distribution of the computed solution.
- Returns:
Attribute
_x. If None, no solution computed.- Return type:
numpy.ndarray
- property t
Accessor of the time instant of the computed solution.
- Returns:
Attribut
_t. If None, no solution computed.- Return type:
float or numpy.ndarray
- plot(show_h: bool = False, show_u: bool = False, show_surface: bool = False, linestyles: list[str] = None, x_unit: str = 'm', h_unit: str = 'm', u_unit: str = 'm/s', show_plot: bool = True, figsize: tuple = None) matplotlib.axes._axes.Axes
Plot the simulation results.
- Parameters:
show_h (bool, optional) – If True, plot the flow height (
_h) curve.show_u (bool, optional) – If True, plot the flow velocity (
_u) curve.show_surface (bool, optional) – If True, plot the slop of the surface.
linestyles (list[str], optional) – List of linestyle to applie to the graph, must have the same since as the numbre of curve to plot or it will not be taken into account (-1), by default None. If None, copper colormap will be applied.
x_unit (str) – Space unit.
h_unit (str) – Height unit.
u_unit (str) – Velocity unit.
show_plot (bool, optional) – If True, show the resulting plot. By default True.
figsize (tuple, optional) – Size of the wanted plot, by default None.
- Returns:
Resulting plot.
- Return type:
matplotlib.axes._axes.Axes
- Raises:
- class tilupy.analytic_sol.Ritter_dry(h_0: float, x_0: float = 0)
Bases:
Depth_resultDam-break solution on a dry domain using shallow water theory.
This class implements the 1D analytical Ritter’s solution of an ideal dam break on a dry domain. The dam break is instantaneous, over an horizontal and flat surface with no friction. It computes the flow height (took verticaly) and velocity over space and time, based on the equation implemanted in SWASHES, based on Ritter’s equation.
Ritter, A., 1892, Die Fortpflanzung der Wasserwellen, Zeitschrift des Vereines Deutscher Ingenieure, vol. 36(33), p. 947-954.
- Parameters:
x_0 (float, optional) – Initial dam location (position along x-axis), by default 0.
h_0 (float) – Initial water depth to the left of the dam.
- _x0
Initial dam location (position along x-axis).
- Type:
float
- _h0
Initial water depth to the left of the dam.
- Type:
float
- xa(t: float) float
Position of the rarefaction wave front (left-most edge) :
\[x_A(t) = x_0 - t \sqrt{g h_0}\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the front edge of the rarefaction wave.
- Return type:
float
- xb(t: float) float
Position of the contact discontinuity:
\[x_B(t) = x_0 + 2 t \sqrt{g h_0}\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the contact wave (end of rarefaction).
- Return type:
float
- compute_h(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Compute the flow height h(x, t) at given time and positions.
\[\begin{split}h(x, t) = \begin{cases} h_0 & \text{if } x \leq x_A(t), \\\\ \frac{4}{9g} \left( \sqrt{g h_0} - \frac{x - x_0}{2t} \right)^2 & \text{if } x_A(t) < x \leq x_B(t), \\\\ 0 & \text{if } x_B(t) < x, \end{cases}\end{split}\]- Parameters:
x (float or np.ndarray) – Spatial positions.
T (float or nd.ndarray) – Time instant.
Notes
Updates the internal
tilupy.analytic_sol.Depth_result._h,tilupy.analytic_sol.Depth_result._x,tilupy.analytic_sol.Depth_result._tattributes with the computed result.
- compute_u(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Compute the flow velocity u(x, t) at given time and positions.
\[\begin{split}u(x,t) = \begin{cases} 0 & \text{if } x \leq x_A(t), \\\\ \frac{2}{3} \left( \frac{x - x_0}{t} + \sqrt{g h_0} \right) & \text{if } x_A(t) < x \leq x_B(t), \\\\ 0 & \text{if } x_B(t) < x, \end{cases}\end{split}\]- Parameters:
x (float or np.ndarray) – Spatial positions.
T (float or np.ndarray) – Time instant.
Notes
Updates the internal
tilupy.analytic_sol.Depth_result._u,tilupy.analytic_sol.Depth_result._x,tilupy.analytic_sol.Depth_result._tattributes with the computed result.
- class tilupy.analytic_sol.Stoker_SWASHES_wet(x_0: float, h_0: float, h_r: float, h_m: float = None)
Bases:
Depth_resultDam-break solution on a wet domain using shallow water theory.
This class implements the 1D analytical Stoker’s solution of an ideal dam break on a wet domain. The dam break is instantaneous, over an horizontal and flat surface with no friction. It computes the flow height (took verticaly) and velocity over space and time, based on the equation implemanted in SWASHES, based on Stoker’s equation.
Delestre, O., Lucas, C., Ksinant, P.-A., Darboux, F., Laguerre, C., Vo, T.-N.-T., James, F. & Cordier, S., 2013, SWASHES: a compilation of shallow water analytic solutions for hydraulic and environmental studies, International Journal for Numerical Methods in Fluids, v. 72(3), p. 269-300, doi:10.1002/fld.3741.
Stoker, J.J., 1957, Water Waves: The Mathematical Theory with Applications, Pure and Applied Mathematics, vol. 4, Interscience Publishers, New York, USA.
- Parameters:
x_0 (float) – Initial dam location (position along x-axis).
h_0 (float) – Water depth to the left of the dam.
h_r (float) – Water depth to the right of the dam.
h_m (float, optional) – Intermediate height used to compute the critical speed cm. If not provided, it will be computed numerically via the ‘compute_cm()’ method.
- _x0
Initial dam location (position along x-axis).
- Type:
float
- _h0
Water depth to the left of the dam.
- Type:
float
- _hr
Water depth to the right of the dam.
- Type:
float
- _cm
Critical velocity.
- Type:
float
- xa(t: float) float
Position of the rarefaction wave front (left-most edge) :
\[x_A(t) = x_0 - t \sqrt{g h_0}\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the front edge of the rarefaction wave.
- Return type:
float
- xb(t: float) float
Position of the contact discontinuity:
\[x_B(t) = x_0 + t \left( 2 \sqrt{g h_0} - 3 c_m \right)\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the contact wave (end of rarefaction).
- Return type:
float
- xc(t: float) float
Position of the shock wave front (right-most wave):
\[x_C(t) = x_0 + t \cdot \frac{2 c_m^2 \left( \sqrt{g h_0} - c_m \right)}{c_m^2 - g h_r}\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the shock front.
- Return type:
float
- equation_cm(cm) float
Equation of the critical velocity cm:
\[-8.g.hr.cm^{2}.(g.h0 - cm^{2})^{2} + (cm^{2} - g.hr)^{2} . (cm^{2} + g.hr) = 0\]
- compute_cm() None
Solves the non-linear equation to compute the critical velocity
_cm.Uses numerical root-finding to find a valid value of cm that separates the flow regimes. Sets
_cmif a valid solution is found.
- compute_h(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Compute the flow height h(x, t) at given time and positions.
\[\begin{split}h(x, t) = \begin{cases} h_0 & \text{if } x \leq x_A(t), \\\\ \frac{4}{9g} \left( \sqrt{g h_0} - \frac{x - x_0}{2t} \right)^2 & \text{if } x_A(t) < x \leq x_B(t), \\\\ \frac{c_m^2}{g} & \text{if } x_B(t) < x \leq x_C(t), \\\\ h_r & \text{if } x_C(t) < x, \end{cases}\end{split}\]- Parameters:
x (float or np.ndarray) – Spatial positions.
T (float or np.ndarray) – Time instant.
Notes
Updates the internal
tilupy.analytic_sol.Depth_result._h,tilupy.analytic_sol.Depth_result._x,tilupy.analytic_sol.Depth_result._tattributes with the computed result.
- compute_u(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Compute the flow velocity u(x, t) at given time and positions.
\[\begin{split}u(x,t) = \begin{cases} 0 & \text{if } x \leq x_A(t), \\\\ \frac{2}{3} \left( \frac{x - x_0}{t} + \sqrt{g h_0} \right) & \text{if } x_A(t) < x \leq x_B(t), \\\\ 2 \left( \sqrt{g h_0} - c_m \right) & \text{if } x_B(t) < x \leq x_C(t), \\\\ 0 & \text{if } x_C(t) < x, \end{cases}\end{split}\]- Parameters:
x (float or np.ndarray) – Spatial positions.
T (float or np.ndarray) – Time instant.
Notes
Updates the internal
tilupy.analytic_sol.Depth_result._u,tilupy.analytic_sol.Depth_result._x,tilupy.analytic_sol.Depth_result._tattributes with the computed result.
- class tilupy.analytic_sol.Stoker_SARKHOSH_wet(h_0: float, h_r: float)
Bases:
Depth_resultDam-break solution on a wet domain using shallow water theory.
This class implements the 1D analytical Stoker’s solution of an ideal dam break on a wet domain. The dam break is instantaneous, over an horizontal and flat surface with no friction. It computes the flow height (took verticaly) and velocity over space and time, based on the equation implemanted in SWASHES, based on Stoker’s equation.
Sarkhosh, P., 2021, Stoker solution package, version 1.0.0, Zenodo. https://doi.org/10.5281/zenodo.5598374
Stoker, J.J., 1957, Water Waves: The Mathematical Theory with Applications, Pure and Applied Mathematics, vol. 4, Interscience Publishers, New York, USA.
- Parameters:
h_0 (float) – Initial water depth to the left of the dam.
h_r (float) – Initial water depth to the right of the dam.
- _h0
Water depth to the left of the dam.
- Type:
float
- _hr
Water depth to the right of the dam.
- Type:
float
- _cm
Shock front speed.
- Type:
float
- _hm
Height of the shock front.
- Type:
float
- xa(t: float) float
Position of the rarefaction wave front (left-most edge) :
\[x_A(t) = x_0 - t \sqrt{g h_0}\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the front edge of the rarefaction wave.
- Return type:
float
- xb(hm: float, t: float) float
Position of the contact discontinuity:
\[x_B(t) = t \left( 2 \sqrt{g h_0} - 3 \sqrt{g h_m} \right)\]- Parameters:
hm (float) – Height of the shock front.
t (float) – Time instant.
- Returns:
Position of the contact wave (end of rarefaction).
- Return type:
float
- xc(cm: float, t: float) float
Position of the shock wave front (right-most wave):
\[x_C(t) = c_m t\]- Parameters:
cm (float) – Shock front speed.
t (float) – Time instant.
- Returns:
Position of the shock front.
- Return type:
float
- compute_cm() float
Compute the shock front speed using Newton-Raphson’s method to find the solution of:
\[c_m h_r - h_r \left( \sqrt{1 + \frac{8 c_m^2}{g h_r}} - 1 \right) \left( \frac{c_m}{2} - \sqrt{g h_0} + \sqrt{\frac{g h_r}{2} \left( \sqrt{1 + \frac{8 c_m^2}{g h_r}} - 1 \right)} \right) = 0\]- Returns:
Speed of the shock front.
- Return type:
float
- compute_h(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Compute the flow height h(x, t) at given time and positions.
\[\begin{split}h(x, t) = \begin{cases} h_0 & \text{if } x \leq x_A(t), \\\\ \frac{\left( 2 \sqrt{g h_0} - \frac{x}{t} \right)^2}{9 g} & \text{if } x_A(t) < x \leq x_B(t), \\\\ h_m = \frac{1}{2} h_r \left( \sqrt{1 + \frac{8 c_m^2}{g h_r}} - 1 \right) & \text{if } x_B(t) < x \leq x_C(t), \\\\ h_r & \text{if } x_C(t) < x \end{cases}\end{split}\]- Parameters:
x (float or np.ndarray) – Spatial positions.
T (float or np.ndarray) – Time instant.
Notes
Updates the internal
tilupy.analytic_sol.Depth_result._h,tilupy.analytic_sol.Depth_result._x,tilupy.analytic_sol.Depth_result._tattributes with the computed result.
- compute_u(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Compute the flow velocity u(x, t) at given time and positions.
\[\begin{split}u(x,t) = \begin{cases} 0 & \text{if } x \leq x_A(t), \\\\ \frac{2}{3} \left( \frac{x}{t} + \sqrt{g h_0} \right) & \text{if } x_A(t) < x \leq x_B(t), \\\\ 2 \sqrt{g h_0} - 2 \sqrt{g h_m} & \text{if } x_B(t) < x \leq x_C(t), \\\\ 0 & \text{if } x_C(t) < x, \end{cases}\end{split}\]- Parameters:
x (float or np.ndarray) – Spatial positions.
T (float or np.ndarray) – Time instant.
Notes
Updates the internal
tilupy.analytic_sol.Depth_result._u,tilupy.analytic_sol.Depth_result._x,tilupy.analytic_sol.Depth_result._tattributes with the computed result.
- class tilupy.analytic_sol.Mangeney_dry(x_0: float, h_0: float, theta: float, delta: float)
Bases:
Depth_resultDam-break solution on an inclined dry domain with friction using shallow water theory.
This class implements the 1D analytical Stoker’s solution of an ideal dam break on a dry domain. The dam break is instantaneous, over an inclined and flat surface with friction. It computes the flow height (took normal to the surface) and velocity over space and time with an infinitely-long fluid mass on an infinite surface.
Mangeney, A., Heinrich, P., & Roche, R., 2000, Analytical solution for testing debris avalanche numerical models, Pure and Applied Geophysics, vol. 157, p. 1081-1096.
- Parameters:
x_0 (float) – Initial dam location (position along x-axis), by default 0.
h_0 (float) – Initial water depth.
theta (float) – Angle of the surface, in degree.
delta (float) – Dynamic friction angle (20°-40° for debris avalanche), in degree.
- _x0
Initial dam location (position along x-axis).
- Type:
float
- _h0
Initial water depth.
- Type:
float
- _delta
Dynamic friction angle, in radian.
- Type:
float
- _c0
Initial wave propagation speed.
- Type:
float
- _m
Constant horizontal acceleration of the front.
- Type:
float
- xa(t: float) float
Edge of the quiet area:
\[x_A(t) = x_0 + \frac{1}{2}mt^2 - c_0 t\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the edge of the quiet region.
- Return type:
float
- xb(t: float) float
Front of the flow:
\[x_B(t) = x_0 + \frac{1}{2}mt^2 + 2 c_0 t\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the front edge of the fluid.
- Return type:
float
- compute_c0() float
Compute the initial wave propagation speed defined by:
\[c_0 = \sqrt{g h_0 \cos{\theta}}\]- Returns:
Value of the initial wave propagation speed.
- Return type:
float
- compute_m() float
Compute the constant horizontal acceleration of the front defined by:
\[m = g \sin{\theta} - g \cos{\theta} \tan{\delta}\]- Returns:
Value of the constant horizontal acceleration of the front.
- Return type:
float
- compute_h(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Compute the flow height h(x, t) at given time and positions.
\[\begin{split}h(x, t) = \begin{cases} h_0 & \text{if } x \leq x_A(t), \\\\ \frac{1}{9g cos(\theta)} \left(2 c_0 - \frac{x-x_0}{t} + \frac{1}{2} m t \right)^2 & \text{if } x_A(t) < x \leq x_B(t), \\\\ 0 & \text{if } x_B(t) < x, \end{cases}\end{split}\]- Parameters:
x (float or np.ndarray) – Spatial positions.
T (float or np.ndarray) – Time instant.
Notes
Updates the internal
tilupy.analytic_sol.Depth_result._h,tilupy.analytic_sol.Depth_result._x,tilupy.analytic_sol.Depth_result._tattributes with the computed result.
- compute_u(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Compute the flow velocity u(x, t) at given time and positions.
\[\begin{split}u(x,t) = \begin{cases} 0 & \text{if } x \leq x_A(t), \\\\ \frac{2}{3} \left( \frac{x-x_0}{t} + c_0 + mt \right) & \text{if } x_A(t) < x \leq x_B(t), \\\\ 0 & \text{if } x_B(t) < x, \end{cases}\end{split}\]- Parameters:
x (float or np.ndarray) – Spatial positions.
T (float or np.ndarray) – Time instant.
Notes
Updates the internal
tilupy.analytic_sol.Depth_result._u,tilupy.analytic_sol.Depth_result._x,tilupy.analytic_sol.Depth_result._tattributes with the computed result.
- class tilupy.analytic_sol.Dressler_dry(x_0: float, h_0: float, C: float = 40)
Bases:
Depth_resultDam-break solution on a dry domain with friction using shallow water theory.
This class implements the 1D analytical Dressler’s solution of an ideal dam break on a dry domain with friction. The dam break is instantaneous, over an horizontal and flat surface with friction. It computes the flow height (took verticaly) and velocity over space and time, based on the equation implemanted in SWASHES, based on Dressler’s equation.
Dressler, R.F., 1952, Hydraulic resistance effect upon the dam‑break functions, Journal of Research of the National Bureau of Standards, vol. 49(3), p. 217-225.
- Parameters:
x_0 (float) – Initial dam location (position along x-axis).
h_0 (float) – Water depth to the left of the dam.
C (float, optional) – Chézy coefficient, by default 40.
- _x0
Initial dam location (position along x-axis).
- Type:
float
- _h0
Water depth to the left of the dam.
- Type:
float
- _c
Chézy coefficient, by default 40.
- Type:
float, optional
- _xt
Position of the tip area, by default None.
- Type:
float
- _ht
Depth of the tip area, by default None.
- Type:
float, optional
- _ut
Velocity of the tip area, by default None.
- Type:
float, optional
- xa(t: float) float
Position of the rarefaction wave front (left-most edge) :
\[x_A(t) = x_0 - t \sqrt{g h_0}\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the front edge of the rarefaction wave.
- Return type:
float
- xb(t: float) float
Position of the contact discontinuity:
\[x_B(t) = x_0 + 2 t \sqrt{g h_0}\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the contact wave (end of rarefaction).
- Return type:
float
- alpha1(x: float, t: float) float
Correction coefficient for the height:
\[\begin{split}\alpha_1(\xi) = \frac{6}{5(2-\xi)} - \frac{2}{3} + \frac{4 \sqrt{3}}{135} (2-\xi)^{3/2}), \\\\\end{split}\]with \(\xi = \frac{x-x_0}{t\sqrt{g h_0}}\)
- Parameters:
x (float) – Spatial position.
t (float) – Time instant.
- Returns:
Correction coefficient.
- Return type:
float
- alpha2(x: float, t: float) float
Correction coefficient for the velocity:
\[\begin{split}\alpha_2(\xi) = \frac{12}{2-(2-\xi)} - \frac{8}{3} + \frac{8 \sqrt{3}}{189} (2-\xi)^{3/2}) - \frac{108}{7(2 - \xi)}, \\\\\end{split}\]with \(\xi = \frac{x-x_0}{t\sqrt{g h_0}}\)
- Parameters:
x (float) – Spatial position.
t (float) – Time instant.
- Returns:
Correction coefficient.
- Return type:
float
- compute_u(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Call
compute_h().
- compute_h(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Compute the flow height h(x, t) and velocity u(x, t) at given time and positions.
\[\begin{split}h(x, t) = \begin{cases} h_0 & \text{if } x \leq x_A(t), \\\\ \frac{1}{g} \left( \frac{2}{3} \sqrt{g h_0} - \frac{x - x_0}{3t} + \frac{g^{2}}{C^2} \alpha_1 t \right)^2 & \text{if } x_A(t) < x \leq x_t(t), \\\\ \frac{-b-\sqrt{b^2 - 4 a (c-x(t))}}{2 a} & \text{if } x_t(t) < x \leq x_B(t), \\\\ 0 & \text{if } x_B(t) < x, \end{cases}\end{split}\]with \(r = \left. \frac{dx}{dh} \right|_{h = h_t}\), \(c = x_B(t)\), \(a = \frac{r h_t + c - x_t}{h_t^2}\), \(b = r - 2 a h_t\). \(x_t\) and \(h_t\) being the position and the flow depth at the beginning of the tip area.
\[\begin{split}u(x,t) = \begin{cases} 0 & \text{if } x \leq x_A(t), \\\\ u_{co} = \frac{2\sqrt{g h_0}}{3} + \frac{2(x - x_0)}{3t} + \frac{g^2}{C^2} \alpha_2 t & \text{if } x_A(t) < x \leq x_t(t), \\\\ \max_{x \in [x_A(t), x_t(t)]} u_{co}(x, t) & \text{if } x_t(t) < x \leq x_B(t), \\\\ 0 & \text{if } x_B(t) < x, \end{cases}\end{split}\]- Parameters:
x (float | np.ndarray) – Spatial positions.
T (float | np.ndarray) – Time instant.
Notes
Updates the internal
tilupy.analytic_sol.Depth_result._h,tilupy.analytic_sol.Depth_result._u,tilupy.analytic_sol.Depth_result._x,tilupy.analytic_sol.Depth_result._tattributes with the computed result.
- class tilupy.analytic_sol.Chanson_dry(h_0: float, x_0: float, f: float)
Bases:
Depth_resultDam-break solution on a dry domain with friction using shallow water theory.
This class implements the 1D analytical Chanson’s solution of an ideal dam break on a dry domain with friction. The dam break is instantaneous, over an horizontal and flat surface with friction. It computes the flow height (took verticaly) and velocity over space and time, based on the equation implemanted in SWASHES, based on Chanson’s equation.
Chanson, H., 2005, Applications of the Saint-Venant Equations and Method of Characteristics to the Dam Break Wave Problem. https://espace.library.uq.edu.au/view/UQ:9438
- Parameters:
x_0 (float) – Initial dam location (position along x-axis).
h_0 (float) – Water depth to the left of the dam.
f (float) – Darcy friction factor.
- _x0
Initial dam location (position along x-axis).
- Type:
float
- _h0
Water depth to the left of the dam.
- Type:
float
- _f
Darcy friction factor.
- Type:
float, optional
- xa(t: float) float
Position of the rarefaction wave front (left-most edge) :
\[x_A(t) = x_0 - t \sqrt{g h_0}\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the front edge of the rarefaction wave.
- Return type:
float
- xb(t: float) float
Position of the tip of the flow:
\[x_B(t) = x_0 + \left( \frac{3}{2} \frac{U(t)}{\sqrt{g h_0}} - 1 \right) t \sqrt{g h_0}\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the flow tip.
- Return type:
float
- xc(t: float) float
Position of the contact discontinuity:
\[x_C(t) = x_0 + \left( \frac{3}{2} \frac{U(t)}{\sqrt{g h_0}} - 1 \right) t \sqrt{\frac{g}{h_0}} + \frac{4}{f\frac{U(t)^2}{g h_0}} \left( 1 - \frac{U(t)}{2 \sqrt{g h_0}} \right)^4\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the contact wave (wave front).
- Return type:
float
- compute_cf(t: float) float
Compute the celerity of the wave front by resolving:
\[\left( \frac{U}{\sqrt{g h_0}} \right)^3 - 8 \left( 0.75 - \frac{3 f t \sqrt{g}}{8 \sqrt{h_0}} \right) \left( \frac{U}{\sqrt{g h_0}} \right)^2 + 12 \left( \frac{U}{\sqrt{g h_0}} \right) - 8 = 0\]- Parameters:
t (float) – Time instant
- Returns:
Value of the front wave velocity.
- Return type:
float
- compute_h(x: float | numpy.ndarray, T: float | numpy.ndarray) None
Compute the flow height h(x, t) at given time and positions.
\[\begin{split}h(x, t) = \begin{cases} h_0 & \text{if } x \leq x_A(t), \\\\ \frac{4}{9g} \left( \sqrt{g h_0} - \frac{x - x_0}{2t} \right)^2 & \text{if } x_A(t) < x \leq x_B(t), \\\\ \sqrt{\frac{f}{4} \frac{U(t)^2}{g h_0} \frac{x_C(t)-x}{h_0}} & \text{if } x_B(t) < x \leq x_C(t), \\\\ 0 & \text{if } x_C(t) < x \end{cases}\end{split}\]- Parameters:
x (float or np.ndarray) – Spatial positions.
T (float or np.ndarray) – Time instant.
Notes
Updates the internal
tilupy.analytic_sol.Depth_result._h,tilupy.analytic_sol.Depth_result._x,tilupy.analytic_sol.Depth_result._tattributes with the computed result.
- compute_u(x: float | numpy.ndarray, T: float | numpy.ndarray) None
No solution
- class tilupy.analytic_sol.Shape_result(theta: float = 0)
Bases:
abc.ABCAbstract base class representing shape results of a simulated flow.
This class defines a common interface for flow simulation that compute the geometry of the final shape of a flow simulation.
- Parameters:
theta (float, optional) – Angle of the surface, in radian, by default 0.
- _g = 9.81
Gravitational constant.
- Type:
float
- _theta
Angle of the surface, in radian.
- Type:
float
- _x
Spatial coordinates.
- Type:
float or np.ndarray
- _h
Flow height depending on space.
- Type:
np.ndarray
- property h
Accessor of the shape h of the flow.
- Returns:
Attribute
_h. If None, no solution computed.- Return type:
numpy.ndarray
- property x
Accessor of the spatial distribution of the computed solution.
- Returns:
Attribute
_x. If None, no solution computed.- Return type:
numpy.ndarray
- property y
Accessor of the lateral spatial distribution of the computed solution.
- Returns:
Attribute
_y. If None, no solution computed.- Return type:
numpy.ndarray
- class tilupy.analytic_sol.Coussot_shape(rho: float, tau: float, theta: float = 0, h_final: float = 1, H_size: int = 100)
Bases:
Shape_resultShape solution on an inclined dry domain without friction.
This class implements the final shape of a simulated flow. The flow is over an inclined and flat surface without friction with a finite volume of fluid. It computes the spatial coordinates from the flow lenght and height.
Coussot, P., Proust, S., & Ancey, C., 1996, Rheological interpretation of deposits of yield stress fluids, Journal of Non-Newtonian Fluid Mechanics, v. 66(1), p. 55-70, doi:10.1016/0377-0257(96)01474-7.
- Parameters:
rho (float) – Fluid density.
tau (float) – Threshold constraint.
theta (float, optional) – Angle of the surface, in degree, by default 0.
h_final (float, optional) – The final flow depth, by default 1.
H_size (int, optional) – Number of value wanted in the H array, by default 100.
- _rho
Fluid density.
- Type:
float
- _tau
Threshold constraint.
- Type:
float
- _D
Normalized distance of the front from the origin.
- Type:
float or numpy.ndarray
- _H
Normalized fluid depth.
- Type:
float or numpy.ndarray
- _d
Distance of the front from the origin.
- Type:
float or numpy.ndarray
- _h
Fluid depth.
- Type:
float or numpy.ndarray
- _H_size
Number of point in H-axis.
- Type:
int
- h_to_H(h: float) float
Normalize the fluid depth by following:
\[H = \frac{\rho g h \sin(\theta)}{\tau_c}\]If \(\theta = 0\), the expression is:
\[H = \frac{\rho g h}{\tau_c}\]- Parameters:
h (float) – Initial fluid depth.
- Returns:
Normalized fluid depth.
- Return type:
float
- H_to_h(H: float) float
Find the original value of the fluid depth from the normalized one by following:
\[h = \frac{H \tau_c}{\rho g \sin(\theta)}\]If \(\theta = 0\), the expression is:
\[h = \frac{H \tau_c}{\rho g}\]- Parameters:
H (float) – Normalized value of the fluid depth.
- Returns:
True value of the fluid depth.
- Return type:
float
- x_to_X(x: float) float
Normalize the spatial coordinates by following:
\[X = \frac{\rho g x (\sin(\theta))^2}{\tau_c \cos(\theta)}\]If \(\theta = 0\), the expression is:
\[X = \frac{\rho g x}{\tau_c}\]- Parameters:
x (float) – Initial spatial coordinate.
- Returns:
Normalized spatial coordinate.
- Return type:
float
- X_to_x(X: float) float
Find the original value of the spatial coordinates from the normalized one by following:
\[x = \frac{X \tau_c \cos(\theta)}{\rho g (\sin(\theta))^2}\]If \(\theta = 0\), the expression is:
\[x = \frac{X \tau_c}{\rho g}\]- Parameters:
X (float) – Normalized values of the spatial coordinates.
- Returns:
True value of the spatial coordinate.
- Return type:
float
- compute_rheological_test_front_morpho() None
Compute the shape of the frontal lobe from the normalized fluid depth for a rheological test on an inclined surface by following :
\[D = - H - \ln(1 - H)\]If \(\theta = 0\), the expression is:
\[D = \frac{H^2}{2}\]
- compute_rheological_test_lateral_morpho() None
Compute the shape of the lateral lobe from the normalized fluid depth for a rheological test on an inclined surface by following :
\[D = 1 - \sqrt{1 - H^2}\]
- compute_slump_test_hf(h_init: float) float
Compute the final fluid depth for a cylindrical slump test following :
\[\frac{h_f}{h_i} = 1 - \frac{2 \tau_c}{\rho g h_i} \left( 1 - \ln{\frac{2 \tau_c}{\rho g h_i}} \right)\]N. Pashias, D. V. Boger, J. Summers, D. J. Glenister; A fifty cent rheometer for yield stress measurement. J. Rheol. 1 November 1996; 40 (6): 1179-1189. https://doi.org/10.1122/1.550780
- Parameters:
h_init (float) – Initial fluid depth.
- Returns:
Final fluid depth
- Return type:
float
- translate_front(d_final: float) None
Translate the shape of the frontal (or transversal) lobe to the wanted x (or y) coordinate.
- Parameters:
d_final (float) – Final wanted coordinate.
- change_orientation_flow() None
Swap the direction of the result.
Notes
There must not have been any prior translation to use this method.
- interpolate_on_d() None
Interpolate the profile on d-axis.
- class tilupy.analytic_sol.Front_result(h0: float)
Class computing front position of a simulated flow.
This class defines multiple methods for flow simulation that compute the position of the front flow at the specified moment.
- Parameters:
h0 (float) – Initial fluid depth.
- _g = 9.81
Gravitational constant.
- Type:
float
- _h0
Initial fluid depth.
- Type:
float
- _xf
Dictionnary of spatial coordinates of the front flow for each time step (keys).
- Type:
dictionnary
- _labels
Dictionnary of spatial coordinates computation’s method for each time step (keys).
- Type:
dictionnary
- xf_mangeney(t: float, delta: float, theta: float = 0) float
Mangeney’s equation for a dam-break solution over an infinite inclined dry domain with friction and an infinitely-long fluid mass:
\[x_f(t) = \frac{1}{2}mt^2 + 2 c_0 t\]with \(c_0\) the initial wave propagation speed defined by:
\[c_0 = \sqrt{g h_0 \cos{\theta}}\]and \(m\) the constant horizontal acceleration of the front defined by:
\[m = g \sin{\theta} - g \cos{\theta} \tan{\delta}\]Mangeney, A., Heinrich, P., & Roche, R., 2000, Analytical solution for testing debris avalanche numerical models, Pure and Applied Geophysics, vol. 157, p. 1081-1096.
- Parameters:
t (float) – Time instant.
delta (float) – Dynamic friction angle, in degree.
theta (float) – Slope angle, in degree.
- Returns:
Position of the front edge of the fluid.
- Return type:
float
- xf_dressler(t: float) float
Dressler’s equation for a dam-break solution over an infinite inclined dry domain with friction:
\[x_f(t) = 2 t \sqrt{g h_0}\]- Parameters:
t (float) – Time instant.
- Returns:
Position of the front edge of the fluid.
- Return type:
float
- xf_ritter(t: float) float
Ritter’s equation for a dam-break solution over an infinite inclined dry domain without friction:
\[x_f(t) = 2 t \sqrt{g h_0}\]Ritter A. Die Fortpflanzung der Wasserwellen. Zeitschrift des Vereines Deuscher Ingenieure August 1892; 36(33): 947-954.
- Parameters:
t (float) – Time instant.
- Returns:
Position of the front edge of the fluid.
- Return type:
float
- xf_stoker(t: float, hr: float) float
Stoker’s equation for a dam-break solution over an infinite inclined wet domain without friction:
\[x_f(t) =t c_m\]with \(c_m\) the front wave velocity solution of:
\[c_m h_r - h_r \left( \sqrt{1 + \frac{8 c_m^2}{g h_r}} - 1 \right) \left( \frac{c_m}{2} - \sqrt{g h_0} + \sqrt{\frac{g h_r}{2} \left( \sqrt{1 + \frac{8 c_m^2}{g h_r}} - 1 \right)} \right) = 0\]Stoker JJ. Water Waves: The Mathematical Theory with Applications, Pure and Applied Mathematics, Vol. 4. Interscience Publishers: New York, USA, 1957.
Sarkhosh, P. (2021). Stoker solution package (1.0.0). Zenodo. https://doi.org/10.5281/zenodo.5598374
- Parameters:
t (float) – Time instant.
hr (float) – Fluid depth at the right of the dam.
- Returns:
Position of the front edge of the fluid.
- Return type:
float
- xf_chanson(t: float, f: float) float
Chanson’s equation for a dam-break solution over an infinite inclined dry domain with friction:
\[x_f(t) = \left( \frac{3}{2} \frac{U(t)}{\sqrt{g h_0}} - 1 \right) t \sqrt{\frac{g}{h_0}} + \frac{4}{f\frac{U(t)^2}{g h_0}} \left( 1 - \frac{U(t)}{2 \sqrt{g h_0}} \right)^4\]with \(U(t)\) the front wave velocity solution of:
\[\left( \frac{U}{\sqrt{g h_0}} \right)^3 - 8 \left( 0.75 - \frac{3 f t \sqrt{g}}{8 \sqrt{h_0}} \right) \left( \frac{U}{\sqrt{g h_0}} \right)^2 + 12 \left( \frac{U}{\sqrt{g h_0}} \right) - 8 = 0\]Chanson, Hubert. (2005). Analytical Solution of Dam Break Wave with Flow Resistance: Application to Tsunami Surges. 137.
- Parameters:
t (float) – Time instant.
f (float) – Darcy friction coefficient.
- Returns:
Position of the front edge of the fluid.
- Return type:
float
- compute_cf(t: float) float
Compute the celerity of the wave front by resolving:
- Parameters:
t (float) – Time instant
- Returns:
Value of the front wave velocity.
- Return type:
float
- show_fronts_over_methods(x_unit: str = 'm') None
Plot the front distance from the initial position for each method.
- Parameters:
x_unit (str, optional) – X-axis unit, by default “m”
- show_fronts_over_time(x_unit: str = 'm') None
Plot the front distance from the initial position over time for each method.
- Parameters:
x_unit (str, optional) – X-axis unit, by default “m”