ModelSpace
All Classes Namespaces Functions Variables Enumerations Pages
sphericalharmonicsutils.h
1/******************************************************************************
2* Copyright (c) ATTX LLC 2024. All Rights Reserved.
3*
4* This software and associated documentation (the "Software") are the
5* proprietary and confidential information of ATTX, LLC. The Software is
6* furnished under a license agreement between ATTX and the user organization
7* and may be used or copied only in accordance with the terms of the agreement.
8* Refer to 'license/attx_license.adoc' for standard license terms.
9*
10* EXPORT CONTROL NOTICE: THIS SOFTWARE MAY INCLUDE CONTENT CONTROLLED UNDER THE
11* INTERNATIONAL TRAFFIC IN ARMS REGULATIONS (ITAR) OR THE EXPORT ADMINISTRATION
12* REGULATIONS (EAR99). No part of the Software may be used, reproduced, or
13* transmitted in any form or by any means, for any purpose, without the express
14* written permission of ATTX, LLC.
15******************************************************************************/
16/* ADDITIONAL COPYRIGHT NOTICE
17 * The header accompanying the original file from which these functions are derived
18 is as follows: */
19/* This file is distributed with 42, */
20/* the (mostly harmless) spacecraft dynamics simulation */
21/* created by Eric Stoneking of NASA Goddard Space Flight Center */
22/* Copyright 2010 United States Government */
23/* as represented by the Administrator */
24/* of the National Aeronautics and Space Administration. */
25/* No copyright is claimed in the United States */
26/* under Title 17, U.S. Code. */
27/* All Other Rights Reserved. */
28#ifndef UTILS_SPHERICALHARMONICSUTILS_H
29#define UTILS_SPHERICALHARMONICSUTILS_H
30
31#include <string>
32
33#include "core/CartesianVector.hpp"
34#include "core/macros.h"
35
36namespace modelspace {
37 // Constants for spherical harmonics
38 const int NMAX_SPHERICAL_HARMONICS = 18;
39
40 /**
41 * @brief Reads gravitational coefficients from a file and stores them in provided matrices.
42 *
43 * This function opens a specified file containing spherical harmonic gravitational
44 * coefficients and populates the provided matrices `cbar` and `sbar` with these values
45 * up to degree and order `NMAX_SPHERICAL_HARMONICS`. The function checks for the file's
46 * existence before attempting to read from it. This code is derived from the open-source
47 * NASA 42 simulation.
48 *
49 * @param[in] filename Path to the file containing the gravitational coefficients.
50 * @param[out] cbar Matrix to store cosine harmonic coefficients.
51 * @param[out] sbar Matrix to store sine harmonic coefficients.
52 * @param[out] n Max degree loaded from egm96 file.
53 * @param[out] m Max order loaded from egm96 file.
54 *
55 * @return Status code indicating success or failure.
56 * @retval NO_ERROR If the file is successfully read and coefficients are loaded.
57 * @retval ERROR_FILE_NOT_FOUND If the specified file does not exist.
58 *
59 * @note This function was adapted from open-source code in the NASA 42 simulation.
60 *
61 * @warning The function assumes that the file format matches the expected format with
62 * columns for degree, order, and the respective coefficients. Incorrect file
63 * formatting may lead to unexpected results.
64 */
65 int readGravityCoefficientsFile(const std::string &filename,
66 double cbar[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1],
67 double sbar[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1],
68 long int &n, long int &m);
69
70 /**
71 * @brief Applies Neumann normalization to spherical harmonic coefficients.
72 *
73 * This function transforms spherical harmonic coefficients from their original format
74 * (typically EGM96 style normalization) to Neumann normalization. The input matrices `cbar`
75 * and `sbar` contain the cosine and sine coefficients, which are normalized and stored
76 * in the output matrices `c` and `s` for use in gravitational calculations.
77 *
78 * The Neumann normalization ensures that the coefficients are scaled to a standard
79 * that provides orthonormality and consistency across spherical harmonics terms.
80 *
81 * @param[in] cbar Matrix of original cosine coefficients.
82 * @param[in] sbar Matrix of original sine coefficients.
83 * @param[out] c Matrix to store Neumann-normalized cosine coefficients.
84 * @param[out] s Matrix to store Neumann-normalized sine coefficients.
85 *
86 * @note This function was adapted from open-source code in the NASA 42 simulation.
87 *
88 * @warning This function assumes the matrices are dimensioned to support up to degree
89 * and order `NMAX_SPHERICAL_HARMONICS`. Applying the normalization to
90 * coefficients beyond this range may result in undefined behavior.
91 */
92 void neumannNormalization(double cbar[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1],
93 double sbar[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1],
94 double c[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1],
95 double s[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1]);
96
97 /**
98 * @brief Computes the gradient of the gravitational potential using spherical harmonics.
99 *
100 * This function calculates the radial, latitudinal, and longitudinal components of the
101 * gravitational potential gradient up to degree N and order M using normalized Legendre
102 * coefficients. The result is stored in gradV. The maximum size is 18x18
103 *
104 * @param[in] N Maximum degree of spherical harmonics.
105 * @param[in] M Maximum order of spherical harmonics.
106 * @param[in] r Radial distance from planet's center.
107 * @param[in] phi Longitude (angle from prime meridian in radians).
108 * @param[in] theta Colatitude (angle from the north pole in radians).
109 * @param[in] Re Reference radius, typically planet's radius in meters.
110 * @param[in] K Scaling constant, typically planet's gravitational parameter divided by Re.
111 * @param[in] C Normalized Legendre coefficients (cosine terms).
112 * @param[in] S Normalized Legendre coefficients (sine terms).
113 * @param[out] grad_v Gradient of the gravitational potential [radial, latitudinal, longitudinal].
114 * @return Error code corresponding to success/failure
115 *
116 * @note This function was adapted from open-source code in the NASA 42 simulation.
117 */
118 int sphericalHarmonics(long N, long M, double r, double phi,
119 double theta, double Re, double K,
120 double C[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1],
121 double S[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1],
122 CartesianVector3D &grad_v);
123
124 /**
125 * @brief Calculates the associated Legendre functions and their derivatives.
126 *
127 * Computes the associated Legendre functions up to degree N and order M
128 * with Neumann normalization. This includes calculating scaled derivatives (sdP),
129 * which handle singularities by multiplying the derivative by sqrt(1 - x^2).
130 *
131 * @param[in] N Maximum degree of the Legendre functions.
132 * @param[in] M Maximum order of the Legendre functions.
133 * @param[in] x Argument for Legendre functions, typically cos(theta).
134 * @param[out] P Matrix of associated Legendre functions up to degree N and order M.
135 * @param[out] sd_p Matrix of scaled derivatives of Legendre functions.
136 * @return Error code corresponding to success/failure
137 *
138 * @note This function was adapted from open-source code in the NASA 42 simulation.
139 * @warning The derivatives sdP[n][1] are singular at x = ±1.
140 */
141 int legendre(long N, long M, double x,
142 double P[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1],
143 double sd_p[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1]);
144
145 /**
146 * @brief Computes the factorial of a non-negative integer.
147 *
148 * This helper function calculates the factorial of a given integer n (n!).
149 * It is used to normalize spherical harmonics coefficients.
150 *
151 * @param[in] n Non-negative integer for which to calculate the factorial.
152 * @return Factorial of the input integer.
153 *
154 * @note This function was adapted from open-source code in the NASA 42 simulation.
155 */
156 double fact(long n);
157
158 /**
159 * @brief Computes the double factorial of an odd integer.
160 *
161 * Calculates the product of all odd numbers up to and including the specified integer n.
162 * This is used in the Legendre function for terms involving odd factorials.
163 *
164 * @param[in] n Non-negative odd integer for which to calculate the double factorial.
165 * @return Double factorial of the input integer.
166 *
167 * @note This function was adapted from open-source code in the NASA 42 simulation.
168 */
169 double oddfact(long n);
170}
171
172#endif
#define CartesianVector3D
Definition macros.h:54
Class to propagate CR3BP dynamics in characteristic units.
Definition ConfigurationWriter.cpp:18
int sphericalHarmonics(long N, long M, double r, double phi, double theta, double Re, double K, double C[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1], double S[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1], clockwerk::CartesianVector< double, 3 > &grad_v)
Computes the gradient of the gravitational potential using spherical harmonics.
Definition sphericalharmonicsutils.cpp:76
double fact(long n)
Computes the factorial of a non-negative integer.
Definition sphericalharmonicsutils.cpp:183
double oddfact(long n)
Computes the double factorial of an odd integer.
Definition sphericalharmonicsutils.cpp:192
void neumannNormalization(double cbar[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1], double sbar[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1], double c[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1], double s[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1])
Applies Neumann normalization to spherical harmonic coefficients.
Definition sphericalharmonicsutils.cpp:58
int readGravityCoefficientsFile(const std::string &filename, double cbar[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1], double sbar[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1], long int &n, long int &m)
Reads gravitational coefficients from a file and stores them in provided matrices.
Definition sphericalharmonicsutils.cpp:37
int legendre(long N, long M, double x, double P[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1], double sd_p[NMAX_SPHERICAL_HARMONICS+1][NMAX_SPHERICAL_HARMONICS+1])
Calculates the associated Legendre functions and their derivatives.
Definition sphericalharmonicsutils.cpp:139