shipmmg package
Submodules
shipmmg.draw_obj module
draw_obj.
Base class for drawing object not only for shipmmg.
- class shipmmg.draw_obj.DrawObj(ax)[source]
Bases:
object
DrawObj class.
General class for drawing object by using matplotlib.animation. Multiple ships can be drawn by using this class.
- draw_obj_with_angle(center_x_list, center_y_list, shape_list, angle_list, obj='ship')[source]
Draw square image with angle.
- Parameters
center_x_list (List[float]) – list of the center x position of the square
center_y_list (List[float]) – list of the center y position of the square
shape_list (List[float]) – list of the square’s shape(length/2, width/2)
angle_list (List[float]) – list of in radians
(str (obj) – optional): object type, ‘ship’ or ‘square’
- Returns
List of Image
- Return type
Image
shipmmg.kt module
kt.
KT (1DOF) simulation code
\(T dr = -r + K \delta\)
- class shipmmg.kt.KTParams(K: float, T: float)[source]
Bases:
object
Dataclass for setting KT parameters of KT simulation.
- K
One of parameters in KT model. [1/s]
- Type
float
- T
One of parameters in KT model. [s]
- Type
float
- K: float
- T: float
- shipmmg.kt.simulate(K: float, T: float, time_list: List[float], δ_list: List[float], r0: float = 0.0, method: str = 'RK45', t_eval=None, events=None, vectorized=False, **options)[source]
KT simulation.
KT simulation by following equation of motion.
\(T dr = -r + K \delta\)
- Parameters
K (float) – parameter K of KT model.
T (float) – parameter T of KT model.
time_list (list[float]) – time list of simulation.
δ_list (list[float]) – rudder angle list of simulation.
r0 (float, optional) – rate of turn [rad/s] in initial condition (time_list[0]). Defaults to 0.0.
method (str, optional) –
Integration method to use in scipy.integrate.solve_ivp():
- ”RK45” (default):
Explicit Runge-Kutta method of order 5(4). The error is controlled assuming accuracy of the fourth-order method, but steps are taken using the fifth-order accurate formula (local extrapolation is done). A quartic interpolation polynomial is used for the dense output. Can be applied in the complex domain.
- ”RK23”:
Explicit Runge-Kutta method of order 3(2). The error is controlled assuming accuracy of the second-order method, but steps are taken using the third-order accurate formula (local extrapolation is done). A cubic Hermite polynomial is used for the dense output. Can be applied in the complex domain.
- ”DOP853”:
Explicit Runge-Kutta method of order 8. Python implementation of the “DOP853” algorithm originally written in Fortran. A 7-th order interpolation polynomial accurate to 7-th order is used for the dense output. Can be applied in the complex domain.
- ”Radau”:
Implicit Runge-Kutta method of the Radau IIA family of order 5. The error is controlled with a third-order accurate embedded formula. A cubic polynomial which satisfies the collocation conditions is used for the dense output.
- ”BDF”:
Implicit multi-step variable-order (1 to 5) method based on a backward differentiation formula for the derivative approximation. A quasi-constant step scheme is used and accuracy is enhanced using the NDF modification. Can be applied in the complex domain.
- ”LSODA”:
Adams/BDF method with automatic stiffness detection and switching. This is a wrapper of the Fortran solver from ODEPACK.
t_eval (array_like or None, optional) – Times at which to store the computed solution, must be sorted and lie within t_span. If None (default), use points selected by the solver.
events (callable, or list of callables, optional) –
Events to track. If None (default), no events will be tracked. Each event occurs at the zeros of a continuous function of time and state. Each function must have the signature event(t, y) and return a float. The solver will find an accurate value of t at which event(t, y(t)) = 0 using a root-finding algorithm. By default, all zeros will be found. The solver looks for a sign change over each step, so if multiple zero crossings occur within one step, events may be missed. Additionally each event function might have the following attributes:
- terminal (bool, optional):
Whether to terminate integration if this event occurs. Implicitly False if not assigned.
- direction (float, optional):
Direction of a zero crossing. If direction is positive, event will only trigger when going from negative to positive, and vice versa if direction is negative. If 0, then either direction will trigger event. Implicitly 0 if not assigned.
You can assign attributes like event.terminal = True to any function in Python.
vectorized (bool, optional) – Whether fun is implemented in a vectorized fashion. Default is False.
options –
Options passed to a chosen solver. All options available for already implemented solvers are listed in scipy.integrate.solve_ivp():
- Returns
- t (ndarray, shape (n_points,)):
Time points.
- y (ndarray, shape (n_points,)):
Values of the solution at t.
- sol (OdeSolution):
Found solution as OdeSolution instance from KT simulation.
- t_events (list of ndarray or None):
Contains for each event type a list of arrays at which an event of that type event was detected. None if events was None.
- y_events (list of ndarray or None):
For each value of t_events, the corresponding value of the solution. None if events was None.
- nfev (int):
Number of evaluations of the right-hand side.
- njev (int):
Number of evaluations of the jacobian.
- nlu (int):
Number of LU decomposition.
- status (int):
- Reason for algorithm termination:
-1: Integration step failed.
0: The solver successfully reached the end of tspan.
1: A termination event occurred.
- message (string):
Human-readable description of the termination reason.
- success (bool):
True if the solver reached the interval end or a termination event occurred (status >= 0).
- Return type
Bunch object with the following fields defined
Examples
>>> K = 0.15 >>> T = 60.0 >>> duration = 300 >>> num_of_sampling = 3000 >>> time_list = np.linspace(0.00, duration, num_of_sampling) >>> δ_list = 35 * np.pi / 180 * np.sin(3.0 * np.pi / Ts * time_list) >>> r0 = 0.0 >>> sol = simulate_kt(K, T, time_list, δ_list, r0) >>> result = sol.sol(time_list)
- shipmmg.kt.simulate_kt(kt_params: shipmmg.kt.KTParams, time_list: List[float], δ_list: List[float], r0: float = 0.0, method: str = 'RK45', t_eval=None, events=None, vectorized=False, **options)[source]
KT simulation.
KT simulation by following equation of motion.
\(T dr = -r + K \delta\)
- Parameters
kt_params (KTParams) – KT parameters.
time_list (list[float]) – time list of simulation.
δ_list (list[float]) – rudder angle list of simulation.
r0 (float, optional) – rate of turn [rad/s] in initial condition (time_list[0]). Defaults to 0.0.
method (str, optional) –
Integration method to use in scipy.integrate.solve_ivp():
- ”RK45” (default):
Explicit Runge-Kutta method of order 5(4). The error is controlled assuming accuracy of the fourth-order method, but steps are taken using the fifth-order accurate formula (local extrapolation is done). A quartic interpolation polynomial is used for the dense output. Can be applied in the complex domain.
- ”RK23”:
Explicit Runge-Kutta method of order 3(2). The error is controlled assuming accuracy of the second-order method, but steps are taken using the third-order accurate formula (local extrapolation is done). A cubic Hermite polynomial is used for the dense output. Can be applied in the complex domain.
- ”DOP853”:
Explicit Runge-Kutta method of order 8. Python implementation of the “DOP853” algorithm originally written in Fortran. A 7-th order interpolation polynomial accurate to 7-th order is used for the dense output. Can be applied in the complex domain.
- ”Radau”:
Implicit Runge-Kutta method of the Radau IIA family of order 5. The error is controlled with a third-order accurate embedded formula. A cubic polynomial which satisfies the collocation conditions is used for the dense output.
- ”BDF”:
Implicit multi-step variable-order (1 to 5) method based on a backward differentiation formula for the derivative approximation. A quasi-constant step scheme is used and accuracy is enhanced using the NDF modification. Can be applied in the complex domain.
- ”LSODA”:
Adams/BDF method with automatic stiffness detection and switching. This is a wrapper of the Fortran solver from ODEPACK.
t_eval (array_like or None, optional) – Times at which to store the computed solution, must be sorted and lie within t_span. If None (default), use points selected by the solver.
events (callable, or list of callables, optional) –
Events to track. If None (default), no events will be tracked. Each event occurs at the zeros of a continuous function of time and state. Each function must have the signature event(t, y) and return a float. The solver will find an accurate value of t at which event(t, y(t)) = 0 using a root-finding algorithm. By default, all zeros will be found. The solver looks for a sign change over each step, so if multiple zero crossings occur within one step, events may be missed. Additionally each event function might have the following attributes:
- terminal (bool, optional):
Whether to terminate integration if this event occurs. Implicitly False if not assigned.
- direction (float, optional):
Direction of a zero crossing. If direction is positive, event will only trigger when going from negative to positive, and vice versa if direction is negative. If 0, then either direction will trigger event. Implicitly 0 if not assigned.
You can assign attributes like event.terminal = True to any function in Python.
vectorized (bool, optional) – Whether fun is implemented in a vectorized fashion. Default is False.
options –
Options passed to a chosen solver. All options available for already implemented solvers are listed in scipy.integrate.solve_ivp():
- Returns
- t (ndarray, shape (n_points,)):
Time points.
- y (ndarray, shape (n_points,)):
Values of the solution at t.
- sol (OdeSolution):
Found solution as OdeSolution instance from KT simulation.
- t_events (list of ndarray or None):
Contains for each event type a list of arrays at which an event of that type event was detected. None if events was None.
- y_events (list of ndarray or None):
For each value of t_events, the corresponding value of the solution. None if events was None.
- nfev (int):
Number of evaluations of the right-hand side.
- njev (int):
Number of evaluations of the jacobian.
- nlu (int):
Number of LU decomposition.
- status (int):
- Reason for algorithm termination:
-1: Integration step failed.
0: The solver successfully reached the end of tspan.
1: A termination event occurred.
- message (string):
Human-readable description of the termination reason.
- success (bool):
True if the solver reached the interval end or a termination event occurred (status >= 0).
- Return type
Bunch object with the following fields defined
Examples
>>> kt_params = KTParams(K=0.15, T=60.0) >>> duration = 300 >>> num_of_sampling = 3000 >>> time_list = np.linspace(0.00, duration, num_of_sampling) >>> δ_list = 35 * np.pi / 180 * np.sin(3.0 * np.pi / Ts * time_list) >>> r0 = 0.0 >>> sol = simulate_kt(kt_params, time_list, δ_list, r0) >>> result = sol.sol(time_list)
- shipmmg.kt.zigzag_test_kt(kt_params: shipmmg.kt.KTParams, target_δ_rad: float, target_ψ_rad_deviation: float, time_list: List[float], δ0: float = 0.0, δ_rad_rate: float = 0.017453292519943295, r0: float = 0.0, ψ0: float = 0.0, method: str = 'RK45', t_eval=None, events=None, vectorized=False, **options)[source]
Zig-zag test simulation.
- Parameters
kt_params (KTParams) – KT parameters.
target_δ_rad (float) – target absolute value of rudder angle.
target_ψ_rad_deviation (float) – target absolute value of psi deviation from ψ0[rad].
time_list (list[float]) – time list of simulation.
δ0 (float) – Initial rudder angle [rad]. Defaults to 0.0.
δ_rad_rate (float) – Initial rudder angle rate [rad/s]. Defaults to 1.0.
r0 (float, optional) – rate of turn [rad/s] in initial condition (time_list[0]). Defaults to 0.0.
ψ0 (float, optional) – Inital azimuth [rad] in initial condition (time_list[0]).. Defaults to 0.0.
method (str, optional) –
Integration method to use in scipy.integrate.solve_ivp():
- ”RK45” (default):
Explicit Runge-Kutta method of order 5(4). The error is controlled assuming accuracy of the fourth-order method, but steps are taken using the fifth-order accurate formula (local extrapolation is done). A quartic interpolation polynomial is used for the dense output. Can be applied in the complex domain.
- ”RK23”:
Explicit Runge-Kutta method of order 3(2). The error is controlled assuming accuracy of the second-order method, but steps are taken using the third-order accurate formula (local extrapolation is done). A cubic Hermite polynomial is used for the dense output. Can be applied in the complex domain.
- ”DOP853”:
Explicit Runge-Kutta method of order 8. Python implementation of the “DOP853” algorithm originally written in Fortran. A 7-th order interpolation polynomial accurate to 7-th order is used for the dense output. Can be applied in the complex domain.
- ”Radau”:
Implicit Runge-Kutta method of the Radau IIA family of order 5. The error is controlled with a third-order accurate embedded formula. A cubic polynomial which satisfies the collocation conditions is used for the dense output.
- ”BDF”:
Implicit multi-step variable-order (1 to 5) method based on a backward differentiation formula for the derivative approximation. A quasi-constant step scheme is used and accuracy is enhanced using the NDF modification. Can be applied in the complex domain.
- ”LSODA”:
Adams/BDF method with automatic stiffness detection and switching. This is a wrapper of the Fortran solver from ODEPACK.
t_eval (array_like or None, optional) – Times at which to store the computed solution, must be sorted and lie within t_span. If None (default), use points selected by the solver.
events (callable, or list of callables, optional) –
Events to track. If None (default), no events will be tracked. Each event occurs at the zeros of a continuous function of time and state. Each function must have the signature event(t, y) and return a float. The solver will find an accurate value of t at which event(t, y(t)) = 0 using a root-finding algorithm. By default, all zeros will be found. The solver looks for a sign change over each step, so if multiple zero crossings occur within one step, events may be missed. Additionally each event function might have the following attributes:
- terminal (bool, optional):
Whether to terminate integration if this event occurs. Implicitly False if not assigned.
- direction (float, optional):
Direction of a zero crossing. If direction is positive, event will only trigger when going from negative to positive, and vice versa if direction is negative. If 0, then either direction will trigger event. Implicitly 0 if not assigned.
You can assign attributes like event.terminal = True to any function in Python.
vectorized (bool, optional) – Whether fun is implemented in a vectorized fashion. Default is False.
options –
Options passed to a chosen solver. All options available for already implemented solvers are listed in scipy.integrate.solve_ivp():
- Returns
list of rudder angle. final_r_list (list[float])) : list of rate of turn.
- Return type
final_δ_list (list[float]))
shipmmg.mmg_3dof module
mmg_3dof.
MMG (3DOF) simulation code
- class shipmmg.mmg_3dof.Mmg3DofBasicParams(L_pp: float, B: float, d: float, x_G: float, D_p: float, m: float, I_zG: float, A_R: float, η: float, m_x: float, m_y: float, J_z: float, f_α: float, ε: float, t_R: float, x_R: float, a_H: float, x_H: float, γ_R_minus: float, γ_R_plus: float, l_R: float, κ: float, t_P: float, w_P0: float, x_P: float)[source]
Bases:
object
Dataclass for setting basic parameters of MMG 3DOF.
- L_pp
Ship length between perpendiculars [m]
- Type
float
- B
Ship breadth [m]
- Type
float
- d
Ship draft [m]
- Type
float
- x_G
Longitudinal coordinate of center of gravity of ship [m]
- Type
float
- D_p
Propeller diameter [m]
- Type
float
- m
Ship mass [kg]
- Type
float
- I_zG
Moment of inertia of ship around center of gravity
- Type
float
- A_R
Profile area of movable part of mariner rudder [m^2]
- Type
float
- η
Ratio of propeller diameter to rudder span (=D_p/HR)
- Type
float
- m_x
Added masses of x axis direction [kg]
- Type
float
- m_y
Added masses of y axis direction [kg]
- Type
float
- J_z
Added moment of inertia
- Type
float
- f_α
Rudder lift gradient coefficient
- Type
float
- ϵ
Ratio of wake fraction at propeller and rudder positions
- Type
float
- t_R
Steering resistance deduction factor
- Type
float
- x_R
Longitudinal coordinate of rudder position.
- Type
float
- a_H
Rudder force increase factor
- Type
float
- x_H
Longitudinal coordinate of acting point of the additional lateral force component induced by steering
- Type
float
- γ_R_minus
Flow straightening coefficient if βR < 0
- Type
float
- γ_R_plus
Flow straightening coefficient if βR > 0
- Type
float
- l_R
Effective longitudinal coordinate of rudder position in formula of βR
- Type
float
- κ
An experimental constant for expressing uR
- Type
float
- t_P
Thrust deduction factor
- Type
float
- w_P0
Wake coefficient at propeller position in straight moving
- Type
float
- x_P
Effective Longitudinal coordinate of propeller position in formula of βP
- Type
float
Note
For more information, please see the following articles.
Yasukawa, H., Yoshimura, Y. (2015) Introduction of MMG standard method for ship maneuvering predictions. J Mar Sci Technol 20, 37-52 https://doi.org/10.1007/s00773-014-0293-y
- A_R: float
- B: float
- D_p: float
- I_zG: float
- J_z: float
- L_pp: float
- a_H: float
- d: float
- f_α: float
- l_R: float
- m: float
- m_x: float
- m_y: float
- t_P: float
- t_R: float
- w_P0: float
- x_G: float
- x_H: float
- x_P: float
- x_R: float
- γ_R_minus: float
- γ_R_plus: float
- ε: float
- η: float
- κ: float
- class shipmmg.mmg_3dof.Mmg3DofManeuveringParams(k_0: float, k_1: float, k_2: float, R_0_dash: float, X_vv_dash: float, X_vr_dash: float, X_rr_dash: float, X_vvvv_dash: float, Y_v_dash: float, Y_r_dash: float, Y_vvv_dash: float, Y_vvr_dash: float, Y_vrr_dash: float, Y_rrr_dash: float, N_v_dash: float, N_r_dash: float, N_vvv_dash: float, N_vvr_dash: float, N_vrr_dash: float, N_rrr_dash: float)[source]
Bases:
object
Dataclass for setting maneuvering parameters of MMG 3ODF.
- k_0
One of manuevering parameters of coefficients representing K_T
- Type
float
- k_1
One of manuevering parameters of coefficients representing K_T
- Type
float
- k_2
One of manuevering parameters of coefficients representing K_T
- Type
float
- R_0_dash
One of manuevering parameters of Ship resistance coefficient in straight moving
- Type
float
- X_vv_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- X_vr_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- X_rr_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- X_vvvv_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- Y_v_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- Y_r_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- Y_vvv_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- Y_vvr_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- Y_vrr_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- Y_rrr_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- N_v_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- N_r_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- N_vvv_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- N_vvr_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- N_vrr_dash
One of manuevering parameters of MMG 3DOF
- Type
float
- N_rrr_dash
One of manuevering parameters of MMG 3DOF
- Type
float
Note
For more information, please see the following articles.
Yasukawa, H., Yoshimura, Y. (2015) Introduction of MMG standard method for ship maneuvering predictions. J Mar Sci Technol 20, 37-52 https://doi.org/10.1007/s00773-014-0293-y
- N_r_dash: float
- N_rrr_dash: float
- N_v_dash: float
- N_vrr_dash: float
- N_vvr_dash: float
- N_vvv_dash: float
- R_0_dash: float
- X_rr_dash: float
- X_vr_dash: float
- X_vv_dash: float
- X_vvvv_dash: float
- Y_r_dash: float
- Y_rrr_dash: float
- Y_v_dash: float
- Y_vrr_dash: float
- Y_vvr_dash: float
- Y_vvv_dash: float
- k_0: float
- k_1: float
- k_2: float
- shipmmg.mmg_3dof.get_sub_values_from_simulation_result(u_list: List[float], v_list: List[float], r_list: List[float], δ_list: List[float], npm_list: List[float], basic_params: shipmmg.mmg_3dof.Mmg3DofBasicParams, maneuvering_params: shipmmg.mmg_3dof.Mmg3DofManeuveringParams, ρ: float = 1025.0, return_all_vals: bool = False)[source]
Get sub values of MMG calculation from simulation result.
- Parameters
u_list (List[float]) – u list of MMG simulation result.
v_list (List[float]) – v list of MMG simulation result.
r_list (List[float]) – r list of MMG simulation result.
δ_list (List[float]) – δ list of MMG simulation result.
npm_list (List[float]) – npm list of MMG simulation result.
basic_params (Mmg3DofBasicParams) – u of MMG simulation result.
maneuvering_params (Mmg3DofManeuveringParams) – u of MMG simulation result.
ρ (float, optional) – seawater density [kg/m^3] Defaults to 1025.0.
return_all_vals (bool, optional) – Whether all sub values are returned or not. Defaults to false.
- Returns
List of X_H X_R_list (List[float]): List of X_R X_P_list (List[float]): List of X_P Y_H_list (List[float]): List of Y_H Y_R_list (List[float]): List of Y_R N_H_list (List[float]): List of N_H N_R_list (List[float]): List of N_R U_list (List[float], optional): List of U if return_all_vals is True β_list (List[float], optional): List of β if return_all_vals is True v_dash_list (List[float], optional): List of v_dash if return_all_vals is True r_dash_list (List[float], optional): List of r_dash if return_all_vals is True w_P_list (List[float], optional): List of w_P if return_all_vals is True J_list (List[float], optional): List of J if return_all_vals is True K_T_list (List[float], optional): List of K_T if return_all_vals is True v_R_list (List[float], optional): List of v_R if return_all_vals is True u_R_list (List[float], optional): List of u_R if return_all_vals is True U_R_list (List[float], optional): List of U_R if return_all_vals is True α_R_list (List[float], optional): List of α_R if return_all_vals is True F_N_list (List[float], optional): List of F_N if return_all_vals is True
- Return type
X_H_list (List[float])
- shipmmg.mmg_3dof.simulate(L_pp: float, B: float, d: float, x_G: float, D_p: float, m: float, I_zG: float, A_R: float, η: float, m_x: float, m_y: float, J_z: float, f_α: float, ε: float, t_R: float, x_R: float, a_H: float, x_H: float, γ_R_minus: float, γ_R_plus: float, l_R: float, κ: float, t_P: float, w_P0: float, x_P: float, k_0: float, k_1: float, k_2: float, R_0_dash: float, X_vv_dash: float, X_vr_dash: float, X_rr_dash: float, X_vvvv_dash: float, Y_v_dash: float, Y_r_dash: float, Y_vvv_dash: float, Y_vvr_dash: float, Y_vrr_dash: float, Y_rrr_dash: float, N_v_dash: float, N_r_dash: float, N_vvv_dash: float, N_vvr_dash: float, N_vrr_dash: float, N_rrr_dash: float, time_list: List[float], δ_list: List[float], npm_list: List[float], u0: float = 0.0, v0: float = 0.0, r0: float = 0.0, x0: float = 0.0, y0: float = 0.0, ψ0: float = 0.0, ρ: float = 1025.0, method: str = 'RK45', t_eval=None, events=None, vectorized=False, **options)[source]
MMG 3DOF simulation.
MMG 3DOF simulation by follwoing equation of motion.
\[ \begin{align}\begin{aligned}m (\dot{u}-vr)&=-m_x\dot{u}+m_yvr+X_H+X_P+X_R\\m (\dot{v}+ur)&=-m_y\dot{v}+m_xur+Y_H+Y_R\\I_{zG}\dot{r}&=-J_Z\dot{r}+N_H+N_R\end{aligned}\end{align} \]- Parameters
L_pp (float) – Ship length between perpendiculars [m]
B (float) – Ship breadth [m]
d (float) – Ship draft [m]
x_G (float) – Longitudinal coordinate of center of gravity of ship
D_p (float) – Propeller diameter [m]
m (float) – Ship mass [kg]
I_zG (float) – Moment of inertia of ship around center of gravity
A_R (float) – Profile area of movable part of mariner rudder [m^2]
η (float) – Ratio of propeller diameter to rudder span (=D_p/HR)
m_x (float) – Added masses of x axis direction [kg]
m_y (float) – Added masses of y axis direction [kg]
J_z (float) – Added moment of inertia
f_α (float) – Rudder lift gradient coefficient
ϵ (float) – Ratio of wake fraction at propeller and rudder positions
t_R (float) – Steering resistance deduction factor
x_R (float) – Longitudinal coordinate of rudder position
a_H (float) – Rudder force increase factor
x_H (float) – Longitudinal coordinate of acting point of the additional lateral force component induced by steering
γ_R_minus (float) – Flow straightening coefficient if βR < 0
γ_R_plus (float) – Flow straightening coefficient if βR > 0
l_R (float) – Effective longitudinal coordinate of rudder position in formula of βR
κ (float) – An experimental constant for expressing uR
t_P (float) – Thrust deduction factor
w_P0 (float) – Wake coefficient at propeller position in straight moving
x_P (float) – Effective Longitudinal coordinate of propeller position in formula of βP
k_0 (float) – One of manuevering parameters of coefficients representing K_T
k_1 (float) – One of manuevering parameters of coefficients representing K_T
k_2 (float) – One of manuevering parameters of coefficients representing K_T
R_0_dash (float) – One of manuevering parameters of MMG 3DOF
X_vv_dash (float) – One of manuevering parameters of MMG 3DOF
X_vr_dash (float) – One of manuevering parameters of MMG 3DOF
X_rr_dash (float) – One of manuevering parameters of MMG 3DOF
X_vvvv_dash (float) – One of manuevering parameters of MMG 3DOF
Y_v_dash (float) – One of manuevering parameters of MMG 3DOF
Y_r_dash (float) – One of manuevering parameters of MMG 3DOF
Y_vvv_dash (float) – One of manuevering parameters of MMG 3DOF
Y_vvr_dash (float) – One of manuevering parameters of MMG 3DOF
Y_vrr_dash (float) – One of manuevering parameters of MMG 3DOF
Y_rrr_dash (float) – One of manuevering parameters of MMG 3DOF
N_v_dash (float) – One of manuevering parameters of MMG 3DOF
N_r_dash (float) – One of manuevering parameters of MMG 3DOF
N_vvv_dash (float) – One of manuevering parameters of MMG 3DOF
N_vvr_dash (float) – One of manuevering parameters of MMG 3DOF
N_vrr_dash (float) – One of manuevering parameters of MMG 3DOF
N_rrr_dash (float) – One of manuevering parameters of MMG 3DOF
time_list (list[float]) – time list of simulation.
δ_list (list[float]) – rudder angle list of simulation.
npm_list (List[float]) – npm list of simulation.
u0 (float, optional) – axial velocity [m/s] in initial condition (time_list[0]). Defaults to 0.0.
v0 (float, optional) – lateral velocity [m/s] in initial condition (time_list[0]). Defaults to 0.0.
r0 (float, optional) – rate of turn [rad/s] in initial condition (time_list[0]). Defaults to 0.0.
x0 (float, optional) – x (surge) position [m] in initial condition (time_list[0]). Defaults to 0.0.
y0 (float, optional) – y (sway) position [m] in initial condition (time_list[0]). Defaults to 0.0.
ψ0 (float, optional) – ψ (yaw) azimuth [rad] in initial condition (time_list[0]). Defaults to 0.0.
ρ (float, optional) – seawater density [kg/m^3] Defaults to 1025.0.
method (str, optional) –
Integration method to use in scipy.integrate.solve_ivp():
- ”RK45” (default):
Explicit Runge-Kutta method of order 5(4). The error is controlled assuming accuracy of the fourth-order method, but steps are taken using the fifth-order accurate formula (local extrapolation is done). A quartic interpolation polynomial is used for the dense output. Can be applied in the complex domain.
- ”RK23”:
Explicit Runge-Kutta method of order 3(2). The error is controlled assuming accuracy of the second-order method, but steps are taken using the third-order accurate formula (local extrapolation is done). A cubic Hermite polynomial is used for the dense output. Can be applied in the complex domain.
- ”DOP853”:
Explicit Runge-Kutta method of order 8. Python implementation of the “DOP853” algorithm originally written in Fortran. A 7-th order interpolation polynomial accurate to 7-th order is used for the dense output. Can be applied in the complex domain.
- ”Radau”:
Implicit Runge-Kutta method of the Radau IIA family of order 5. The error is controlled with a third-order accurate embedded formula. A cubic polynomial which satisfies the collocation conditions is used for the dense output.
- ”BDF”:
Implicit multi-step variable-order (1 to 5) method based on a backward differentiation formula for the derivative approximation. A quasi-constant step scheme is used and accuracy is enhanced using the NDF modification. Can be applied in the complex domain.
- ”LSODA”:
Adams/BDF method with automatic stiffness detection and switching. This is a wrapper of the Fortran solver from ODEPACK.
t_eval (array_like or None, optional) – Times at which to store the computed solution, must be sorted and lie within t_span. If None (default), use points selected by the solver.
events (callable, or list of callables, optional) –
Events to track. If None (default), no events will be tracked. Each event occurs at the zeros of a continuous function of time and state. Each function must have the signature event(t, y) and return a float. The solver will find an accurate value of t at which event(t, y(t)) = 0 using a root-finding algorithm. By default, all zeros will be found. The solver looks for a sign change over each step, so if multiple zero crossings occur within one step, events may be missed. Additionally each event function might have the following attributes:
- terminal (bool, optional):
Whether to terminate integration if this event occurs. Implicitly False if not assigned.
- direction (float, optional):
Direction of a zero crossing. If direction is positive, event will only trigger when going from negative to positive, and vice versa if direction is negative. If 0, then either direction will trigger event. Implicitly 0 if not assigned.
You can assign attributes like event.terminal = True to any function in Python.
vectorized (bool, optional) – Whether fun is implemented in a vectorized fashion. Default is False.
options –
Options passed to a chosen solver. All options available for already implemented solvers are listed in scipy.integrate.solve_ivp():
- Returns
- t (ndarray, shape (n_points,)):
Time points.
- y (ndarray, shape (n_points,)):
Values of the solution at t.
- sol (OdeSolution):
Found solution as OdeSolution instance from MMG 3DOF simulation.
- t_events (list of ndarray or None):
Contains for each event type a list of arrays at which an event of that type event was detected. None if events was None.
- y_events (list of ndarray or None):
For each value of t_events, the corresponding value of the solution. None if events was None.
- nfev (int):
Number of evaluations of the right-hand side.
- njev (int):
Number of evaluations of the jacobian.
- nlu (int):
Number of LU decomposition.
- status (int):
- Reason for algorithm termination:
-1: Integration step failed.
0: The solver successfully reached the end of tspan.
1: A termination event occurred.
- message (string):
Human-readable description of the termination reason.
- success (bool):
True if the solver reached the interval end or a termination event occurred (status >= 0).
- Return type
Bunch object with the following fields defined
Examples
>>> duration = 200 # [s] >>> sampling = 2000 >>> time_list = np.linspace(0.00, duration, sampling) >>> δ_list = np.full(len(time_list), 35.0 * np.pi / 180.0) >>> npm_list = np.full(len(time_list), 17.95) >>> L_pp=7.00 >>> B=1.27 >>> d=0.46 >>> x_G=0.25 >>> D_p=0.216 >>> m=3.27*1.025 >>> I_zG=m*((0.25 * L_pp) ** 2) >>> A_R=0.0539 >>> η=D_p/0.345 >>> m_x=0.022*(0.5 * ρ * (L_pp ** 2) * d) >>> m_y=0.223*(0.5 * ρ * (L_pp ** 2) * d) >>> J_z=0.011*(0.5 * ρ * (L_pp ** 4) * d) >>> f_α=2.747 >>> ϵ=1.09 >>> t_R=0.387 >>> x_R=-0.500*L_pp >>> a_H=0.312 >>> x_H=-0.464*L_pp >>> γ_R=0.395 >>> l_R=-0.710 >>> κ=0.50 >>> t_P=0.220 >>> w_P0=0.40 >>> x_P=-0.650 >>> k_0 = 0.2931 >>> k_1 = -0.2753 >>> k_2 = -0.1385 >>> R_0_dash = 0.022 >>> X_vv_dash = -0.040 >>> X_vr_dash = 0.002 >>> X_rr_dash = 0.011 >>> X_vvvv_dash = 0.771 >>> Y_v_dash = -0.315 >>> Y_r_dash = 0.083 >>> Y_vvv_dash = -1.607 >>> Y_vvr_dash = 0.379 >>> Y_vrr_dash = -0.391 >>> Y_rrr_dash = 0.008 >>> N_v_dash = -0.137 >>> N_r_dash = -0.049 >>> N_vvv_dash = -0.030 >>> N_vvr_dash = -0.294 >>> N_vrr_dash = 0.055 >>> N_rrr_dash = -0.013 >>> sol = simulate_mmg_3dof( >>> L_pp=L_pp, >>> B=B, >>> d=d, >>> x_G=x_G, >>> D_p=D_p, >>> m=m, >>> I_zG=I_zG, >>> A_R=A_R, >>> η=η, >>> m_x=m_x, >>> m_y=m_y, >>> J_z=J_z, >>> f_α=f_α, >>> ϵ=ϵ, >>> t_R=t_R, >>> x_R=x_R, >>> a_H=a_H, >>> x_H=x_H, >>> γ_R=γ_R, >>> l_R=l_R, >>> κ=κ, >>> t_P=t_P, >>> w_P0=w_P0, >>> x_P=x_P, >>> k_0=k_0, >>> k_1=k_1, >>> k_2=k_2, >>> X_0=X_0, >>> X_ββ=X_ββ, >>> X_βγ=X_βγ, >>> X_γγ=X_γγ, >>> X_vvvv_dash=X_vvvv_dash, >>> Y_β=Y_β, >>> Y_γ=Y_γ, >>> Y_βββ=Y_βββ, >>> Y_vvr_dash=Y_vvr_dash, >>> Y_vrr_dash=Y_vrr_dash, >>> Y_rrr_dash=Y_rrr_dash, >>> N_β=N_β, >>> N_γ=N_γ, >>> N_vvv_dash=N_vvv_dash, >>> N_vvr_dash=N_vvr_dash, >>> N_vrr_dash=N_vrr_dash, >>> N_rrr_dash=N_rrr_dash, >>> time_list, >>> δ_rad_list, >>> npm_list, >>> u0=2.29 * 0.512, >>> ) >>> result = sol.sol(time_list)
Note
For more information, please see the following articles.
Yasukawa, H., Yoshimura, Y. (2015) Introduction of MMG standard method for ship maneuvering predictions. J Mar Sci Technol 20, 37-52 https://doi.org/10.1007/s00773-014-0293-y
- shipmmg.mmg_3dof.simulate_mmg_3dof(basic_params: shipmmg.mmg_3dof.Mmg3DofBasicParams, maneuvering_params: shipmmg.mmg_3dof.Mmg3DofManeuveringParams, time_list: List[float], δ_list: List[float], npm_list: List[float], u0: float = 0.0, v0: float = 0.0, r0: float = 0.0, x0: float = 0.0, y0: float = 0.0, ψ0: float = 0.0, ρ: float = 1025.0, method: str = 'RK45', t_eval=None, events=None, vectorized=False, **options)[source]
MMG 3DOF simulation.
MMG 3DOF simulation by follwoing equation of motion.
\[ \begin{align}\begin{aligned}m (\dot{u}-vr)&=-m_x\dot{u}+m_yvr+X_H+X_P+X_R\\m (\dot{v}+ur)&=-m_y\dot{v}+m_xur+Y_H+Y_R\\I_{zG}\dot{r}&=-J_Z\dot{r}+N_H+N_R\end{aligned}\end{align} \]- Parameters
basic_params (Mmg3DofBasicParams) – Basic paramters for MMG 3DOF simulation.
maneuvering_params (Mmg3DofManeuveringParams) – Maneuvering parameters for MMG 3DOF simulation.
time_list (list[float]) – time list of simulation.
δ_list (list[float]) – rudder angle list of simulation.
npm_list (List[float]) – npm list of simulation.
u0 (float, optional) – axial velocity [m/s] in initial condition (time_list[0]). Defaults to 0.0.
v0 (float, optional) – lateral velocity [m/s] in initial condition (time_list[0]). Defaults to 0.0.
r0 (float, optional) – rate of turn [rad/s] in initial condition (time_list[0]). Defaults to 0.0.
x0 (float, optional) – x (surge) position [m] in initial condition (time_list[0]). Defaults to 0.0.
y0 (float, optional) – y (sway) position [m] in initial condition (time_list[0]). Defaults to 0.0.
ψ0 (float, optional) – ψ (yaw) azimuth [rad] in initial condition (time_list[0]). Defaults to 0.0.
ρ (float, optional) – seawater density [kg/m^3] Defaults to 1025.0.
method (str, optional) –
Integration method to use in scipy.integrate.solve_ivp():
- ”RK45” (default):
Explicit Runge-Kutta method of order 5(4). The error is controlled assuming accuracy of the fourth-order method, but steps are taken using the fifth-order accurate formula (local extrapolation is done). A quartic interpolation polynomial is used for the dense output. Can be applied in the complex domain.
- ”RK23”:
Explicit Runge-Kutta method of order 3(2). The error is controlled assuming accuracy of the second-order method, but steps are taken using the third-order accurate formula (local extrapolation is done). A cubic Hermite polynomial is used for the dense output. Can be applied in the complex domain.
- ”DOP853”:
Explicit Runge-Kutta method of order 8. Python implementation of the “DOP853” algorithm originally written in Fortran. A 7-th order interpolation polynomial accurate to 7-th order is used for the dense output. Can be applied in the complex domain.
- ”Radau”:
Implicit Runge-Kutta method of the Radau IIA family of order 5. The error is controlled with a third-order accurate embedded formula. A cubic polynomial which satisfies the collocation conditions is used for the dense output.
- ”BDF”:
Implicit multi-step variable-order (1 to 5) method based on a backward differentiation formula for the derivative approximation. A quasi-constant step scheme is used and accuracy is enhanced using the NDF modification. Can be applied in the complex domain.
- ”LSODA”:
Adams/BDF method with automatic stiffness detection and switching. This is a wrapper of the Fortran solver from ODEPACK.
t_eval (array_like or None, optional) – Times at which to store the computed solution, must be sorted and lie within t_span. If None (default), use points selected by the solver.
events (callable, or list of callables, optional) –
Events to track. If None (default), no events will be tracked. Each event occurs at the zeros of a continuous function of time and state. Each function must have the signature event(t, y) and return a float. The solver will find an accurate value of t at which event(t, y(t)) = 0 using a root-finding algorithm. By default, all zeros will be found. The solver looks for a sign change over each step, so if multiple zero crossings occur within one step, events may be missed. Additionally each event function might have the following attributes:
- terminal (bool, optional):
Whether to terminate integration if this event occurs. Implicitly False if not assigned.
- direction (float, optional):
Direction of a zero crossing. If direction is positive, event will only trigger when going from negative to positive, and vice versa if direction is negative. If 0, then either direction will trigger event. Implicitly 0 if not assigned.
You can assign attributes like event.terminal = True to any function in Python.
vectorized (bool, optional) – Whether fun is implemented in a vectorized fashion. Default is False.
options –
Options passed to a chosen solver. All options available for already implemented solvers are listed in scipy.integrate.solve_ivp():
- Returns
- t (ndarray, shape (n_points,)):
Time points.
- y (ndarray, shape (n_points,)):
Values of the solution at t.
- sol (OdeSolution):
Found solution as OdeSolution instance from MMG 3DOF simulation.
- t_events (list of ndarray or None):
Contains for each event type a list of arrays at which an event of that type event was detected. None if events was None.
- y_events (list of ndarray or None):
For each value of t_events, the corresponding value of the solution. None if events was None.
- nfev (int):
Number of evaluations of the right-hand side.
- njev (int):
Number of evaluations of the jacobian.
- nlu (int):
Number of LU decomposition.
- status (int):
- Reason for algorithm termination:
-1: Integration step failed.
0: The solver successfully reached the end of tspan.
1: A termination event occurred.
- message (string):
Human-readable description of the termination reason.
- success (bool):
True if the solver reached the interval end or a termination event occurred (status >= 0).
- Return type
Bunch object with the following fields defined
Examples
>>> duration = 200 # [s] >>> sampling = 2000 >>> time_list = np.linspace(0.00, duration, sampling) >>> δ_list = np.full(len(time_list), 35.0 * np.pi / 180.0) >>> npm_list = np.full(len(time_list), 20.338) >>> basic_params = Mmg3DofBasicParams( >>> L_pp=7.00, >>> B=1.27, >>> d=0.46, >>> x_G=0.25, >>> D_p=0.216, >>> m=3.27*1.025, >>> I_zG=m*((0.25 * L_pp) ** 2), >>> A_R=0.0539, >>> η=D_p/0.345, >>> m_x=0.022*(0.5 * ρ * (L_pp ** 2) * d), >>> m_y=0.223*(0.5 * ρ * (L_pp ** 2) * d), >>> J_z=0.011*(0.5 * ρ * (L_pp ** 4) * d), >>> f_α=2.747, >>> ϵ=1.09, >>> t_R=0.387, >>> x_R=-0.500*L_pp, >>> a_H=0.312, >>> x_H=-0.464*L_pp, >>> γ_R_minus=0.395, >>> γ_R_plus=0.640, >>> l_R=-0.710, >>> κ=0.50, >>> t_P=0.220, >>> w_P0=0.40, >>> x_P=-0.650, >>> ) >>> maneuvering_params = Mmg3DofManeuveringParams( >>> k_0 = 0.2931, >>> k_1 = -0.2753, >>> k_2 = -0.1385, >>> R_0_dash = 0.022, >>> X_vv_dash = -0.040, >>> X_vr_dash = 0.002, >>> X_rr_dash = 0.011, >>> X_vvvv_dash = 0.771, >>> Y_v_dash = -0.315, >>> Y_r_dash = 0.083, >>> Y_vvv_dash = -1.607, >>> Y_vvr_dash = 0.379, >>> Y_vrr_dash = -0.391, >>> Y_rrr_dash = 0.008, >>> N_v_dash = -0.137, >>> N_r_dash = -0.049, >>> N_vvv_dash = -0.030, >>> N_vvr_dash = -0.294, >>> N_vrr_dash = 0.055, >>> N_rrr_dash = -0.013, >>> ) >>> sol = simulate_mmg_3dof( >>> basic_params, >>> maneuvering_params, >>> time_list, >>> δ_rad_list, >>> npm_list, >>> u0=2.29 * 0.512, >>> ) >>> result = sol.sol(time_list)
Note
For more information, please see the following articles.
Yasukawa, H., Yoshimura, Y. (2015) Introduction of MMG standard method for ship maneuvering predictions. J Mar Sci Technol 20, 37-52 https://doi.org/10.1007/s00773-014-0293-y
- shipmmg.mmg_3dof.zigzag_test_mmg_3dof(basic_params: shipmmg.mmg_3dof.Mmg3DofBasicParams, maneuvering_params: shipmmg.mmg_3dof.Mmg3DofManeuveringParams, target_δ_rad: float, target_ψ_rad_deviation: float, time_list: List[float], npm_list: List[float], δ0: float = 0.0, δ_rad_rate: float = 0.017453292519943295, u0: float = 0.0, v0: float = 0.0, r0: float = 0.0, x0: float = 0.0, y0: float = 0.0, ψ0: float = 0.0, ρ: float = 1025.0, method: str = 'RK45', t_eval=None, events=None, vectorized=False, **options)[source]
Zig-zag test simulation.
- Parameters
basic_params (Mmg3DofBasicParams) – Basic paramters for MMG 3DOF simulation.
maneuvering_params (Mmg3DofManeuveringParams) – Maneuvering parameters for MMG 3DOF simulation.
target_δ_rad (float) – target absolute value of rudder angle.
target_ψ_rad_deviation (float) – target absolute value of psi deviation from ψ0[rad].
time_list (list[float]) – time list of simulation.
npm_list (List[float]) – npm list of simulation.
δ0 (float) – Initial rudder angle [rad]. Defaults to 0.0.
δ_rad_rate (float) – Initial rudder angle rate [rad/s]. Defaults to 1.0.
u0 (float, optional) – axial velocity [m/s] in initial condition (time_list[0]). Defaults to 0.0.
v0 (float, optional) – lateral velocity [m/s] in initial condition (time_list[0]). Defaults to 0.0.
r0 (float, optional) – rate of turn [rad/s] in initial condition (time_list[0]). Defaults to 0.0.
x0 (float, optional) – x (surge) position [m] in initial condition (time_list[0]). Defaults to 0.0.
y0 (float, optional) – y (sway) position [m] in initial condition (time_list[0]). Defaults to 0.0.
ψ0 (float, optional) – Inital azimuth [rad] in initial condition (time_list[0]).. Defaults to 0.0.
ρ (float, optional) – seawater density [kg/m^3] Defaults to 1025.0.
method (str, optional) –
Integration method to use in scipy.integrate.solve_ivp():
- ”RK45” (default):
Explicit Runge-Kutta method of order 5(4). The error is controlled assuming accuracy of the fourth-order method, but steps are taken using the fifth-order accurate formula (local extrapolation is done). A quartic interpolation polynomial is used for the dense output. Can be applied in the complex domain.
- ”RK23”:
Explicit Runge-Kutta method of order 3(2). The error is controlled assuming accuracy of the second-order method, but steps are taken using the third-order accurate formula (local extrapolation is done). A cubic Hermite polynomial is used for the dense output. Can be applied in the complex domain.
- ”DOP853”:
Explicit Runge-Kutta method of order 8. Python implementation of the “DOP853” algorithm originally written in Fortran. A 7-th order interpolation polynomial accurate to 7-th order is used for the dense output. Can be applied in the complex domain.
- ”Radau”:
Implicit Runge-Kutta method of the Radau IIA family of order 5. The error is controlled with a third-order accurate embedded formula. A cubic polynomial which satisfies the collocation conditions is used for the dense output.
- ”BDF”:
Implicit multi-step variable-order (1 to 5) method based on a backward differentiation formula for the derivative approximation. A quasi-constant step scheme is used and accuracy is enhanced using the NDF modification. Can be applied in the complex domain.
- ”LSODA”:
Adams/BDF method with automatic stiffness detection and switching. This is a wrapper of the Fortran solver from ODEPACK.
t_eval (array_like or None, optional) – Times at which to store the computed solution, must be sorted and lie within t_span. If None (default), use points selected by the solver.
events (callable, or list of callables, optional) –
Events to track. If None (default), no events will be tracked. Each event occurs at the zeros of a continuous function of time and state. Each function must have the signature event(t, y) and return a float. The solver will find an accurate value of t at which event(t, y(t)) = 0 using a root-finding algorithm. By default, all zeros will be found. The solver looks for a sign change over each step, so if multiple zero crossings occur within one step, events may be missed. Additionally each event function might have the following attributes:
- terminal (bool, optional):
Whether to terminate integration if this event occurs. Implicitly False if not assigned.
- direction (float, optional):
Direction of a zero crossing. If direction is positive, event will only trigger when going from negative to positive, and vice versa if direction is negative. If 0, then either direction will trigger event. Implicitly 0 if not assigned.
You can assign attributes like event.terminal = True to any function in Python.
vectorized (bool, optional) – Whether fun is implemented in a vectorized fashion. Default is False.
options –
Options passed to a chosen solver. All options available for already implemented solvers are listed in scipy.integrate.solve_ivp():
- Returns
list of rudder angle. final_u_list (list[float])) : list of surge velocity. final_v_list (list[float])) : list of sway velocity. final_r_list (list[float])) : list of rate of turn.
- Return type
final_δ_list (list[float]))
shipmmg.ship_obj_3dof module
ship_obj_3dof.
Ship class for drawing the simulation results and estimating maneuvering parameters.
- class shipmmg.ship_obj_3dof.ShipObj3dof(L: float, B: float, time: List[float] = <factory>, u: List[float] = <factory>, v: List[float] = <factory>, r: List[float] = <factory>, x: List[float] = <factory>, y: List[float] = <factory>, psi: List[float] = <factory>, δ: List[float] = <factory>, npm: List[float] = <factory>)[source]
Bases:
object
Ship 3DOF class just for drawing.
- L
ship length [m]
- Type
float
- B
ship breath [m]
- Type
float
- time
Time list of simulation result.
- Type
list[float]
- u
List of axial velocity [m/s] in simulation result.
- Type
list[float]
- v
List of lateral velocity [m/s] in simulation result.
- Type
list[float]
- r
List of rate of turn [rad/s] in simulation result.
- Type
list[float]
- x
List of position of X axis [m] in simulation result.
- Type
list[float]
- y
List of position of Y axis [m/s] in simulation result.
- Type
list[float]
- psi
List of azimuth [rad] in simulation result.
- Type
list[float]
- δ
rudder angle list of simulation.
- Type
list[float]
- npm
npm list of simulation.
- Type
List[float]
- B: float
- L: float
- draw_chart(x_index: str, y_index: str, xlabel: Optional[str] = None, ylabel: Optional[str] = None, num: Optional[int] = None, figsize: List[float] = [6.4, 4.8], dpi: float = 100.0, facecolor: Optional[str] = None, edgecolor: Optional[str] = None, frameon: bool = True, FigureClass: matplotlib.figure.Figure = <class 'matplotlib.figure.Figure'>, clear: bool = False, tight_layout: bool = False, constrained_layout: bool = False, fmt: Optional[str] = None, save_fig_path: Optional[str] = None, **kwargs) matplotlib.figure.Figure [source]
Draw chart.
- Parameters
x_index (string) – Index value of X axis.
y_index (string) – Index value of Y axis.
xlabel (string, optional) – Label of X axis. Defaults to None.
ylabel (string, optional) – Label of Y axis. Defaults to None.
num (int or str, optional) – A unique identifier for the figure. If a figure with that identifier already exists, this figure is made active and returned. An integer refers to the Figure.number attribute, a string refers to the figure label. If there is no figure with the identifier or num is not given, a new figure is created, made active and returned. If num is an int, it will be used for the Figure.number attribute. Otherwise, an auto-generated integer value is used (starting at 1 and incremented for each new figure). If num is a string, the figure label and the window title is set to this value. Default to None.
figsize ((float, float), optional) – Width, height in inches. Default to [6.4, 4.8]
dpi (float, optional) – The resolution of the figure in dots-per-inch. Default to 100.0.
facecolor (str, optional) – The background color.
edgecolor (str, optional) – The border color.
frameon (bool, optional) – If False, suppress drawing the figure frame. Defaults to True.
FigureClass (subclass of matplotlib.figure.Figure, optional) – Optionally use a custom Figure instance. Defaults to matplotlib.figure.Figure.
clear (bool, optional) – If True and the figure already exists, then it is cleared. Defaults to False.
tight_layout (bool, optional) – If False use subplotpars. If True adjust subplot parameters using tight_layout with default padding. When providing a dict containing the keys pad, w_pad, h_pad, and rect, the default tight_layout paddings will be overridden. Defaults to False.
constrained_layout (bool, optional) – If True use constrained layout to adjust positioning of plot elements. Like tight_layout, but designed to be more flexible. See Constrained Layout Guide for examples. (Note: does not work with add_subplot or subplot2grid.) Defaults to False.
fmt (str, optional) – A format string, e.g. ‘ro’ for red circles. See the Notes section for a full description of the format strings. Format strings are just an abbreviation for quickly setting basic line properties. All of these and more can also be controlled by keyword arguments. This argument cannot be passed as keyword. Defaults to None.
save_fig_path (str, optional) – Path of saving figure. Defaults to None.
**kwargs (matplotlib.lines.Line2D properties, optional) –
kwargs are used to specify properties like a line label (for auto legends), linewidth, antialiasing, marker face color. You can show the detailed information at `matplotlib.lines.Line2D
- Returns
Figure
- Return type
matplotlib.pyplot.Figure
Examples
>>> ship.draw_chart("time", "r", xlabel="time [sec]", >>> ylabel=r"$u$" + " [rad/s]",save_fig_path='test.png')
- draw_gif(dimensionless: bool = False, aspect_equal: bool = True, frate: int = 10, interval: int = 100, num: Optional[int] = None, figsize: List[float] = [6.4, 4.8], dpi: float = 100.0, facecolor: Optional[str] = None, edgecolor: Optional[str] = None, frameon: bool = True, FigureClass: matplotlib.figure.Figure = <class 'matplotlib.figure.Figure'>, clear: bool = False, tight_layout: bool = False, constrained_layout: bool = False, fmt: str = '--k', save_fig_path: Optional[str] = None, **kwargs) matplotlib.figure.Figure [source]
Draw GIF of ship trajectory.
- Parameters
dimensionless (bool, optional) – drawing with dimensionless by using L or not. Defaults to False
aspect_equal (bool, optional) – Set equal of figure aspect or not. Defaults to True.
frate (int, optional) – One of the parameter of frames in matplotlib.FuncAnimation(). frames expresses source of data to pass func and each frame of the animation. frames = int (len(time) / frate) Defaults to 10.
interval (int, optional) – Delay between frames in milliseconds. Defaults to 100.
num (int or str, optional) – A unique identifier for the figure. If a figure with that identifier already exists, this figure is made active and returned. An integer refers to the Figure.number attribute, a string refers to the figure label. If there is no figure with the identifier or num is not given, a new figure is created, made active and returned. If num is an int, it will be used for the Figure.number attribute. Otherwise, an auto-generated integer value is used (starting at 1 and incremented for each new figure). If num is a string, the figure label and the window title is set to this value. Default to None.
figsize ((float, float), optional) – Width, height in inches. Default to [6.4, 4.8]
dpi (float, optional) – The resolution of the figure in dots-per-inch. Default to 100.0.
facecolor (str, optional) – The background color.
edgecolor (str, optional) – The border color.
frameon (bool, optional) – If False, suppress drawing the figure frame. Defaults to True.
FigureClass (subclass of matplotlib.figure.Figure, optional) – Optionally use a custom Figure instance. Defaults to matplotlib.figure.Figure.
clear (bool, optional) – If True and the figure already exists, then it is cleared. Defaults to False.
tight_layout (bool, optional) – If False use subplotpars. If True adjust subplot parameters using tight_layout with default padding. When providing a dict containing the keys pad, w_pad, h_pad, and rect, the default tight_layout paddings will be overridden. Defaults to False.
constrained_layout (bool, optional) – If True use constrained layout to adjust positioning of plot elements. Like tight_layout, but designed to be more flexible. See Constrained Layout Guide for examples. (Note: does not work with add_subplot or subplot2grid.) Defaults to False.
fmt (str, optional) – A format string, e.g. ‘ro’ for red circles. See the Notes section for a full description of the format strings. Format strings are just an abbreviation for quickly setting basic line properties. All of these and more can also be controlled by keyword arguments. This argument cannot be passed as keyword. Defaults to “–k”.
save_fig_path (str, optional) – Path of saving figure. Defaults to None.
**kwargs (matplotlib.lines.Line2D properties, optional) –
kwargs are used to specify properties like a line label (for auto legends), linewidth, antialiasing, marker face color. You can show the detailed information at `matplotlib.lines.Line2D
Examples
>>> ship.draw_gif(save_fig_path='test.gif')
- draw_multi_x_chart(x_index_list: List[str], y_index: str, xlabel: Optional[str] = None, ylabel: Optional[str] = None, num: Optional[int] = None, figsize: List[float] = [6.4, 4.8], dpi: float = 100.0, facecolor: Optional[str] = None, edgecolor: Optional[str] = None, frameon: bool = True, FigureClass: matplotlib.figure.Figure = <class 'matplotlib.figure.Figure'>, clear: bool = False, tight_layout: bool = False, constrained_layout: bool = False, fmt: Optional[str] = None, save_fig_path: Optional[str] = None, **kwargs) matplotlib.figure.Figure [source]
Draw chart of multiple Y variables.
- Parameters
x_index_list (List[string]) – List of index value of X axis.
y_index (string) – Index value of Y axis.
xlabel (string, optional) – Label of X axis. Defaults to None.
ylabel (string, optional) – Label of Y axis. Defaults to None.
num (int or str, optional) – A unique identifier for the figure. If a figure with that identifier already exists, this figure is made active and returned. An integer refers to the Figure.number attribute, a string refers to the figure label. If there is no figure with the identifier or num is not given, a new figure is created, made active and returned. If num is an int, it will be used for the Figure.number attribute. Otherwise, an auto-generated integer value is used (starting at 1 and incremented for each new figure). If num is a string, the figure label and the window title is set to this value. Default to None.
figsize ((float, float), optional) – Width, height in inches. Default to [6.4, 4.8]
dpi (float, optional) – The resolution of the figure in dots-per-inch. Default to 100.0.
facecolor (str, optional) – The background color.
edgecolor (str, optional) – The border color.
frameon (bool, optional) – If False, suppress drawing the figure frame. Defaults to True.
FigureClass (subclass of matplotlib.figure.Figure, optional) – Optionally use a custom Figure instance. Defaults to matplotlib.figure.Figure.
clear (bool, optional) – If True and the figure already exists, then it is cleared. Defaults to False.
tight_layout (bool, optional) – If False use subplotpars. If True adjust subplot parameters using tight_layout with default padding. When providing a dict containing the keys pad, w_pad, h_pad, and rect, the default tight_layout paddings will be overridden. Defaults to False.
constrained_layout (bool, optional) – If True use constrained layout to adjust positioning of plot elements. Like tight_layout, but designed to be more flexible. See Constrained Layout Guide for examples. (Note: does not work with add_subplot or subplot2grid.) Defaults to False.
fmt (str, optional) – A format string, e.g. ‘ro’ for red circles. See the Notes section for a full description of the format strings. Format strings are just an abbreviation for quickly setting basic line properties. All of these and more can also be controlled by keyword arguments. This argument cannot be passed as keyword. Defaults to None.
save_fig_path (str, optional) – Path of saving figure. Defaults to None.
**kwargs (matplotlib.lines.Line2D properties, optional) –
kwargs are used to specify properties like a line label (for auto legends), linewidth, antialiasing, marker face color. You can show the detailed information at `matplotlib.lines.Line2D
- Returns
Figure
- Return type
matplotlib.pyplot.Figure
Examples
>>> ship.draw_chart("time", "r", xlabel="time [sec]", >>> ylabel=r"$u$" + " [rad/s]",save_fig_path='test.png')
- draw_multi_y_chart(x_index: str, y_index_list: List[str], xlabel: Optional[str] = None, ylabel: Optional[str] = None, num: Optional[int] = None, figsize: List[float] = [6.4, 4.8], dpi: float = 100.0, facecolor: Optional[str] = None, edgecolor: Optional[str] = None, frameon: bool = True, FigureClass: matplotlib.figure.Figure = <class 'matplotlib.figure.Figure'>, clear: bool = False, tight_layout: bool = False, constrained_layout: bool = False, fmt: Optional[str] = None, save_fig_path: Optional[str] = None, **kwargs) matplotlib.figure.Figure [source]
Draw chart of multiple Y variables.
- Parameters
x_index (string) – Index value of X axis.
y_index_list (List[string]) – List of index value of Y axis.
xlabel (string, optional) – Label of X axis. Defaults to None.
ylabel (string, optional) – Label of Y axis. Defaults to None.
num (int or str, optional) – A unique identifier for the figure. If a figure with that identifier already exists, this figure is made active and returned. An integer refers to the Figure.number attribute, a string refers to the figure label. If there is no figure with the identifier or num is not given, a new figure is created, made active and returned. If num is an int, it will be used for the Figure.number attribute. Otherwise, an auto-generated integer value is used (starting at 1 and incremented for each new figure). If num is a string, the figure label and the window title is set to this value. Default to None.
figsize ((float, float), optional) – Width, height in inches. Default to [6.4, 4.8]
dpi (float, optional) – The resolution of the figure in dots-per-inch. Default to 100.0.
facecolor (str, optional) – The background color.
edgecolor (str, optional) – The border color.
frameon (bool, optional) – If False, suppress drawing the figure frame. Defaults to True.
FigureClass (subclass of matplotlib.figure.Figure, optional) – Optionally use a custom Figure instance. Defaults to matplotlib.figure.Figure.
clear (bool, optional) – If True and the figure already exists, then it is cleared. Defaults to False.
tight_layout (bool, optional) – If False use subplotpars. If True adjust subplot parameters using tight_layout with default padding. When providing a dict containing the keys pad, w_pad, h_pad, and rect, the default tight_layout paddings will be overridden. Defaults to False.
constrained_layout (bool, optional) – If True use constrained layout to adjust positioning of plot elements. Like tight_layout, but designed to be more flexible. See Constrained Layout Guide for examples. (Note: does not work with add_subplot or subplot2grid.) Defaults to False.
fmt (str, optional) – A format string, e.g. ‘ro’ for red circles. See the Notes section for a full description of the format strings. Format strings are just an abbreviation for quickly setting basic line properties. All of these and more can also be controlled by keyword arguments. This argument cannot be passed as keyword. Defaults to None.
save_fig_path (str, optional) – Path of saving figure. Defaults to None.
**kwargs (matplotlib.lines.Line2D properties, optional) –
kwargs are used to specify properties like a line label (for auto legends), linewidth, antialiasing, marker face color. You can show the detailed information at `matplotlib.lines.Line2D
- Returns
Figure
- Return type
matplotlib.pyplot.Figure
Examples
>>> ship.draw_chart("time", "r", xlabel="time [sec]", >>> ylabel=r"$u$" + " [rad/s]",save_fig_path='test.png')
- draw_xy_trajectory(dimensionless: bool = False, aspect_equal: bool = True, num: Optional[int] = None, figsize: List[float] = [6.4, 4.8], dpi: float = 100.0, fmt: Optional[str] = None, facecolor: Optional[str] = None, edgecolor: Optional[str] = None, frameon: bool = True, FigureClass: matplotlib.figure.Figure = <class 'matplotlib.figure.Figure'>, clear: bool = False, tight_layout: bool = False, constrained_layout: bool = False, save_fig_path: Optional[str] = None, **kwargs) matplotlib.figure.Figure [source]
Draw trajectry(x,y).
- Parameters
dimensionless (bool, optional) – drawing with dimensionless by using L or not. Defaults to False
aspect_equal (bool, optional) – Set equal of figure aspect or not. Defaults to True.
num (int or str, optional) – A unique identifier for the figure. If a figure with that identifier already exists, this figure is made active and returned. An integer refers to the Figure.number attribute, a string refers to the figure label. If there is no figure with the identifier or num is not given, a new figure is created, made active and returned. If num is an int, it will be used for the Figure.number attribute. Otherwise, an auto-generated integer value is used (starting at 1 and incremented for each new figure). If num is a string, the figure label and the window title is set to this value. Default to None.
figsize ((float, float), optional) – Width, height in inches. Default to [6.4, 4.8]
dpi (float, optional) – The resolution of the figure in dots-per-inch. Default to 100.0.
figsize – Width, height in inches. Default to [6.4, 4.8]
dpi – The resolution of the figure in dots-per-inch. Default to 100.0
facecolor (str, optional) – The background color.
edgecolor (str, optional) – The border color.
frameon (bool, optional) – If False, suppress drawing the figure frame. Defaults to True.
FigureClass (subclass of matplotlib.figure.Figure, optional) – Optionally use a custom Figure instance. Defaults to matplotlib.figure.Figure.
clear (bool, optional) – If True and the figure already exists, then it is cleared. Defaults to False.
tight_layout (bool, optional) – If False use subplotpars. If True adjust subplot parameters using tight_layout with default padding. When providing a dict containing the keys pad, w_pad, h_pad, and rect, the default tight_layout paddings will be overridden. Defaults to False.
constrained_layout (bool, optional) – If True use constrained layout to adjust positioning of plot elements. Like tight_layout, but designed to be more flexible. See Constrained Layout Guide for examples. (Note: does not work with add_subplot or subplot2grid.) Defaults to False.
fmt (str, optional) – A format string, e.g. ‘ro’ for red circles. See the Notes section for a full description of the format strings. Format strings are just an abbreviation for quickly setting basic line properties. All of these and more can also be controlled by keyword arguments. This argument cannot be passed as keyword. Defaults to None.
save_fig_path (str, optional) – Path of saving figure. Defaults to None.
**kwargs (matplotlib.lines.Line2D properties, optional) –
kwargs are used to specify properties like a line label (for auto legends), linewidth, antialiasing, marker face color. You can show the detailed information at `matplotlib.lines.Line2D
- Returns
Figure
- Return type
matplotlib.pyplot.Figure
Examples
>>> ship.draw_xy_trajectory(save_fig_path="test.png")
- load_simulation_result(time: List[float], u: List[float], v: List[float], r: List[float], x0: float = 0.0, y0: float = 0.0, psi0: float = 0.0)[source]
Load simulation result (time, u, v, r).
By running this, x, y and psi of this class are registered automatically.
- Parameters
time (list[float]) – Time list of simulation result.
u (list[float]) – List of axial velocity [m/s] in simulation result.
v (list[float]) – List of lateral velocity [m/s] in simulation result.
r (list[float]) – List of rate of turn [rad/s] in simulation result.
x0 (float, optional) – Inital position of X axis [m]. Defaults to 0.0.
y0 (float, optional) – Inital position of Y axis [m/s]. Defaults to 0.0.
psi0 (float, optional) – Inital azimuth [rad]. Defaults to 0.0.
Examples
>>> time_list = np.linspace(0.00, duration, num_of_sampling) >>> delta_list = np.full(len(time_list), 10 * np.pi / 180) >>> kt_params = KTParams(K=0.15, T=60.0) >>> result = kt.simulate_kt(kt_params, time_list, delta_list) >>> u_list = np.full(len(time_list), 20 * (1852.0 / 3600)) >>> v_list = np.zeros(len(time_list)) >>> r_list = result[0] >>> ship = ShipObj3dof(L = 180, B = 20) >>> ship.load_simulation_result(time_list, u_list, v_list, r_list) >>> print(ship.x, ship.y, ship.psi)
- npm: List[float]
- psi: List[float]
- r: List[float]
- register_simulation_result(time: List[float], u: List[float], v: List[float], r: List[float], x: List[float], y: List[float], psi: List[float])[source]
register simulation result (time, u, v, r, x, y, psi).
- Parameters
time (list[float]) – Time list of simulation result.
u (list[float]) – List of axial velocity [m/s] in simulation result.
v (list[float]) – List of lateral velocity [m/s] in simulation result.
r (list[float]) – List of rate of turn [rad/s] in simulation result.
x (list[float]) – List of position of X axis [m].
y (list[float]) – List of position of Y axis [m/s].
psi (list[float]) – List of inital azimuth [rad].
- time: List[float]
- u: List[float]
- v: List[float]
- x: List[float]
- y: List[float]
- δ: List[float]