environment#

Functionalities of environment objects.

This module provides functionalities for environment objects. Specifically, it contains classes and functions that perform computations related to environment models of natural and artificial bodies. Much of the functionality in this module concerns classes stored inside Body objects, a list of which is in turn stored in a SystemOfBodies object. Note that the classes in this module are rarely created manually, but are instead created by the functionality in the environment_setup module.

Functions#

save_vehicle_mesh_to_file(...[, ...])

Function to save the mesh used for a hypersonic local inclination analysis to a file.

save_vehicle_mesh_to_file(local_inclination_analysis_object: tudatpy.kernel.numerical_simulation.environment.HypersonicLocalInclinationAnalysis, output_directory: str, output_file_prefix: str = '') None#

Function to save the mesh used for a hypersonic local inclination analysis to a file.

Function to save the mesh used for a hypersonic local inclination analysis to a file. This function saves two files to the specified directory, with filenames: “ShapeFile.dat” and “SurfaceNormalFile.dat”, where these files names may be prefixed by an optional string (see below). The first of these files contains four columns defining the surface points that define mesh, with Column 0: point index; Column 1: x-position of point; Column 1: y-position of point; Column 2: z-position of point. The second file contains four columns with Column 0: point index; Column 1: x-component of surface normal; Column 1: y-position of surface normal; Column 2: z-position of surface normal.

Parameters:
  • local_inclination_analysis_object (HypersonicLocalInclinationAnalysis) – Object used to calculate the aerodynamics of the vehicle

  • output_directory (str) – Directory to which the files are to be saved

  • output_file_prefix (str, default='') – Optional prefix of output file names

Enumerations#

AerodynamicsReferenceFrames

Enumeration of reference frame identifiers typical for aerodynamic calculations.

AerodynamicCoefficientFrames

Enumeration of reference frames used for definition of aerodynamic coefficients.

AerodynamicCoefficientsIndependentVariables

Enumeration of the independent variables that can be used to compute aerodynamic coefficients.

class AerodynamicsReferenceFrames#

Enumeration of reference frame identifiers typical for aerodynamic calculations.

Enumeration of reference frame identifiers typical for aerodynamic calculations. Note that the frames are also defined in the absence of any aerodynamic forces and/or atmosphere. They define frames of a body w.r.t. a central body, with the details given by Mooij (1994). The chain of frames starts from the inertial frame, to the frame fixed to the central body (corotating), to the vertical frame (defined by the body’s relative position), the trajectory and aerodynamic frames (defined by the body’s relative velocity) and finally the body’s own body-fixed frame.

Members:

inertial_frame :

The global orientation (which is by definition inertial).

corotating_frame :

The body-fixed frame of the central body.

vertical_frame :

Frame with z-axis pointing towards origin of central body, the x-axis lies in the meridian plane and points towards the central-body-fixed z-axis (the y-axis completes the frame).

trajectory_frame :

The (airspeed-based) trajectory frame has the x-axis in the direction of the velocity vector relative to the atmosphere (airspeed-based velocity vector), z-axis lies in the vertical plane and points downwards (the y-axis completes the frame).

aerodynamic_frame :

The (airspeed-based) aerodynamic frame has the x-axis in the direction of the velocity vector relative to the atmosphere (airspeed-based velocity vector), z-axis co-linear with the aerodynamic lift vector, pointing in the opposite direction (the y-axis completes the frame)..

body_frame :

The body-fixed frame of the body itself.

property name#
class AerodynamicCoefficientFrames#

Enumeration of reference frames used for definition of aerodynamic coefficients.

Enumeration of reference frames used for definition of aerodynamic coefficients. There is a partial overlap between this enum and the AerodynamicsReferenceFrames. This enum combines a subset of those frames (which are typically used for aerodynamic coefficient definition), and a swap in sign. For instance, aerodynamic force coefficients are often defined positive along negative axes of the aerodynamic frame (drag, side force and lift coefficients)

Members:

positive_body_fixed_frame_coefficients :

The coefficients are defined in the body-fixed frame, with the directions the same as the body-fixed axes. For aerodynamic forces and moments, this results in the typical \(C_{x}, C_{y}, C_{y}\) (force) and \(C_{l}, C_{m}, C_{n}\) (moment) coefficients

negative_body_fixed_frame_coefficients :

Same as positive_body_fixed_frame_coefficients, but opposite in direction (so axes along negative body-fixed frame axes)

positive_aerodynamic_frame_coefficients :

Same as negative_aerodynamic_frame_coefficients, but opposite in direction (so axes along positive aerodynamic frame axes)

negative_aerodynamic_frame_coefficients :

The coefficients are defined in aerodynamic frame, with the directions the same as the negative axes. For aerodynamic forces, this results in the typical \(C_{D}, C_{S}, C_{D}\) force coefficients

property name#
class AerodynamicCoefficientsIndependentVariables#

Enumeration of the independent variables that can be used to compute aerodynamic coefficients.

Members:

mach_number_dependent :

Mach number of the propagated vehicle.

angle_of_attack_dependent :

Angle of attack of the propagated vehicle.

sideslip_angle_dependent :

Sideslip angle of the propagated vehicle.

altitude_dependent :

Altitude of the propagated vehicle.

time_dependent :

Current simulation epoch.

control_surface_deflection_dependent :

Angle of deflection of the control surface of the propagated vehicle.

undefined_independent_variable :

Can be used for a custom coefficient interface with other variables, at the expense of being able to use the FlightConditions class to automatically updates the aerodynamic coefficients during propagation.

property name#

Classes#

AerodynamicCoefficientInterface

Base class for computing the current aerodynamic coefficients of the body

HypersonicLocalInclinationAnalysis

FlightConditions

Object that calculates various state-derived quantities typically relevant for flight dynamics.

AtmosphericFlightConditions

Object that calculates various state-derived quantities typically relevant for flight dynamics, for flight in an atmosphere.

AerodynamicAngleCalculator

Object to calculate (aerodynamic) orientation angles, and frame transformations, from current vehicle state.

RotationalEphemeris

Object that stores the rotational state of the bodies.

VehicleSystems

Object used to store physical (hardware) properties of a vehicle.

Ephemeris

Object that computes the state of a body as a function of time

EngineModel

Body

Object that stores the environment properties and current state of a single body.

SystemOfBodies

Object that contains a set of Body objects and associated frame information.

class AerodynamicCoefficientInterface#

Base class for computing the current aerodynamic coefficients of the body

Base class for computing the current aerodynamic coefficients of the body. The implementation of the computation depends on the choice of aerodynamic coefficient model (see aerodynamic_coefficients for available options). During the propagation, this object is automatically updated to the current state by the AtmosphericFlightConditions object. The user may override the current aerodynamic coefficients when using, for instance, a custom aerodynamic guidance model (see here for an example). using the member functions of this class.

current_control_surface_force_coefficient_increment(self: tudatpy.kernel.numerical_simulation.environment.AerodynamicCoefficientInterface, control_surface_name: str) numpy.ndarray[numpy.float64[3, 1]]#

Function to get the contribution from a single control surface to the aerodynamic force coefficient, as compute by last call to update_full_coefficients()

Parameters:

control_surface_name (str) – The name of the control surface for which the contribution is to be retrieved

Returns:

Contribution from the requested control surface to the aerodynamic force coefficient

Return type:

np.ndarray

current_control_surface_moment_coefficient_increment(self: tudatpy.kernel.numerical_simulation.environment.AerodynamicCoefficientInterface, control_surface_name: str) numpy.ndarray[numpy.float64[3, 1]]#

Function to get the contribution from a single control surface to the aerodynamic moment coefficients, as compute by last call to update_full_coefficients()

Parameters:

control_surface_name (str) – The name of the control surface for which the contribution is to be retrieved

Returns:

Contribution from the requested control surface to the aerodynamic moment coefficients

Return type:

np.ndarray

set_control_surface_increments(self: tudatpy.kernel.numerical_simulation.environment.AerodynamicCoefficientInterface, control_surface_list: dict[str, tudat::aerodynamics::ControlSurfaceIncrementAerodynamicInterface]) None#

No documentation found.

update_coefficients(self: tudatpy.kernel.numerical_simulation.environment.AerodynamicCoefficientInterface, independent_variables: list[float], time: float) None#

Function to update the aerodynamic coefficients of the body only

Function to update the aerodynamic coefficients of the body only (without the control surface contribution), based on the current state. This function may be called by the user, but will set only the current_force_coefficients and current_moment_coefficients (while leaving the current_control_surface_free_force_coefficients and current_control_surface_free_moment_coefficients unchanged)

Parameters:
  • independent_variables (list[float]) – List of inputs from which the aerodynamic coefficients are to be computed, with each entry corresponding to the value of the physical variable defined by the independent_variable_names attribute.

  • time (float) – Current time (in seconds since J2000)

Returns:

Contribution from the requested control surface to the aerodynamic moment coefficients

Return type:

np.ndarray

update_full_coefficients(self: tudatpy.kernel.numerical_simulation.environment.AerodynamicCoefficientInterface, independent_variables: list[float], control_surface_independent_variables: dict[str, list[float]], time: float, check_force_contribution: bool = True) None#

Function to update the aerodynamic coefficients, from both the body and its control surfaces

Function to update the aerodynamic coefficients of both the body and its control surfaces, based on the current state. This function will call the update_coefficients() function to update the body coefficients. This function may be called by the user, and will set the following attributes: current_force_coefficients, current_moment_coefficients , current_control_surface_free_force_coefficients and current_control_surface_free_moment_coefficients. In addition, it will modify the coefficients returned by the current_control_surface_force_coefficient_increment() and current_control_surface_moment_coefficient_increment() functions

Parameters:
  • independent_variables (list[float]) – List of inputs from which the aerodynamic coefficients of the body are to be computed, with each entry corresponding to the value of the physical variable defined by the independent_variable_names attribute.

  • control_surface_independent_variables (dict[str,list[float]]) – List of inputs from which the control surface aerodynamic coefficients are to be computed (with dictionary key the control surface name), with each entry corresponding to the value of the physical variable defined by the control_surface_independent_variable_names attribute.

  • time (float) – Current time (in seconds since J2000)

  • check_force_contribution (bool, default = True) – Boolean that determines if the force contribution to the aerodynamic moments should be added. Note that this input is only used if the add_force_contribution_to_moments attribute is set to True.

property control_surface_independent_variable_names#

read-only

List of independent variables from which the aerodynamic coefficients of each control surface are computed, with dictionary key being the control surface name (e.g. required input to update_full_coefficients() function).

Type:

dict[str,list[AerodynamicCoefficientsIndependentVariables]]

property current_coefficients#

read-only

Concatenation of current_force_coefficients and current_moment_coefficients

Type:

np.ndarray

property current_control_surface_free_force_coefficients#

read-only

Same as current_force_coefficients, but without contribution (if any) from control surfaces

Type:

np.ndarray

property current_control_surface_free_moment_coefficients#

read-only

Same as current_moment_coefficients, but without contribution (if any) from control surfaces

Type:

np.ndarray

property current_force_coefficients#

read-only

The current aerodynamic force coefficients, in the frame defined by the force_coefficient_frame attribute, as computed by the last call to the update_coefficients() function.

Type:

np.ndarray

property current_moment_coefficients#

read-only

The current aerodynamic moment coefficients, in the frame defined by the moment_coefficient_frame attribute, as computed by the last call to the update_coefficients() function.

Type:

np.ndarray

property force_coefficient_frame#

read-only

Reference frame in which the current_force_coefficients are defined

Type:

AerodynamicCoefficientFrames

property independent_variable_names#

read-only

List of independent variables from which the aerodynamic coefficients are computed (e.g. required input to update_coefficients() function).

Type:

list[AerodynamicCoefficientsIndependentVariables]

property moment_coefficient_frame#

read-only

Reference frame in which the current_moment_coefficients are defined

Type:

AerodynamicCoefficientFrames

property reference_area#

read-only

The aerodynamic reference area \(A\) of the coefficients

Type:

float

class HypersonicLocalInclinationAnalysis#
__init__(self: tudatpy.kernel.numerical_simulation.environment.HypersonicLocalInclinationAnalysis, independent_variable_points: list[list[float]], body_shape: tudatpy.kernel.math.geometry.SurfaceGeometry, number_of_lines: list[int], number_of_points: list[int], invert_orders: list[bool], selected_methods: list[list[int]], reference_area: float, reference_length: float, moment_reference_point: numpy.ndarray[numpy.float64[3, 1]], save_pressure_coefficients: bool = False) None#

Class constructor, taking the shape of the vehicle, and various analysis options as input.

Parameters:
  • independent_variable_points (list[list[float]]) – List containing three lists, with each sublist containing the data points of each of the independent variables for the coefficient generation. The physical meaning of each of the three independent variables is: 0 = mach number, 1 = angle of attack, 2 = angle of sideslip. Each of the subvectors must be sorted in ascending order.

  • body_shape (SurfaceGeometry) – Class that defines the shape of the vehicle as a continuous surface. The local inclination analysis discretizes the surface of the vehicle into quadrilateral panels, defined by the other inputs to this constructor. In case the tudat.geometry.SurfaceGeometry object is made up of multiple sub-shapes, different settings may be used for each

  • number_of_lines (List[ float ]) – Number of discretization points in the first independent surface variable of each of the subparts of body_shape. The size of this list should match the number of parts of which the body_shape is composed. The first independent variable of a subpart typically runs along the longitudinal vehicle direction

  • number_of_points (List[ float ]) – Number of discretization points in the second independent surface variable of each of the subparts of body_shape. The size of this list should match the number of parts of which the body_shape is composed. The first independent variable of a subpart typically runs along the lateral vehicle direction

  • invert_orders (List[ bool ]) – Booleans to denote whether the surface normals of the panels of each discretized body_shape subpart are to be inverted (i.e. inward-facing->outward facing or vice versa). The size of this list should match the number of parts of which the body_shape is composed.

  • selected_methods (List[ List[ int ] ]) – Double list of selected local inclination methods, the first index (outer list) represents compression or expansion (0 and 1), the second index (inner list) denotes the vehicle part index. The size of this inner list should match the number of parts of which the body_shape is composed. The int defining the method type is interpreted as follows. For the compression methods, the following are available: * 0: Newtonian Method. * 1: Modified Newtonian. * 2 and 3: not available at this moment. * 4: Tangent-wedge method. * 5: Tangent-cone method. * 6: Modified Dahlem-Buck method. * 7: VanDyke unified pressure method. * 8: Smyth Delta Wing method. * 9: Hankey flat surface method The expansion method has the following options: * 0: Vacuum Pressure coefficient method. * 1: Zero Pressure function. * 4: High Mach base pressure method. * 3 or 5: Prandtl-Meyer method. * 6: ACM empirical pressure coefficient.

  • reference_area (float) – Reference area used to non-dimensionalize aerodynamic forces and moments.

  • moment_reference_point (numpy.ndarray) – Reference point wrt which aerodynamic moments are calculated.

  • save_pressure_coefficients (Bool) – Boolean denoting whether to save the pressure coefficients that are computed to files

clear_data(self: tudatpy.kernel.numerical_simulation.environment.HypersonicLocalInclinationAnalysis) None#
class FlightConditions#

Object that calculates various state-derived quantities typically relevant for flight dynamics.

Object that calculates various state-derived quantities typically relevant for flight dynamics, such as latitude, longitude, altitude, etc. It also contains an AerodynamicAngleCalculator that computes derived angles (flight path, heading angle, etc.). This object is limited to non-atmospheric flight. For flight through Body objects endowed with an atmosphere model, the derived class AtmosphericFlightConditions is used. This object is stored inside a Body object, and represents the flight conditions of a single body w.r.t. a single central body.

update_conditions(self: tudatpy.kernel.numerical_simulation.environment.FlightConditions, current_time: float) None#
property aerodynamic_angle_calculator#

read-only

The object that is responsible for computing various relevant flight dynamics angles and frame rotations.

Type:

AerodynamicAngleCalculator

property altitude#

read-only

The current time, at which this object was last updated

Type:

float

property body_centered_body_fixed_state#

read-only

Cartesian translational state, expressed in a frame centered at, and fixed to, the central body. Note that, due to the rotation of the central body, the norm of the body-fixed, body-centered, velocity differs from the norm of the inertial body-centered velocity.

Type:

numpy.ndarray

property geodetic_latitude#

read-only

The body-fixed geographic latitude of the body w.r.t. its central body.

Type:

float

property latitude#

read-only

The body-fixed geographic latitude of the body w.r.t. its central body.

Type:

float

property longitude#

read-only

The body-fixed longitude of the body w.r.t. its central body.

Type:

float

property time#

read-only

The current time, at which this object was last updated

Type:

float

class AtmosphericFlightConditions#

Object that calculates various state-derived quantities typically relevant for flight dynamics, for flight in an atmosphere.

Object that calculates various state-derived quantities typically relevant for flight dynamics, for flight in an atmosphere, such as latitude, longitude, altitude, density, Mach number etc. It also contains an AerodynamicAngleCalculator that computes derived angles (flight path, heading angle, etc.). This object is derived from FlightConditions, which performs computations for non-atmospheric flight only. This object is stored inside a Body object, and represents the flight conditions of a single body w.r.t. a single central body.

property aero_coefficient_independent_variables#

read-only

List of current values of independent variables of aerodynamic coefficients. This list is only defined if the body has an AerodynamicCoefficientInterface that has dependencies on environmental variables (e.g. Mach number, angle of attack, etc.).

Type:

numpy.ndarray

property aerodynamic_coefficient_interface#

read-only

Object extracted from the same Body object as this AtmosphericFlightConditions object, which defines the aerodynamic coefficients.

Type:

AerodynamicCoefficientInterface

property airspeed#

read-only

The airspeed of the body w.r.t. the atmosphere.

Type:

float

property airspeed_velocity#

read-only

The velocity vector of the body w.r.t. the freestream atmosphere (e.g. vectorial counterpart of airspeed).

Type:

numpy.ndarray

property control_surface_aero_coefficient_independent_variables#

read-only

List of lists current values of independent variables of aerodynamic coefficients for control surfaces. The outer list defines the control surface, the inner list the values of the independent variables. This list is only defined if the body has an AerodynamicCoefficientInterface with control surfaces that have dependencies on environmental variables (e.g. Mach number, angle of attack, etc.).

Type:

numpy.ndarray

property density#

read-only

The freestream atmospheric density at the body’s current location.

Type:

float

property dynamic_pressure#

read-only

The freestream atmospheric dynamic pressure at the body’s current location.

Type:

float

property mach_number#

read-only

The freestream Mach number of the body.

Type:

float

property pressure#

read-only

The freestream atmospheric static pressure at the body’s current location.

Type:

float

property speed_of_sound#

read-only

The freestream atmospheric speed of sound at the body’s current location.

Type:

float

property temperature#

read-only

The freestream atmospheric temperature at the body’s current location.

Type:

float

class AerodynamicAngleCalculator#

Object to calculate (aerodynamic) orientation angles, and frame transformations, from current vehicle state.

Object to calculate (aerodynamic) orientation angles (list given by the AerodynamicsReferenceFrameAngles enum) and transformations between frames (list given by the AerodynamicsReferenceFrames enum) from current vehicle state.

get_angle(self: tudatpy.kernel.numerical_simulation.environment.AerodynamicAngleCalculator, angle_type: tudatpy.kernel.numerical_simulation.environment.AerodynamicsReferenceFrameAngles) float#

Function to get a single orientation angle

Function to get a single orientation angle. This function is meant to be used only during a numerical propagation, in particular for the definition of a custom (e.g. guidance) model.

Parameters:

original_frame (AerodynamicsReferenceFrameAngles) – The identifier for the angle that is to be returnd

Returns:

Value of requested angle

Return type:

double

get_rotation_matrix_between_frames(self: tudatpy.kernel.numerical_simulation.environment.AerodynamicAngleCalculator, original_frame: tudatpy.kernel.numerical_simulation.environment.AerodynamicsReferenceFrames, target_frame: tudatpy.kernel.numerical_simulation.environment.AerodynamicsReferenceFrames) numpy.ndarray[numpy.float64[3, 3]]#

Function to get the rotation matrix between two frames.

Function to get the rotation matrix between two frames. This function is meant to be used only during a numerical propagation, in particular for the definition of a custom (e.g. guidance) model.

Parameters:
Returns:

Rotation matrix \(\mathbf{R}^{B/A}\) from frame \(A\) to frame B

Return type:

np.ndarray

set_body_orientation_angle_functions(self: tudatpy.kernel.numerical_simulation.environment.AerodynamicAngleCalculator, angle_of_attack_function: Callable[[], float] = None, angle_of_sideslip_function: Callable[[], float] = None, bank_angle_function: Callable[[], float] = None, angle_update_function: Callable[[float], None] = None, silence_warnings: bool = False) None#
set_body_orientation_angles(self: tudatpy.kernel.numerical_simulation.environment.AerodynamicAngleCalculator, angle_of_attack: float = nan, angle_of_sideslip: float = nan, bank_angle: float = nan, silence_warnings: bool = False) None#
class RotationalEphemeris#

Object that stores the rotational state of the bodies.

Object that stores the rotational state of the bodies. This object can be used to calculate rotation matrices, which are used to transform coordinates between reference frames.

angular_velocity_in_body_fixed_frame(self: tudatpy.kernel.numerical_simulation.environment.RotationalEphemeris, time: float) numpy.ndarray[numpy.float64[3, 1]]#

Function to get the body’s angular velocity vector, expressed in the body-fixed frame.

Function to get the body’s angular velocity vector \(\boldsymbol{\omega}^{(B)}\), expressed in the body-fixed frame \(B\). The calculation of the angular velocity depends on the specific rotation model that has been defined, either from an a priori definition (see rotation_model submodule) or from processing the results of propagation of the rotational equations of motion. Note that when numerically propagating rotational dynamics, this angular velocity vector is typically directly defined in the last three entries of the state vector.

Parameters:

current_time (float) – The time (in seconds since epoch J2000, TDB time scale) at which the angular velocity vector is evaluated

Returns:

Angular velocity vector of the body \(\boldsymbol{\omega}^{(B)}\) expressed in the body-fixed frame \(B\)

Return type:

np.ndarray

angular_velocity_in_inertial_frame(self: tudatpy.kernel.numerical_simulation.environment.RotationalEphemeris, time: float) numpy.ndarray[numpy.float64[3, 1]]#

Function to get the body’s angular velocity vector, expressed in the inertial frame.

Function to get the body’s angular velocity vector \(\boldsymbol{\omega}^{(I)}\), expressed in the body-fixed frame \(I\). This quantity is computed from \(\mathbf{R}^{I/B}\boldsymbol{\omega}^{(B)}\), see the angular_velocity_in_body_fixed_frame and body_fixed_to_inertial_rotation functions.

Parameters:

current_time (float) – The time (in seconds since epoch J2000, TDB time scale) at which the angular velocity vector is evaluated

Returns:

Angular velocity vector of the body \(\boldsymbol{\omega}^{(B)}\) expressed in the body-fixed frame \(B\)

Return type:

np.ndarray

body_fixed_to_inertial_rotation(self: tudatpy.kernel.numerical_simulation.environment.RotationalEphemeris, time: float) numpy.ndarray[numpy.float64[3, 3]]#

Function to get rotation matrix from body-fixed frame to inertial frame over time.

Function to get rotation matrix from body-fixed (target) frame to inertial (base) frame over time. The calculation of this rotation matrix depends on the specific rotation model that has been defined, either from an a priori definition (see rotation_model submodule) or from processing the results of propagation of the rotational equations of motion.

Parameters:

current_time (float) – The time (in seconds since epoch J2000, TDB time scale) at which the rotation matrix is evaluated

Returns:

Rotation matrix \(\mathbf{R}^{I/B}\) from body-fixed frame \(B\) to inertial frame I

Return type:

np.ndarray

inertial_to_body_fixed_rotation(self: tudatpy.kernel.numerical_simulation.environment.RotationalEphemeris, time: float) numpy.ndarray[numpy.float64[3, 3]]#

Function to get rotation matrix from inertial frame to body-fixed frame over time.

Function computes the inverse (equal to transpose) rotation of the body_fixed_to_inertial_rotation function.

Parameters:

current_time (float) – The time (in seconds since epoch J2000, TDB time scale) at which the rotation matrix is evaluated

time_derivative_body_fixed_to_inertial_rotation(self: tudatpy.kernel.numerical_simulation.environment.RotationalEphemeris, time: float) numpy.ndarray[numpy.float64[3, 3]]#

Function to get time derivative of rotation matrix from body-fixed frame to inertial frame over time.

Function to get time derivative of rotation matrix from body-fixed frame to inertial frame over time (see body_fixed_to_inertial_rotation), denoted \(\dot{\mathbf{R}}^{(I/B)}\),

Parameters:

current_time (float) – The time (in seconds since epoch J2000, TDB time scale) at which the rotation matrix derivative is evaluated

Returns:

Rotation matrix \(\dot{\mathbf{R}}^{I/B}\) from body-fixed frame \(B\) to inertial frame I

Return type:

np.ndarray

time_derivative_inertial_to_body_fixed_rotation(self: tudatpy.kernel.numerical_simulation.environment.RotationalEphemeris, time: float) numpy.ndarray[numpy.float64[3, 3]]#

Function to get time derivative of rotation matrix from inertial frame to body-fixed frame over time.

Function to get time derivative of rotation matrix from inertial frame to body-fixed frame over time (see inertial_to_body_fixed_rotation), denoted \(\dot{\mathbf{R}}^{(B/I)}\),

Parameters:

current_time (float) – The time (in seconds since epoch J2000, TDB time scale) at which the rotation matrix derivative is evaluated

Returns:

Rotation matrix \(\dot{\mathbf{R}}^{B/I}\) from inertial frame I to body-fixed frame \(B\)

Return type:

np.ndarray

property body_fixed_frame_name#

read-only

The identifier of the body-fixed frame, used in other parts of the simulation to identify it.

Type:

str

property inertial_frame_name#

read-only

The identifier of the inertial frame, used in other parts of the simulation to identify it.

Type:

str

class VehicleSystems#

Object used to store physical (hardware) properties of a vehicle.

get_control_surface_deflection(self: tudatpy.kernel.numerical_simulation.environment.VehicleSystems, control_surface_id: str) float#

No documentation found.

get_engine_model(self: tudatpy.kernel.numerical_simulation.environment.VehicleSystems, engine_name: str) tudat::system_models::EngineModel#

Function to retrieve an engine model from the vehicle

Parameters:

engine_name (str) – The identifier for the engine model that is to be retrieved

Returns:

Model for the engine that is requested

Return type:

EngineModel

set_control_surface_deflection(self: tudatpy.kernel.numerical_simulation.environment.VehicleSystems, control_surface_id: str, deflection_angle: float) None#

Function to set the current deflection of an aerodynamic control surface.

Function to set the current deflection of an aerodynamic control surface, identified by its name. To set the control surface deflection, the control surface has to exist. A control surface is created whenever control surfaces are defined in a body’s aerodynamic coefficient interface.

Parameters:
  • control_surface_id (str) – The identified (name) of the given control surface

  • deflection_angle (float) – The deflection (in radians) that the control surface is to be set to. This will typically influence the aerodynamic coefficients of the vehicle

set_transponder_turnaround_ratio(self: tudatpy.kernel.numerical_simulation.environment.VehicleSystems, transponder_ratio_per_uplink_and_downlink_frequency_band: dict[tuple[tudat::observation_models::FrequencyBands, tudat::observation_models::FrequencyBands], float]) None#

No documentation found.

class Ephemeris#

Object that computes the state of a body as a function of time

Object (typically stored inside a Body object) that computes the state of a body as a function of time, both outside of a propagation, and during a propagation if the given body’s translational state is not propagated. Note that this object computes the state w.r.t. its own origin (defined by frame_origin), which need not be the same as the global frame origin of the environment.

cartesian_position(self: tudatpy.kernel.numerical_simulation.environment.Ephemeris, current_time: float) numpy.ndarray[numpy.float64[3, 1]]#

As cartesian_state, but only the three position components

Parameters:

current_time (float) – Time (in seconds since J2000 in TDB time scale) at which the state is to be computed.

Returns:

Requested Cartesian position

Return type:

np.ndarray

cartesian_state(self: tudatpy.kernel.numerical_simulation.environment.Ephemeris, current_time: float) numpy.ndarray[numpy.float64[6, 1]]#

This function returns the Cartesian state (position and velocity) at the given time, w.r.t. the frame_origin.

Parameters:

current_time (float) – Time (in seconds since J2000 in TDB time scale) at which the state is to be computed.

Returns:

Requested Cartesian state

Return type:

np.ndarray

cartesian_velocity(self: tudatpy.kernel.numerical_simulation.environment.Ephemeris, current_time: float) numpy.ndarray[numpy.float64[3, 1]]#

As cartesian_state, but only the three velocity components

Parameters:

current_time (float) – Time (in seconds since J2000 in TDB time scale) at which the state is to be computed.

Returns:

Requested Cartesian velocity

Return type:

np.ndarray

property frame_orientation#

read-only

Name of the frame orientation w.r.t which this object provides its states

Type:

str

property frame_origin#

read-only

Name of the reference body/point w.r.t. which this object provides its states

Type:

str

class EngineModel#
class Body#

Object that stores the environment properties and current state of a single body.

Object that stores the environment properties and current state of a single celestial body (natural or artificial). Each separate environment model (gravity field, ephemeris, etc.) is stored as a member object in this class. During each time step, the Body gets updated to the current time/propagated state, and the current properties, in as much as they are time-dependent, can be extracted from this object

get_ground_station(self: tudatpy.kernel.numerical_simulation.environment.Body, station_name: str) tudatpy.kernel.numerical_simulation.environment.GroundStation#

This function extracts a ground station object from the body.

This function extracts a ground station object, for a station of a given name, from the body. If no station of this name exists, an exception is thrown

Parameters:

station_name (str) – Name of the ground station that is to be retrieved.

Returns:

Ground station object of the station of requested name

Return type:

GroundStation

set_constant_mass(self: tudatpy.kernel.numerical_simulation.environment.Body, mass: float) None#
state_in_base_frame_from_ephemeris(self: tudatpy.kernel.numerical_simulation.environment.Body, time: float) numpy.ndarray[numpy.float64[6, 1]]#
property aerodynamic_coefficient_interface#

Object defining the aerodynamic coefficients of a vehicle (force-only, or force and moment) as a function of any number of independent variables. Depending on the selected type of model, the type of this attribute is of type AerodynamicCoefficientInterface, or a derived class thereof.

Type:

AerodynamicCoefficientInterface

property atmosphere_model#

Atmosphere model of this body, used to calculate density, temperature, etc. at a given state/time. Depending on the selected type of model, the type of this attribute is of type AtmosphereModel, or a derived class thereof.

Type:

AtmosphereModel

property body_fixed_angular_velocity#

read-only

Angular velocity vector of the body, expressed in body-fixed frame (see inertial_to_body_fixed_frame).

Type:

numpy.ndarray

property body_fixed_to_inertial_frame#

read-only

The rotation from this Body’s body-fixed frame to inertial frame (see inertial_to_body_fixed_frame).

Type:

numpy.ndarray

property body_fixed_to_inertial_frame_derivative#

read-only

Time derivative of rotation matrix from this Body’s body-fixed frame to inertial frame (see inertial_to_body_fixed_frame).

Type:

numpy.ndarray

property ephemeris#

Ephemeris model of this body, used to calculate its current state as a function of time. Depending on the selected type of model, the type of this attribute is of type Ephemeris, or a derived class thereof.

Type:

Ephemeris

property flight_conditions#

Object used to calculated and store the current flight conditions of a vehicle (altitude, latitude, longitude, flight-path angle, etc.) w.r.t. a central body. In case the central body contains an atmosphere, this object also stores current local density, Mach number, etc. This object is typically used for aerodynamic accelerations, guidance models or other central-body-related custom models.

Type:

FlightConditions

property gravitational_parameter#

read-only

Attribute of convenience, equivalent to .gravity_field_model.gravitational_parameter

Type:

float

property gravity_field_model#

Gravity field model of this body, used to define the exterior gravitational potential, and its gradient(s). Depending on the selected type of model, the type of this attribute is of type GravityFieldModel, or a derived class thereof.

Type:

GravityFieldModel

property ground_station_list#

Dictionary of all ground stations that exist in the body, with dictionary key being the name of the station, and the ground station object the key of the dictionary.

Type:

dict[str,GroundStation]

property inertia_tensor#

The current inertia tensor \(\mathbf{I}\) of the vehicle, as used in the calculation of (for instance) the reponse to torques. This attribute is a shorthand for accessing the inertia tensor as computed/stored in the rigid_body_properties attribute. For certain types of rigid-body properties, this attribute cannot be used to (re)set the current mass.

Unlike the attributes containing the state, orientation, angular velocity of the Body, this attribute may be used to retrieve the state during the propagation and to define the mass of a vehicle.

Type:

np.ndarray

property inertial_angular_velocity#

read-only

Angular velocity vector of the body, expressed in inertial frame (see inertial_to_body_fixed_frame).

Type:

numpy.ndarray

property inertial_to_body_fixed_frame#

read-only

The rotation from inertial frame (with global frame orientation) to this Body’s body-fixed frame. The rotation is always returned here as a rotation matrix. If the body’s rotational state is numerically propagated, this property gets extracted from the propagated state vector. If it is not propagated, the state is extracted from this body’s rotational ephemeris.

Note

This function is only valid during the numerical propagation if any aspects of the dynamics or dependent variables require the body’s rotational state.

Type:

numpy.ndarray

property inertial_to_body_fixed_frame_derivative#

read-only

Time derivative of rotation matrix from inertial frame to this Body’s body-fixed frame (see inertial_to_body_fixed_frame).

Type:

numpy.ndarray

property mass#

The current mass \(m\) of the vehicle, as used in the calculation of non-conservative acceleration. This attribute is a shorthand for accessing the mass as computed/stored in the rigid_body_properties attribute. For certain types of rigid-body properties, this attribute cannot be used to (re)set the current mass. If the body has no rigid_body_properties, and this function is used to set a mass, a new object is automatically created, with settings analogous to the the constant_rigid_body_properties() setting.

Unlike the attributes containing the state, orientation, angular velocity of the Body, this attribute may be used to retrieve the state during the propagation and to define the mass of a vehicle.

Type:

float

property position#

read-only

The translational position of the Body, as set during the current step of the numerical propagation (see state).

Type:

numpy.ndarray

property rigid_body_properties#

No documentation found.

property rotation_model#

Object defining the orientation of a body, used to calculate the rotation to/from a body-fixed frame (and its derivate). Depending on the selected type of model, the type of this attribute is of type RotationalEphemeris, or a derived class thereof.

Type:

RotationalEphemeris

property shape_model#

Shape model of this body, used to define the exterior shape of the body, for instance for the calculation of vehicle’s altitude. Depending on the selected type of model, the type of this attribute is of type BodyShapeModel, or a derived class thereof.

Type:

BodyShapeModel

property state#

read-only

The translational state of the Body, as set during the current step of the numerical propagation. The translational state stored here is always in Cartesian elements, w.r.t. the global frame origin, with axes along the global frame orientation. If the body’s translational state is numerically propagated, this property gets extracted from the propagated state vector. If it is not propagated, the state is extracted from this body’s ephemeris. In both cases, any required state transformations are automatically applied. Note that this function is only valid during the numerical propagation if any aspects of the dynamics or dependent variables require the body’s state.

Type:

numpy.ndarray

property system_models#

Object used to store physical (hardware) properties of a vehicle, such as engines, control surfaces, etc. This object is typically created automatically whenever such a hardware model needs to be assigned to a vehicle.

Type:

VehicleSystems

property velocity#

read-only

The translational velocity of the Body, as set during the current step of the numerical propagation (see state).

Type:

numpy.ndarray

class SystemOfBodies#

Object that contains a set of Body objects and associated frame information.

Object that contains a set of Body objects and associated frame information. This object stored the entire environment for a typical Tudat numerical simulation, and is fundamental for the overall Tudat architecture.

add_body(self: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies, body_to_add: tudatpy.kernel.numerical_simulation.environment.Body, body_name: str, process_body: bool = 1) None#

This function adds an existing body, which the user has separately created, to the SystemOfBodies.

Parameters:
  • body_to_add (Body) – Body object that is to be added.

  • body_name (numpy.ndarray) – Name of the Body that is to be added.

  • process_body (bool, default=True) –

    Variable that defines whether this new Body will have its global frame origin/orientation set to conform to rest of the environment.

    Warning

    Only in very rare cases should this variable be anything other than True. Users are recommended to keep this default value intact.

create_empty_body(self: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies, body_name: str, process_body: bool = 1) None#

This function creates a new empty body.

This function creates a new empty body, and adds it to the SystemOfBodies. Since the body is empty, it will not have any environment models defined. These must all be added manually by a user.

Parameters:
  • body_name (string) – Name of the Body that is to be added

  • process_body (bool, default=True) –

    Variable that defines whether this new Body will have its global frame origin/orientation set to conform to rest of the environment.

    Warning

    Only in very rare cases should this variable be anything other than True. Users are recommended to keep this default value intact.

Examples

This function is often used early on in the environment creation segment of a simulation, following the creation of a SystemOfBodies from the default settings for celestial bodies.

# Define string names for bodies to be created from default.
bodies_to_create = ["Sun", "Earth", "Moon", "Mars", "Venus"]

# Use "Earth"/"J2000" as global frame origin and orientation.
global_frame_origin = "Earth"
global_frame_orientation = "J2000"

# Create default body settings, usually from `spice`.
body_settings = environment_setup.get_default_body_settings(
    bodies_to_create,
    global_frame_origin,
    global_frame_orientation)

# Create system of selected celestial bodies
bodies = environment_setup.create_system_of_bodies(body_settings)

# Create vehicle objects.
bodies.create_empty_body("Delfi-C3")
does_body_exist(self: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies, body_name: str) bool#

No documentation found.

get(self: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies, body_name: str) tudatpy.kernel.numerical_simulation.environment.Body#

This function extracts a single Body object from the SystemOfBodies.

Parameters:

body_name (str) – Name of the Body that is to be retrieved.

Returns:

Body object of the requested name

Return type:

Body

get_body(self: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies, body_name: str) tudatpy.kernel.numerical_simulation.environment.Body#

Deprecated version of get()

global_frame_orientation(self: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies) str#

No documentation found.

global_frame_origin(self: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies) str#

No documentation found.

list_of_bodies(self: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies) list[str]#

No documentation found.

remove_body(self: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies, body_name: str) None#

This function removes an existing body from the SystemOfBodies.

Warning

This function does not necessarily delete the Body object, it only removes it from this object. If any existing models in the simulation refer to this Body, it will persist in memory.

Parameters:

body_name (numpy.ndarray) – Name of the Body that is to be removed.