LRO Albedo and Thermal IR (Simple)

This example demonstrates how to compute albedo and thermal infrared accelerations using the object-oriented interface of the pyRTX library.

In this simplified case, we use single uniform values for:

  • Planet emissivity

  • Planet albedo

  • Planet surface temperature

This provides a quick way to compute planetary radiation pressure effects without requiring detailed surface property maps.

lro_alb_ir_simple.py
  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 a single value for the emissivity, albedo and temperature for the planet.
  9
 10### ------------------------------------------------------------------------------------------------------- ###
 11### IMPORTS
 12
 13import numpy as np
 14import spiceypy as sp
 15import matplotlib.pyplot as plt
 16import timeit
 17
 18from pyRTX.classes.Spacecraft import Spacecraft
 19from pyRTX.classes.Planet import Planet
 20from pyRTX.classes.Radiation import Albedo, Emissivity
 21from pyRTX.classes.Precompute import Precompute
 22from pyRTX.classes.LookUpTable import LookUpTable
 23from pyRTX.core.analysis_utils import epochRange2
 24
 25from numpy import floor, mod
 26
 27### ------------------------------------------------------------------------------------------------------- ###
 28### INPUTS
 29
 30# NOTE: before running this script you should generate the input lutfile running the 
 31# example 'compute_lut.py' using type = 'accel'.
 32
 33ref_epc		=  "2010 may 10 09:25:00"
 34duration    =  10000  									  # seconds
 35timestep    =  100
 36METAKR      = '../example_data/LRO/metakernel_lro.tm'     # metakernel
 37obj_path    = '../example_data/LRO/'				      # folder with .obj files
 38lutfile     = 'luts/lro_accel_lut.nc'					  # lookup table file
 39base_flux   =  1361.5
 40ref_radius  =  1737.4
 41n_cores     =  10
 42sc_mass		=  2000  									  # can be a float, int or xarray [kg]
 43
 44### ------------------------------------------------------------------------------------------------------- ###
 45### OBJECTS DEFINITION
 46
 47# Time initialization
 48tic = timeit.default_timer()
 49
 50# Load the metakernel containing references to the necessary SPICE frames
 51sp.furnsh(METAKR)
 52
 53# Define epochs
 54epc_et0  = sp.str2et( ref_epc )
 55epc_et1  = epc_et0 + duration
 56epochs   = epochRange2(startEpoch = epc_et0, endEpoch = epc_et1, step = timestep)
 57 
 58# Define the Spacecraft Object (Refer to the class documentation for further details)
 59lro = Spacecraft( name = 'LRO',
 60                 
 61				  base_frame = 'LRO_SC_BUS', 					     # Name of the spacecraft body-fixed frame
 62      
 63                  mass = sc_mass,
 64      
 65				  spacecraft_model = {						         # Define a spacecraft model
 66                          
 67					'LRO_BUS': { 
 68							 'file' : obj_path + 'bus_rotated.obj',	 # .obj file of the spacecraft component
 69							 'frame_type': 'Spice',				     # type of frame (can be 'Spice' or 'UD'
 70							 'frame_name': 'LRO_SC_BUS',			 # Name of the frame
 71							 'center': [0.0,0.0,0.0],			     # Origin of the component
 72							 'diffuse': 0.1,				         # Diffuse reflect. coefficient
 73							 'specular': 0.3,				         # Specular reflect. coefficient
 74							 },
 75
 76					'LRO_SA': {	
 77							'file': obj_path + 'SA_recentred.obj',
 78							'frame_type': 'Spice',
 79							'frame_name': 'LRO_SA',
 80							'center': [-1,-1.1, -0.1],
 81							'diffuse': 0,
 82							'specular': 0.3,
 83							},
 84
 85
 86					'LRO_HGA': { 	
 87							'file': obj_path + 'HGA_recentred.obj',
 88							'frame_type': 'Spice',
 89							'frame_name': 'LRO_HGA',
 90							'center':[-0.99,    -0.3,  -3.1],
 91							'diffuse': 0.2,
 92							'specular': 0.1,
 93							},
 94					}
 95					)
 96
 97# Define the Moon object
 98moon = Planet(  fromFile = None,
 99                radius = ref_radius,
100                name = 'Moon',
101                bodyFrame = 'MOON_PA',
102                sunFixedFrame = 'GSE_MOON',
103                units = 'km',
104                subdivs = 5,
105                )
106
107# Set the albedo and emissivity values
108# Here we use dummy values and assume that
109# albedo and emissivity are constant over the whole planet
110# and set a different dayside and nightside temperature
111# pyRTX supports also gridded values for albedo and emissivity (see examples/lro_alb_ir_grid.py)
112
113moon.albedo = 0.3
114moon.emissivity  = 0.8 
115moon.dayside_temperature = 300
116moon.nightside_temperature = 200
117
118# Load the Look up table
119LUT  = LookUpTable(lutfile)
120
121# Precomputation object
122prec = Precompute(epochs = epochs,)
123prec.precomputePlanetaryRadiation(lro, moon, LUT.moving_frames, correction='CN')
124prec.dump()
125
126# Create the albedo object
127albedo = Albedo(lro, LUT, moon, precomputation  = prec, baseflux  = base_flux,)
128
129# Create the thermal infrared object
130thermal_ir = Emissivity(lro, LUT, moon, precomputation  = prec, baseflux  = base_flux,)
131
132### ------------------------------------------------------------------------------------------------------- ###
133### COMPUTATIONS
134
135# Both the albedo and emissivity objects have a .compute() method.
136# This method returns the normalized fluxes, direction and albedo values
137# for each of the planet faces contributing to the computation
138# The general formula for computing the acceleration of an elementary face is:
139#
140# acc_i = L * albedo_value/mass * norm_flux
141#
142# where L is the normalized optical response of the spacecraft which can be computed
143# with raytracing setting a unitary radiance of the impacting rays
144# Due to the very high number of rays involved in these computations
145# it is mandatory to precompute L in the form of a lookup table.
146# Here we import the lookup table for LRO which can be computed using the example script
147# 'compute_lut.py'
148
149# Parallel computations
150alb_accel = albedo.compute(epochs, n_cores = n_cores)[0] * 1e3 
151ir_accel  = thermal_ir.compute(epochs, n_cores = n_cores)[0] * 1e3 
152
153# Always unload the SPICE kernels
154sp.unload(METAKR)
155
156### ... Elapsed time
157toc = timeit.default_timer()
158time_min = int(floor((toc-tic)/60))
159time_sec = int(mod((toc-tic), 60))
160print("")
161print("\t Elapsed time: %d min, %d sec" %(time_min, time_sec))
162print("")
163
164### ------------------------------------------------------------------------------------------------------- ###
165### PLOT
166
167epochs  = [float( epc - epc_et0 )/3600 for epc in epochs]
168
169# ALBEDO 
170fig, ax = plt.subplots(3, 1, figsize=(14,8), sharex = True)
171
172ax[0].plot(epochs, alb_accel[:,0], linewidth = 2, color = "tab:blue")
173ax[0].set_ylabel('X [m/s^2]')
174ax[1].plot(epochs, alb_accel[:,1], linewidth = 2, color = "tab:blue")
175ax[1].set_ylabel('Y [m/s^2]')
176ax[2].plot(epochs, alb_accel[:,2], linewidth = 2, color = "tab:blue")
177ax[2].set_ylabel('Z [m/s^2]')
178ax[2].set_xlabel('Hours from CA')
179fig.suptitle('Albedo in S/C body frame')
180plt.tight_layout()
181
182# IR
183fig, ax = plt.subplots(3, 1, figsize=(14,8), sharex = True)
184
185ax[0].plot(epochs, ir_accel[:,0], linewidth = 2, color = "tab:blue")
186ax[0].set_ylabel('X [m/s^2]')
187ax[1].plot(epochs, ir_accel[:,1], linewidth = 2, color = "tab:blue")
188ax[1].set_ylabel('Y [m/s^2]')
189ax[2].plot(epochs, ir_accel[:,2], linewidth = 2, color = "tab:blue")
190ax[2].set_ylabel('Z [m/s^2]')
191ax[2].set_xlabel('Hours from CA')
192fig.suptitle('IR in S/C body frame')
193plt.tight_layout()
194
195plt.show()
196
197### ------------------------------------------------------------------------------------------------------- ###