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.

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.

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.

shadowFactorAtdA is optional to be calculated with eclipseCase being True or can be assigned directly by the user with eclipseCase 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 shadow 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 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)

Typedefs

typedef class Config instConfig_t

albedo instrument configuration class

class Config
#include <albedo.h>

albedo instrument configuration class

Public Functions

Config()
~Config()

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 SysModel
#include <albedo.h>

albedo class

Public Functions

Albedo()

Albedo module constructor

Returns

void

~Albedo()

Albedo module destructor

Returns

void

void UpdateState(uint64_t CurrentSimNanos)

updates the state

Read Messages, calculate albedo then write it out

Returns

void

void Reset(uint64_t CurrentSimNanos)

resets the module

This method resets the module

Returns

void

void addInstrumentConfig(instConfig_t configMsg)

adds instrument configuration (overloaded function)

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

Returns

void

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)

Returns

void

void addPlanetandAlbedoAverageModel(Message<SpicePlanetStateMsgPayload> *msg)

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

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

Returns

void

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)

Returns

void

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

Returns

void

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

vector of output messages for albedo data

ReadFunctor<SpicePlanetStateMsgPayload> sunPositionInMsg

input message name for sun data

ReadFunctor<SCStatesMsgPayload> spacecraftStateInMsg

input message name for spacecraft data

std::vector<ReadFunctor<SpicePlanetStateMsgPayload>> planetInMsgs

vector of planet data input data

BSKLogger bskLogger

BSK Logging

int defaultNumLat

[-] default number of latitude grid points

int defaultNumLon

[-] default number of longitude grid points

Eigen::Vector3d nHat_B_default

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

double fov_default

[rad] default value for instrument’s field of view half angle

bool eclipseCase

consider eclipse at dA, if true

double shadowFactorAtdA

[-] shadow factor at incremental area

double altitudeRateLimit

[-] rate limit of the instrument’s altitude to the planet’s radius for albedo calculations

Private Functions

void readMessages()

reads the inpt messages

This method reads the messages and saves the values to member attributes

Returns

void

void writeMessages(uint64_t CurrentSimNanos)

writes the outpus messages

This method writes the output albedo messages

Returns

void

void getPlanetRadius(std::string planetSpiceName)

gets the planet’s radius

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

Returns

void

void evaluateAlbedoModel(int idx)

evaluates the ALB model

This method evaluates the albedo model

Returns

void

void computeAlbedo(int idx, int instIdx, SpicePlanetStateMsgPayload planetMsg, bool AlbArray, double outData[])

computes the albedo at instrument’s location

This method calculates the albedo at instrument

Returns

void

double computeEclipseAtdA(double Rplanet, Eigen::Vector3d r_dAP_N, Eigen::Vector3d r_SP_N)

computes the shadow factor at dA

This method computes eclipse at the incremental area if eclipseCase is defined true

Returns

double

Private Members

std::vector<std::string> dataPaths

string with the path to the ALB coefficient folder

std::vector<std::string> fileNames

file names containing the ALB coefficients

std::vector<std::string> modelNames

albedo model names

std::vector<double> ALB_avgs

[-] albedo average value vector for each planet defined

std::vector<double> REQ_planets
std::vector<double> RP_planets

[m] equatorial and polar radius of the planets

Eigen::MatrixXd albLon
Eigen::MatrixXd albLat

sunlit area seen by the instroment fov

double scLon
double scLat

[deg, deg] spaceccraft footprint

std::vector<int> numLats
std::vector<int> numLons

[-] vector of latitude and longitude number

double albedoAtInstrument

[-] total albedo at instrument location

double albedoAtInstrumentMax

[-] max total albedo at instrument location

double SfluxAtInstrument

[W/m2] solar flux at instrument’s position

double AfluxAtInstrumentMax

[W/m2] max albedo flux at instrument’s position

double AfluxAtInstrument

[W/m2] albedo flux at instrument’s position

std::vector<Eigen::Vector3d> r_IB_Bs

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

std::vector<Eigen::Vector3d> nHat_Bs

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

std::vector<double> fovs

[rad] vector containing instrument’s field of view half angle

std::vector<Eigen::Vector4d> albOutData

the vector that keeps the albedo output data

std::map<int, std::vector<double>> gdlat

[rad] geodetic latitude

std::map<int, std::vector<double>> gdlon

[rad] geodetic longitude

std::map<int, double> latDiff

[rad] latitude difference between grid points

std::map<int, double> lonDiff

[rad] longitude difference between grid points

std::map<int, std::vector<std::vector<double>>> ALB

[-] ALB coefficients

bool readFile

defines if there is a need for reading an albedo model file or not

std::vector<bool> albArray

defines if the albedo data is formatted as array or not

Eigen::Vector3d r_PN_N

[m] planet position (inertial)

Eigen::Vector3d r_SN_N

[m] sun position (inertial)

Eigen::MRPd sigma_BN

[-] Current spaceraft MRPs (inertial)

Eigen::Vector3d r_BN_N

[m] s/c position (inertial)

Eigen::Vector3d nHat_N

[-] Unit normal vector of the instrument (inertial)

Eigen::Vector3d rHat_PI_N

[-] direction vector from instrument to planet

std::vector<SpicePlanetStateMsgPayload> planetMsgData

vector of incoming planet message states

SpicePlanetStateMsgPayload sunMsgData

sun message data

SCStatesMsgPayload scStatesMsgData

spacecraft message data