LRO Albedo and Thermal IR (Complex)
Overview
This example demonstrates advanced computation of albedo and thermal infrared accelerations using the object-oriented interface of the pyRTX library.
Key Features
Spatially-varying albedo: Uses a grid of albedo values across the planet surface
Spatially-varying temperature: Incorporates temperature variations with location
Digital elevation model: Represents planetary topography for accurate shadowing
OBJ file format: Loads detailed shape models from standard mesh files
This advanced approach is essential for high-precision orbit determination and propagation where planetary radiation pressure effects need to be accurately modeled.
When to Use This Approach
Use this complex method when:
High-fidelity force modeling is required
Surface property data is available
Topographic effects are significant
Comparing with the simple uniform model shows non-negligible differences
Code
1### ------------------------------------------------------------------------------------------------------- ###
2
3# Example purpose:
4# Show the object-oriented interface of the pyRTX library
5#
6# Example case:
7# Show how to compute albedo and thermal-ir accelerations
8# Here we use an albedo and temperature grid for the planet. We also represent the shape of the planet
9# with a digital elevation model stored in a .obj file.
10
11### ------------------------------------------------------------------------------------------------------- ###
12### IMPORTS
13
14import numpy as np
15import spiceypy as sp
16import xarray as xr
17import matplotlib.pyplot as plt
18import timeit
19
20from pyRTX.classes.Spacecraft import Spacecraft
21from pyRTX.classes.Planet import Planet, TemperatureGrid, AlbedoGrid
22from pyRTX.classes.Radiation import Albedo, Emissivity
23from pyRTX.classes.Precompute import Precompute
24from pyRTX.classes.LookUpTable import LookUpTable
25from pyRTX.core.analysis_utils import epochRange2
26
27from numpy import floor, mod
28
29### ------------------------------------------------------------------------------------------------------- ###
30### INPUTS
31
32# NOTE: before running this script you should generate the input lutfile running the
33# example 'compute_lut.py' using type = 'accel'.
34
35ref_epc = "2010 may 10 09:25:00"
36duration = 10000 # seconds
37timestep = 100
38METAKR = '../example_data/LRO/metakernel_lro.tm' # metakernel
39obj_path = '../example_data/LRO/' # folder with .obj files
40lutfile = 'luts/lro_accel_lut.nc' # lookup table file
41base_flux = 1361.5
42ref_radius = 1737.4
43n_cores = 10
44
45# The spacecraft mass can be a float, int or a xarray with times and values [kg]
46sc_mass = xr.open_dataset('mass/lro_mass.nc')
47sc_mass.load()
48sc_mass.close()
49
50# Planet representation (set None for representing the planet as a simple sphere.
51# Load an obj for representing it with digital elevation models)
52fromFile = 'moon_obj/fib_1e4_v2.obj'
53
54# Albedo grids
55# Reference: https://ode.rsl.wustl.edu/moon/pagehelp/Content/Missions_Instruments/LRO/LOLA/GDR/GDRDAM.htm
56alb_grid = 'grids/bond_albedo.npy'
57alb_lon = 'grids/ldam_4_lon.npy'
58alb_lat = 'grids/ldam_4_lat.npy'
59
60# Temperature grids
61# Reference: https://doi.org/10.1016/j.icarus.2016.08.012
62temp_grid = 'grids/temp.npy'
63temp_lon = 'grids/temp_lon.npy'
64temp_lat = 'grids/temp_lat.npy'
65
66### ------------------------------------------------------------------------------------------------------- ###
67### OBJECTS DEFINITION
68
69# Time initialization
70tic = timeit.default_timer()
71
72# Load the metakernel containing references to the necessary SPICE frames
73sp.furnsh(METAKR)
74
75# Define epochs
76epc_et0 = sp.str2et( ref_epc )
77epc_et1 = epc_et0 + duration
78epochs = epochRange2(startEpoch = epc_et0, endEpoch = epc_et1, step = timestep)
79
80# Define the Spacecraft Object (Refer to the class documentation for further details)
81lro = Spacecraft( name = 'LRO',
82
83 base_frame = 'LRO_SC_BUS', # Name of the spacecraft body-fixed frame
84
85 mass = sc_mass,
86
87 spacecraft_model = { # Define a spacecraft model
88
89 'LRO_BUS': {
90 'file' : obj_path + 'bus_rotated.obj', # .obj file of the spacecraft component
91 'frame_type': 'Spice', # type of frame (can be 'Spice' or 'UD'
92 'frame_name': 'LRO_SC_BUS', # Name of the frame
93 'center': [0.0,0.0,0.0], # Origin of the component
94 'diffuse': 0.1, # Diffuse reflect. coefficient
95 'specular': 0.3, # Specular reflect. coefficient
96 },
97
98 'LRO_SA': {
99 'file': obj_path + 'SA_recentred.obj',
100 'frame_type': 'Spice',
101 'frame_name': 'LRO_SA',
102 'center': [-1,-1.1, -0.1],
103 'diffuse': 0,
104 'specular': 0.3,
105 },
106
107
108 'LRO_HGA': {
109 'file': obj_path + 'HGA_recentred.obj',
110 'frame_type': 'Spice',
111 'frame_name': 'LRO_HGA',
112 'center':[-0.99, -0.3, -3.1],
113 'diffuse': 0.2,
114 'specular': 0.1,
115 },
116 }
117 )
118
119# Define the Moon object
120moon = Planet( fromFile = fromFile,
121 radius = ref_radius,
122 name = 'Moon',
123 bodyFrame = 'MOON_PA',
124 sunFixedFrame = 'GSE_MOON',
125 units = 'km',
126 subdivs = 5,
127 )
128
129# Define the Albedo grid object
130Lon = np.load(alb_lon)
131Lat = np.load(alb_lat)
132
133ALB = AlbedoGrid(
134 radius = ref_radius,
135 frame = 'MOON_PA',
136 planet_name = 'Moon',
137 from_array = alb_grid,
138 axes = (Lon, Lat),
139 )
140
141# Define the temperature grid object
142Lon = np.load(temp_lon)
143Lat = np.load(temp_lat)
144
145TEMP = TemperatureGrid(
146 radius = ref_radius,
147 frame = 'GSE_MOON',
148 planet_name = 'Moon',
149 from_array = temp_grid,
150 axes = (Lon, Lat),
151 )
152
153# Set thermal properties
154moon.emissivity = 0.9
155moon.albedo = ALB
156moon.gridded_temperature = TEMP
157
158# Load the Look up table
159LUT = LookUpTable(lutfile)
160
161# Precomputation object
162prec = Precompute(epochs = epochs,)
163prec.precomputePlanetaryRadiation(lro, moon, LUT.moving_frames, correction='CN')
164prec.dump()
165
166# Create the albedo object
167albedo = Albedo(lro, LUT, moon, precomputation = prec, baseflux = base_flux,)
168
169# Create the thermal infrared object
170thermal_ir = Emissivity(lro, LUT, moon, precomputation = prec, baseflux = base_flux,)
171
172### ------------------------------------------------------------------------------------------------------- ###
173### COMPUTATIONS
174
175# Both the albedo and emissivity objects have a .compute() method.
176# This method returns the normalized fluxes, direction and albedo values
177# for each of the planet faces contributing to the computation
178# The general formula for computing the acceleration of an elementary face is:
179#
180# acc_i = L * albedo_value/mass * norm_flux
181#
182# where L is the normalized optical response of the spacecraft which can be computed
183# with raytracing setting a unitary radiance of the impacting rays
184# Due to the very high number of rays involved in these computations
185# it is mandatory to precompute L in the form of a lookup table.
186# Here we import the lookup table for LRO which can be computed using the example script
187# 'compute_lut.py'
188
189# Parallel computations
190alb_accel = albedo.compute(epochs, n_cores = n_cores)[0] * 1e3
191ir_accel = thermal_ir.compute(epochs, n_cores = n_cores)[0] * 1e3
192
193# Always unload the SPICE kernels
194sp.unload(METAKR)
195
196### ... Elapsed time
197toc = timeit.default_timer()
198time_min = int(floor((toc-tic)/60))
199time_sec = int(mod((toc-tic), 60))
200print("")
201print("\t Elapsed time: %d min, %d sec" %(time_min, time_sec))
202print("")
203
204### ------------------------------------------------------------------------------------------------------- ###
205### PLOT
206
207epochs = [float( epc - epc_et0 )/3600 for epc in epochs]
208
209# ALBEDO
210fig, ax = plt.subplots(3, 1, figsize=(14,8), sharex = True)
211
212ax[0].plot(epochs, alb_accel[:,0], linewidth = 2, color = "tab:blue")
213ax[0].set_ylabel('X [m/s^2]')
214ax[1].plot(epochs, alb_accel[:,1], linewidth = 2, color = "tab:blue")
215ax[1].set_ylabel('Y [m/s^2]')
216ax[2].plot(epochs, alb_accel[:,2], linewidth = 2, color = "tab:blue")
217ax[2].set_ylabel('Z [m/s^2]')
218ax[2].set_xlabel('Hours from CA')
219fig.suptitle('Albedo in S/C body frame')
220plt.tight_layout()
221
222# IR
223fig, ax = plt.subplots(3, 1, figsize=(14,8), sharex = True)
224
225ax[0].plot(epochs, ir_accel[:,0], linewidth = 2, color = "tab:blue")
226ax[0].set_ylabel('X [m/s^2]')
227ax[1].plot(epochs, ir_accel[:,1], linewidth = 2, color = "tab:blue")
228ax[1].set_ylabel('Y [m/s^2]')
229ax[2].plot(epochs, ir_accel[:,2], linewidth = 2, color = "tab:blue")
230ax[2].set_ylabel('Z [m/s^2]')
231ax[2].set_xlabel('Hours from CA')
232fig.suptitle('IR in S/C body frame')
233plt.tight_layout()
234
235plt.show()
236
237### ------------------------------------------------------------------------------------------------------- ###