C++ Module: albedo

Executive Summary

This document describes how albedo is modeled in the Basilisk software. The purpose of this module is to calculate albedo value at the given instrument locations.

The C++ Module: albedo module is a subclass of C++ Module: planetRadiationBase, which provides the core infrastructure for planetary radiation modeling including patch grid generation, IR flux computation, and eclipse handling.

Message Connection Descriptions

The following table lists all the module input and output messages. The module msg variable names are set by the user from python. The msg type contains a link to the message structure definition, while the description provides information on what this message is used for.

albedo module input and output messages

Module I/O Messages

Msg Variable Name

Msg Type

Description

spacecraftStateInMsg

SCStatesMsgPayload

input message with thruster commands.

sunPositionInMsg

SpicePlanetStateMsgPayload

Sun input message.

planetInMsgs

SpicePlanetStateMsgPayload

vector of planet input messages are set by using either the addPlanetandAlbedoAverageModel() method and addPlanetandAlbedoDataModel() method.

albOutMsgs

AlbedoMsgPayload

vector of albedo output message names. The order matches the order with which instruments are added.

Detailed Module Description

The albedo module is responsible for determining whether or not an incremental area on a planet is sunlit and within the field of view of an instrument and if so, what is the albedo value at the instrument’s location. The module uses the configuration of the instruments, and states of the sun, spacecraft and planets. It calculates the albedo value at the instrument locations as ratio [-] and albedo flux \([W/m^2]\). Albedo ratio of \(0.0\) represents no albedo, and \(1.0\) represents an equal flux with solar flux at the location.

To determine the states of the bodies in question, messages are passed into the code. For the spacecraft, Cartesian vectors provide the position and velocity of its center of mass. For the sun and planets their ephemeris messages are read in. If given multiple instruments, the code iterates through the sensor list and creates an array consisting of albedo values for each instrument. If given multiple planets, the code iterates through the planet list but sum all the values at the end and creates only one value for each instrument.

Determining States

The initial step in the module is to obtain the albedo model, then compute the albedo. In computing albedo, the celestial body state data are obtained and transformed into usable terms. The relationships shown below remove the dependency on relating position to the inertial frame \(N\) and instead utilize the spacecraft body \(B\), instrument body \(I\), sun \(S\), planet \(P\), and incremental area \(dA\) on the planet.

(1)\[\mathbf{r}_{BP} = \mathbf{r}_{BN} - \mathbf{r}_{PN}\]
(2)\[\mathbf{r}_{IB} = \mathbf{r}_{IN} - \mathbf{r}_{BN}\]
(3)\[\mathbf{r}_{IP} = \mathbf{r}_{IB} + \mathbf{r}_{BP}\]
(4)\[\mathbf{r}_{SP} = \mathbf{r}_{SN} - \mathbf{r}_{PN}\]
(5)\[\mathbf{r}_{IdA} = \mathbf{r}_{IP} - \mathbf{r}_{dAP}\]
(6)\[\mathbf{r}_{SdA} = \mathbf{r}_{SP} - \mathbf{r}_{dAP}\]

The previous two equations provide the sun’s and instrument’s position with respect to the incremental area using Eq. (1) - (4) and \(\mathbf{r}_{dAP}\), which is transformed from latitude and longitude of the grid points.

Sunlit Field of View Area

In determining the illuminated area within the instrument’s fov, \(f_1\), \(f_2\) and \(f_3\) are computed as shown below,

(7)\[f_1 = \frac{\mathbf{r}_{dAP}}{| \mathbf{r}_{dAP}|} \cdot \frac{\mathbf{r}_{SdA}}{| \mathbf{r}_{SdA}|}\]
(8)\[f_2 = \frac{\mathbf{r}_{dAP}}{| \mathbf{r}_{dAP}|} \cdot \frac{\mathbf{r}_{IdA}}{| \mathbf{r}_{IdA}|}\]
(9)\[f_3 = \hat{n}_N \cdot \frac{-\mathbf{r}_{IdA}}{| \mathbf{r}_{IdA}|}\]

Here \(\hat{n}_N\) indicates the unit normal vector of the instrument in inertial frame. \(f_1 > 0\) presents the sunlit \(f_2 > 0\) presents the instrument’s maximum fov, \(f_3 > \cos(fov)\) presents the instrument’s specified fov.

Albedo module needs three variables related to instrument’s configuration which are instrument’s misalignment vector with respect to spacecraft’s body frame (\(r_{{IB}_B}\)), unit normal vector of the instrument in spacecraft body frame (\(\hat{n}_B\)), and instrument’s field of view half angle in radian (\(fov\)). These variables can be added to the module using addInstrumentConfig() method. First term for the method is the instrument name. The rest of the terms can be set using the instConfig_t class or variable by variable respectively as: \(fov\), \(\hat{n}_B\), and \(r_{{IB}_B}\).

In the module, for planets that have polar radius, \(RP_{planet}\) and equatorial radius, \(REQ_{planet}\) defined, authalic radius is calculated. By doing this, the sphere is having the same surface area with the reference ellipsoid. If the polar radius is not defined, module uses the equatorial radius.

Albedo Value

Albedo flux ratio can be calculated as,

(10)\[\text{albedoAtInstrument} = ALB \frac{f_1 \cdot f_2 \cdot f_3 \cdot d_{Area}}{\pi \cdot |\mathbf{r}_{IdA}|^2}\]

where \(d_{Area}\) is the area of the incremental area, \(ALB\) is the albedo coefficient. There are albedo models based on an average albedo value and albedo data. The existing data files are placed under Basilisk/supportData/AlbedoData as .csv file format consisting \(ALB\) matrix. The number of rows represent the \(numLat\), number of latitude (between -90 to 90 deg) and columns represent the \(numLon\), number of longitude (between -180 to 180 deg).

The Earth’s albedo data is obtained from CERES instrument as .nd format and converted to .csv format for consistency with 1x1, 5x5, and 10x10 degree resolutions under clear-sky and all-sky conditions.

The Mars’ albedo data is obtained from TES instrument as VICAR format and converted to .csv format for consistency with 1x1, 5x5, and 10x10 degree resolutions.

illuminationFactorAtdA is optional to be calculated when getEclipseCase() returns True or can be assigned directly by the user when getEclipseCase() returns False. It is used as a multiplication term in Eq. (10), if defined. Therefore, when using albedo output on an instrument, it should be used after the illumination factor multiplication of the instrument, if exists.

A limit can be set in order not to compute the albedo for planets too far by \(altitudeRateLimit\) which is the limit for the rate of the instrument’s altitude to the planet’s radius.

Module Assumptions and Limitations

  • Albedo Model: The albedo models based on average value or specified data can be used.

  • Planet Shape: The module uses approximated authalic sphere which has the same surface area with the reference ellipsoid.

  • Planet Radius: The module have a list of planets with specified radius.

User Guide

This section outlines the steps needed to add Albedo module to a sim. First, the albedo module should be imported:

from Basilisk.simulation import albedo
albModule = albedo.Albedo()
albModule.ModelTag = "Albedo_module"

The instruments’ configuration must be added by using,

instConfig = albedo.instConfig_t()
instConfig.fov
instConfig.nHat_B
instConfig.r_IB_B
albModule.addInstrumentConfig(instConfig)

or by using,

albModule.addInstrumentConfig(fov, nHat_B, r_IB_B)

In the first case, if the variables are not defined for some reason and they are empty; then, default values are going to be used as \(fov = 90.\) deg, \(\hat{n}_B = [ 1.0, 0.0, 0.0 ]\), \(r_{{IB}_B} = [ 0.0, 0.0, 0.0 ]\). The default values can be defined by the user as well. Both functions for the instrument configuration has the ability to do a sanity check for \(fov\) being positive and \(\hat{n}_B\) not having all zero elements. Also, \(\hat{n}_B\) is always normalized. Then, the planet and albedo model function must be added.

There are three options based on the albedo model to be used. In all cases the planet input message is provided as an argument and the method makes the albedo modue subscribe to this message. For ALBEDO_AVG case,

albModule.addPlanetandAlbedoAverageModel(planetMsg)

where albedo average value is calculated automatically based on the given planet, and

albModule.addPlanetandAlbedoAverageModel(planetMsg, ALB_avg, numLat, numLon)

where the user can set the albedo average value. Number of latitude/longitude can also be specified or set to a negative, or zero, value to let default values being used instead (defaultNumLat = 180 and defaultNumLon = 360). The default values can be changed by the user as well. For ALBEDO_DATA case,

albModule.addPlanetandAlbedoDataModel(planetMsg, dataPath, fileName)

where the user can define the data path and file name for the albedo data to be used. The model can be added to a task like other simModels.

unitTestSim.AddModelToTask(simTaskName, albModule)

Warning

shadowFactor is deprecated in favor of illuminationFactor. In C/C++ both names work until Dec 31, 2026. After that date, code must use illuminationFactor exclusively.

Warning

eclipseCase will become a protected member, prefer using getter/setter. In C/C++ both names work until May 1, 2027. After that date, code must use getter/setter exclusively.


Typedefs

typedef class Config instConfig_t

albedo instrument configuration class

class Config
#include <albedo.h>

albedo instrument configuration class

Public Functions

Config()

Config constructor

~Config()

Config destructor

Public Members

double fov

[rad] instrument’s field of view half angle

Eigen::Vector3d nHat_B

[-] unit normal of the instrument (spacecraft body)

Eigen::Vector3d r_IB_B

[m] instrument’s misalignment wrt spacecraft’s body frame

class Albedo : public PlanetRadiationBase
#include <albedo.h>

albedo class

Note: spacecraftStateInMsg and sunPositionInMsg are provided by PlanetRadiationBase. The inherited planetInMsg is unused in this class; planets must be registered via addPlanetandAlbedoAverageModel() or addPlanetandAlbedoDataModel().

Public Functions

Albedo()

Albedo module constructor

~Albedo()

Albedo module destructor

void addInstrumentConfig(instConfig_t configMsg)

adds instrument configuration (overloaded function)

Adds the instrument configuration and automatically creates an output message name (overloaded function)

void addInstrumentConfig(double fov, Eigen::Vector3d nHat_B, Eigen::Vector3d r_IB_B)

adds instrument configuration (overloaded function)

Adds the instrument configuration and automatically creates an output message name (overloaded function)

void addPlanetandAlbedoAverageModel(Message<SpicePlanetStateMsgPayload> *msg)

This method adds planet msg and albedo average model name (overloaded function).

Sentinel values in the PlanetGrid cfg signal resolvePlanetEntry() to look up the value from the planet name: REQ_m = -1 → always filled from planet name RP_m = -1 → always filled from planet name albedoAvg = -1 → fill from planet name (ALBEDO_AVG_IMPLICIT case) nLat = 0, nLon = 0 → use base-class defaultNumLat/defaultNumLon

This method subscribes to the planet msg and sets the albedo average model (overloaded function)

void addPlanetandAlbedoAverageModel(Message<SpicePlanetStateMsgPayload> *msg, double ALB_avg, int numLat, int numLon)

This method adds planet name and albedo average model name (overloaded function).

This method subscribes to the planet msg and sets the albedo average model (overloaded function)

void addPlanetandAlbedoDataModel(Message<SpicePlanetStateMsgPayload> *msg, std::string dataPath, std::string fileName)

This method adds planet name and albedo data model.

This method subscribes to the planet msg and sets the albedo data model

double getAlbedoAverage(std::string planetSpiceName)

gets the average albedo value of the specified planet

This method gets the average albedo value of the planet

Returns:

double

Public Members

std::vector<Message<AlbedoMsgPayload>*> albOutMsgs

output messages (one per instrument)

Eigen::Vector3d nHat_B_default = {1.0, 0.0, 0.0}

[-] default instrument unit normal (body frame)

double fov_default = 90.0 * D2R

[rad] default instrument FOV half angle

double illuminationFactorAtdA = 1.0

[-] constant illumination factor (m_eclipseCase=false)

double altitudeRateLimit = -1.0

[-] max altitude/radius ratio; <0 = disabled

Protected Functions

void resolvePlanetEntry(const SpicePlanetStateMsgPayload &planetMsg, PlanetGrid &grid, int idx) override

fills in planet-dependent grid parameters

Fills in planet-dependent grid parameters (REQ_m, RP_m, albedoAvg, nLat, nLon) from the planet name

void customReset(uint64_t nanos) override

validates instrument and planet configuration

This method resets the module; validates that instruments and planets have been configured

void onUpdateBegin(const SCStatesMsgPayload &scMsg, const SpicePlanetStateMsgPayload &sunMsg, uint64_t nanos) override

caches spacecraft/sun state; clears accumulators

This method reads the messages and caches spacecraft/sun state; clears per-instrument accumulators

void evaluatePlanet(const SCStatesMsgPayload &scMsg, const SpicePlanetStateMsgPayload &sunMsg, SpicePlanetStateMsgPayload planetMsg, PlanetEntry &entry, const double S_sun, uint64_t nanos) override

evaluates the albedo model for one planet

This method evaluates the albedo model for one planet: applies FOV and eclipse filtering per instrument

void onUpdateEnd(uint64_t nanos) override

writes the output albedo messages

This method writes the output albedo messages from accumulated per-instrument data

Private Functions

std::pair<double, double> lookupPlanetRadius(const std::string &name)

gets the planet’s equatorial and polar radii in meters

Planet’s equatorial radii and polar radii (if exists) in meters

Returns:

pair {REQ_m, RP_m}; RP_m < 0 for spherical bodies

Private Members

std::vector<Eigen::Vector3d> r_IB_Bs

[m] instrument misalignment, body frame

std::vector<Eigen::Vector3d> nHat_Bs

[-] instrument unit normal, body frame

std::vector<double> fovs

[rad] instrument FOV half angles

Eigen::Vector3d r_BN_N

[m] spacecraft position, N frame

Eigen::Matrix3d dcm_BN

[-] rotation matrix body to inertial

Eigen::Vector3d r_SN_N

[m] sun position, N frame

uint64_t currentNanos = 0

[ns] current simulation time

std::vector<double> albedoAtInstrumentMax

[-] albedo ratio sum over all visible patches

std::vector<double> albedoAtInstrument

[-] FOV-filtered albedo ratio sum

std::vector<double> SfluxAtInstrument

[W/m^2] solar flux at instrument (last-planet value)

std::vector<double> AfluxAtInstrumentMax

[-] Albedo flux ratio at instrument

std::vector<double> AfluxAtInstrument

[W/m^2] Albedo flux at instrument