C++ Module: planetRadiationBase
Executive Summary
Abstract Basilisk base class for modules that compute planetary radiation using the Knocke et al. (1988) [Knocke1988] patch model.
PlanetRadiationBase owns the full simulation loop:
message validation and reading (spacecraft, Sun, per-planet)
authalic planet radius computation from WGS-84/spheroid parameters
patch grid initialisation and distance-corrected solar flux
per-timestep call to
evaluatePlanet()for every registered planet
Derived classes implement a set of virtual hooks rather than overriding
Reset() / UpdateState() directly.
User Guide
Register planets by calling addPlanetEntry(msg, gridConfig) from a derived
addPlanet*() method. The grid config may contain sentinel values (e.g.
REQ_m = -1, albedoAvg = -1, nLat = 0) that
resolvePlanetEntry() fills in during Reset().
|
pure virtual |
Called once per planet per timestep. Call |
|
optional |
Called in |
|
optional |
Called before the planet loop. Clear per-timestep accumulators or cache spacecraft attitude. Default: no-op. |
|
optional |
Called after the planet loop. Write output messages. Default: no-op. |
|
optional |
Called at the end of |
|
optional |
Determine if a patch is in eclipse (umbra/penumbra). Default implementation
applies the Knocke penumbra model when |
Single-planet backward-compatibility path
If no planets are registered via addPlanetEntry() at Reset() time
(i.e. planets_ is empty), the base class automatically builds one entry
from the public scalar config members (nLat, nLon, irFluxMean,
albedoAvg, useAlbedoData, REQ_m, RP_m) and the public
planetInMsg functor. This means EarthRadiationModel users need no
code changes.
Current derived classes
C++ Module: earthRadiationModel — ERP flux and flux-weighted direction vectors. Uses the backward-compat single-planet path; overrides
evaluatePlanet(),onUpdateBegin(), andonUpdateEnd().C++ Module: albedo — per-instrument albedo ratio and flux for multiple planets and multiple instruments. Registers planets via
addPlanetEntry()insideaddPlanetandAlbedoAverageModel()/addPlanetandAlbedoDataModel(); overrides all five hooks.
Known limitations
EarthRadiationModel has no Python-facing addPlanetandAlbedoAverageModel()-style
API. Multi-planet support is wired in the base class (the planets_ loop runs
for all registered entries), but the user must call addPlanetEntry() directly
from C++ or a future Python-facing wrapper is needed.
Message Connection Descriptions
The following messages are defined in PlanetRadiationBase and inherited by all
derived classes.
Msg Variable Name |
Msg Type |
Description |
|---|---|---|
spacecraftStateInMsg |
spacecraft position and attitude in inertial frame N |
|
planetInMsg |
single-planet state (ERM backward-compat path only) |
|
sunPositionInMsg |
Sun position used for solar flux distance correction and patch illumination |
References
Knocke, P. C., Ries, J. C., and Tapley, B. D., Earth Radiation Pressure Effects on Satellites, AIAA/AAS Astrodynamics Conference, 1988. DOI: 10.2514/6.1988-4292
-
struct PatchResult
- #include <planetRadiationBase.h>
Per-patch flux and geometry returned by PlanetGrid::computePatches().
Public Members
-
double albedoFlux
[W/m^2] albedo flux from this patch; 0 if not illuminated
-
double irFlux
[W/m^2] IR flux from this patch; > 0 for all visible patches
-
double rHat_IdA_N[3]
[-] unit vector patch centroid -> satellite, inertial N
-
double solidAngle
[sr] solid angle of patch as seen from satellite
-
double r_dAP_N[3]
[m] patch centroid position relative to planet centre, N frame
-
double albedoFlux
-
struct PlanetGrid
- #include <planetRadiationBase.h>
Initialised patch grid for one planet.
Configure all public members, then call initialize(). After that, call computePatches() every timestep.
Public Functions
-
void initialize(BSKLogger bskLogger)
Builds the patch grid and loads optional CSV albedo data. Call once after setting all configuration members.
-
std::vector<PatchResult> computePatches(const double r_sat_N[3], const double r_sun_N[3], const double r_planet_N[3], double J20002Pfix[3][3], double R_planet, const double S_sun) const
Computes per-patch radiation fluxes for one timestep.
Only patches with cos_sat > 0 (visible from r_sat_N) are returned. Eclipse shadowing is not applied here; use isPatchEclipsed() post-kernel.
- Parameters:
r_sat_N – [m] satellite (or instrument) position, inertial N
r_sun_N – [m] Sun position, inertial N
r_planet_N – [m] planet centre position, inertial N
J20002Pfix – [-] SpicePlanetStateMsgPayload rotation matrix (row-major)
R_planet – [m] authalic planet radius
S_sun – [W/m^2] solar flux at planet, distance-corrected
Public Members
-
int nLat = 180
[-] latitude zones
-
int nLon = 360
[-] longitude zones per band
-
double irFluxMean = 237.0
[W/m^2] mean
-
double albedoAvg = 0.30
[-] Bond albedo (albedoDataPath=””)
-
bool useAlbedoData = false
true = load CSV albedo data
-
std::string albedoDataPath = ""
path to supportData/AlbedoData/
-
std::string albedoDataFile = ""
CSV filename.
-
double REQ_m = -1.0
[m] equatorial radius; -1 = must be set before initialize()
-
double RP_m = -1.0
[m] polar radius; -1 = spherical (REQ used)
Private Functions
-
std::vector<std::vector<double>> loadAlbedoGridFromCsv(BSKLogger bskLogger) const
Private Members
-
std::vector<std::array<double, 3>> patchDirs_P
[-] patch centroid directions in planet frame
-
std::vector<double> normAreas
[sr] normalised patch areas
-
std::vector<double> patchAlbedo
[-] per-patch albedo values
-
double t_aut[3] = {0.0, 0.0, 0.0}
authalic latitude correction coefficients
-
bool hasAuthalic = false
true when oblate spheroid correction is active
-
double latDiff = 0.0
[rad] latitude spacing between grid points
-
double lonDiff = 0.0
[rad] longitude spacing between grid points
-
bool initialized = false
true after initialize() completes successfully
-
void initialize(BSKLogger bskLogger)
-
struct PlanetEntry
- #include <planetRadiationBase.h>
One registered planet: message functor, initialised grid, authalic radius.
Public Members
-
ReadFunctor<SpicePlanetStateMsgPayload> planetMsg
per-planet input message
-
PlanetGrid grid
initialised patch grid
-
double R_planet = 0.0
[m] authalic radius, set in Reset()
-
ReadFunctor<SpicePlanetStateMsgPayload> planetMsg
-
class PlanetRadiationBase : public SysModel
- #include <planetRadiationBase.h>
Abstract base class for Basilisk planetary-radiation environment modules.
PlanetRadiationBase reads spacecraft and sun state messages, iterates over registered planets, and dispatches per-planet evaluation to the derived class via evaluatePlanet(). Concrete subclasses implement evaluatePlanet() and optionally override the other virtual hooks.
Public Functions
-
PlanetRadiationBase() = default
-
~PlanetRadiationBase() override = default
-
void Reset(uint64_t CurrentSimNanos) override
Validates messages, builds planet entries, and initialises patch grids.
- Parameters:
CurrentSimNanos – Current simulation time (ns).
-
void UpdateState(uint64_t CurrentSimNanos) override
Reads inputs, dispatches per-planet evaluation, and writes outputs.
- Parameters:
CurrentSimNanos – Current simulation time (ns).
-
void addPlanetEntry(ReadFunctor<SpicePlanetStateMsgPayload> msg, const PlanetGrid &gridConfig)
Registers a planet. Call from derived addPlanet*() methods.
gridConfig may have sentinel values (e.g. REQ_m=-1, albedoAvg=-1, nLat=0) that resolvePlanetEntry() fills in during Reset().
- Parameters:
msg – Planet state message functor.
gridConfig – Partially-configured patch grid.
-
const bool getEclipseCase()
getter for the eclipseCase flag
-
void setEclipseCase(bool value)
setter for the eclipseCase flag
Public Members
-
int defaultNumLat = 180
[-] default number of latitude grid points
-
int defaultNumLon = 360
[-] default number of longitude grid points
-
double irFluxMean = 237.0
[W/m^2] mean outgoing IR radiation
-
double albedoAvg = 0.30
[-] Bond albedo, used when albedoDataPath is left empty
-
std::string albedoDataPath = ""
path to supportData/AlbedoData/
-
std::string albedoDataFile = ""
CSV filename.
-
double REQ_m = -1.0
[m] equatorial radius; -1 (default) = use WGS-84 Earth value
-
double RP_m = -1.0
[m] polar radius; -1 (default) = use WGS-84 Earth value
-
ReadFunctor<SCStatesMsgPayload> spacecraftStateInMsg
spacecraft position input message
-
ReadFunctor<SpicePlanetStateMsgPayload> planetInMsg
single-planet convenience (ERM path)
-
ReadFunctor<SpicePlanetStateMsgPayload> sunPositionInMsg
Sun position input message.
-
BSKLogger bskLogger
BSK logging.
-
bool eclipseCase = false
true = compute eclipse shadowing at each patch. Deprecated, will be removed after May 1st 2027.
-
std::vector<PlanetEntry> planets
registered planets, built in Reset() or addPlanetEntry()
Protected Functions
-
virtual void evaluatePlanet(const SCStatesMsgPayload &scMsg, const SpicePlanetStateMsgPayload &sunMsg, SpicePlanetStateMsgPayload planetMsg, PlanetEntry &entry, const double S_sun, uint64_t nanos) = 0
Called once per planet per timestep. PURE VIRTUAL.
- Parameters:
scMsg – spacecraft state
sunMsg – sun state
planetMsg – planet state (non-const copy, safe to pass J20002Pfix onward)
entry – planet entry with initialised grid and authalic radius
S_sun – [W/m^2] distance-corrected solar flux at planet
nanos – simulation time
-
inline virtual void resolvePlanetEntry(const SpicePlanetStateMsgPayload &planetMsg, PlanetGrid &grid, int idx)
Called in Reset() for each entry after reading the planet message.
Override to fill in grid members that depend on the planet name, e.g. REQ_m, RP_m, albedoAvg, nLat, nLon. Default: no-op (ERM relies on scalar config set before Reset()).
- Parameters:
planetMsg – Planet state message.
grid – PlanetGrid to fill in.
idx – Index of this planet entry.
-
inline virtual void onUpdateBegin(const SCStatesMsgPayload &scMsg, const SpicePlanetStateMsgPayload &sunMsg, uint64_t nanos)
Called at the start of UpdateState(), before the planet loop. Override to clear per-timestep accumulators or cache the spacecraft attitude. Default: no-op.
- Parameters:
scMsg – Spacecraft state message.
sunMsg – Sun state message.
nanos – Current simulation time (ns).
-
inline virtual void onUpdateEnd(uint64_t nanos)
Called after the planet loop in UpdateState(). Override to write output messages. Default: no-op.
- Parameters:
nanos – Current simulation time (ns).
-
inline virtual void customReset(uint64_t nanos)
Called at the end of Reset(). Override for additional validation (e.g. instrument configuration). Default: no-op.
- Parameters:
nanos – Current simulation time (ns).
-
inline virtual bool allowSinglePlanetFallback() const
Return true to enable the single-planet backward-compat fallback in Reset(): when planets is empty and planetInMsg is linked, a default PlanetGrid is built from the scalar config fields. Default: false.
-
virtual bool isPatchEclipsed(const PatchResult &patch, const SCStatesMsgPayload &scMsg, const SpicePlanetStateMsgPayload &sunMsg)
Determines if a patch is in eclipse (umbra/penumbra). Returns true if the patch should contribute zero albedo flux.
- Parameters:
patch – Patch result to evaluate.
scMsg – Spacecraft state message.
sunMsg – Sun state message.
Protected Attributes
-
Eigen::Vector3d currentR_SP_N = Eigen::Vector3d::Zero()
[m] Sun position relative to current planet centre
-
double currentPlanetRadius = 0.0
[m] current planet authalic radius
-
bool useAlbedoData = false
true = load CSV albedo data
-
bool m_eclipseCase = false
true = compute eclipse shadowing at each patch
Protected Static Functions
-
static double computeIlluminationAtdA(double Rplanet, Eigen::Vector3d r_dAP_N, Eigen::Vector3d r_SP_N)
Computes the umbra/penumbra illumination factor [0, 1] for a surface element at r_dAP_N (relative to planet centre). Returns 0 in full umbra, 1 in full sunlight, smooth intermediate in penumbra.
Compute the fraction of the solar disk visible from a surface patch. Returns: 1.0 -> fully illuminated 0.0 -> full umbra 0.0 < value < 1.0 -> penumbra
The calculation follows the apparent angular radii of the Sun and planet as seen from the patch and computes the overlap area of the two disks.
- Parameters:
Rplanet – [m] planet radius
r_dAP_N – [m] patch position relative to planet centre
r_SP_N – [m] Sun position relative to planet centre
-
PlanetRadiationBase() = default