/******************************************************** * ██████╗ ██████╗████████╗██╗ * ██╔════╝ ██╔════╝╚══██╔══╝██║ * ██║ ███╗██║ ██║ ██║ * ██║ ██║██║ ██║ ██║ * ╚██████╔╝╚██████╗ ██║ ███████╗ * ╚═════╝ ╚═════╝ ╚═╝ ╚══════╝ * Geophysical Computational Tools & Library (GCTL) * * Copyright (c) 2022 Yi Zhang (yizhang-geo@zju.edu.cn) * * GCTL is distributed under a dual licensing scheme. You can redistribute * it and/or modify it under the terms of the GNU Lesser General Public * License as published by the Free Software Foundation, either version 2 * of the License, or (at your option) any later version. You should have * received a copy of the GNU Lesser General Public License along with this * program. If not, see . * * If the terms and conditions of the LGPL v.2. would prevent you from using * the GCTL, please consider the option to obtain a commercial license for a * fee. These licenses are offered by the GCTL's original author. As a rule, * licenses are provided "as-is", unlimited in time for a one time fee. Please * send corresponding requests to: yizhang-geo@zju.edu.cn. Please do not forget * to include some description of your company and the realm of its activities. * Also add information on how to contact you by electronic and paper mail. ******************************************************/ #ifndef _GCTL_GM_DATA_H #define _GCTL_GM_DATA_H #include "gctl_potential_config.h" #ifdef GCTL_POTENTIAL_OPENMP #include "omp.h" #include #endif #include "gctl/core.h" #include "gctl/utility.h" #include "gctl/geometry.h" #include "gctl/maths.h" #include "gctl/algorithms.h" #include "gctl/io.h" namespace gctl { /** * @brief components of potential field */ enum gravitational_field_type_e { GravPot, ///< field potential Vx, ///< x-gradient of V (longitudinal-gradient in the spherical coordinates) Vy, ///< y-gradient of V (latitudinal-gradient in the spherical coordinates) Vz, ///< z-gradient of V (radial-gradient in the spherical coordinates) Txx, ///< x-gradient of Gx Txy, ///< y-gradient of Gx Txz, ///< z-gradient of Gx Tyx, ///< x-gradient of Gy Tyy, ///< y-gradient of Gy Tyz, ///< z-gradient of Gy Tzx, ///< x-gradient of Gz Tzy, ///< y-gradient of Gz Tzz, ///< z-gradient of Gz }; enum magnetic_field_type_e { MagPot, Hax, Hay, Za, Bx, By, Bz, Bxx, Bxy, Bxz, Byx, Byy, Byz, Bzx, Bzy, Bzz, DeltaT, ///< total magnetic strength DeltaTx, ///< x-gradient of total magnetic strength DeltaTy, ///< x-gradient of total magnetic strength DeltaTz, ///< x-gradient of total magnetic strength }; /** * @brief A structure holds all outputs of the standard IGRF or EMM program. * */ struct IGRF_para { int Year, Month, Day; char CSystem, AltiCode; double Altitude, Latitude, Longitude; double D, I, H_nT, X_nT, Y_nT, Z_nT, F_nT, dD_min, dI_min, dH_nT, dX_nT, dY_nT, dZ_nT, dF_nT; int D_deg, D_min, I_deg, I_min; }; /** * @brief A structure holds a block of spherical harmonic coefficients * */ struct shc_data { int ns, ms, ne, me; // start degree, start order, end degree, end order double time; // snapshot array Snm, Cnm; // Snm or gnm, Cnm or hnm double coeff_s(int n, int m){return Snm[n*(n + 1)/2 + m];} double coeff_c(int n, int m){return Cnm[n*(n + 1)/2 + m];} }; /** * @brief Calculate the transform matrix of the vector componments wrt. the spherical coordinates. * * @param op The observation point on a sphere * @return The transform matrix */ tensor transform_matrix(const point3ds &op); /** * @brief Transform localized magnetization vectors on a sphere to the absolute Cartesian coordinates. * * @param abs_b_ Output absolute vectors * @param geo_b_ Local magnetization vectors on a sphere * @param loc_s_ Spherical coordinates of the local vectors */ void geomag_local2Cartesian(array &abs_b_, const array &geo_b_, const array &loc_s_); /** * @brief Read output table of the standard IGRF or EMM program * * @param file Input filename * @param IGRFs Return IGRF parameters * @param head_record Lines of head records * @param ext_f File name extension */ void read_IGRF_table(std::string file, array &IGRFs, int head_record = 1, std::string ext_f = ".txt"); /** * @brief Spherical harmonic coefficients stored in the Swarm level-2 product format. This format * is mainly used for storing spherical harmonic models of the Earth's magnetic field. * * @param file File name * @param SHCs Returned coefficients started from 0/0 to N/M * @param spline_odr Spline order for time dependent model * @param sample_num Sample order for time dependent model */ void read_Swarm_shc(std::string file, array &SHCs, int &spline_odr, int &sample_num); /** * @brief Save spherical harmonic coefficients to the Swarm level-2 product format. * * @param file File name. * @param SHCs coefficients started from 0/0 * @param spline_odr Spline order for time dependent model * @param sample_num Sample order for time dependent model * @param info Information to be saved with the file */ void save_Swarm_shc(std::string file, const array &SHCs, int spline_odr = 1, int sample_num = 1, std::string info = "This file is saved in the Swarm level-2 product format."); /** * @brief Transform magnetic components data to delta_T anomalies data. * * @param Hax Hax component of the magnetic data * @param Hay Hay component of the magnetic data * @param Za Za component of the magnetic data * @param deltaT Output delta_T anomalies data * @param T0_inclina Inclination degree of the earth normal magnetic field. * @param T0_declina Declination degree of the earth normal magnetic field. */ void magnetic_components2deltaT(const _1d_array &Hax, const _1d_array &Hay, const _1d_array &Za, _1d_array &deltaT, double T0_inclina, double T0_declina); /** * @brief Transform magnetic components data to delta_T anomalies data. * * @param Mag_components Magnetic components of the magnetic data * @param deltaT Output delta_T anomalies data * @param T0_inclina Inclination degree of the earth normal magnetic field. * @param T0_declina Declination degree of the earth normal magnetic field. */ void magnetic_components2deltaT(const array &Mag_components, _1d_array &deltaT, double T0_inclina, double T0_declina); /** * @brief Transform magnetic tensor data to delta_T gradient anomalies data. * * @param Mag_tensors Magnetic tensors of the magnetic data * @param deltaTs Output delta_T gradient anomalies data * @param T0_inclina Inclination degree of the earth normal magnetic field. * @param T0_declina Declination degree of the earth normal magnetic field. */ void magnetic_tensors2deltaTs(const array &Mag_tensors, array &deltaTs, double T0_inclina, double T0_declina); /** * @brief Transform magnetic components data to delta_T anomalies data. * * @note note here Mag_components[i].x is the reversed radial component. * Mag_components[i].y is the latitudinal component (south pointing). * to use it in a local cratesian coordinate, we need to reverse it. * Mag_components[i].z is the longtidinal component (east pointing) * * @param Mag_components Magnetic components of the magnetic data * @param deltaT Output delta_T anomalies data * @param T0_inclina Inclination degrees of the earth normal magnetic field. * @param T0_declina Declination degrees of the earth normal magnetic field. */ void magnetic_components2deltaT_sph(const array &Mag_components, const _1d_array &T0_inclina, const _1d_array &T0_declina, _1d_array &deltaT); /** * @brief Transform magnetic components data to delta_T anomalies data. * * @note note here Mag_components.x is the reversed radial component. * Mag_components.y is the latitudinal component (south pointing). * to use it in a local cratesian coordinate, we need to reverse it. * Mag_components.z is the longtidinal component (east pointing) * * @param Mag_components Magnetic components of the magnetic data * @param T0_inclina Inclination degrees of the earth normal magnetic field. * @param T0_declina Declination degrees of the earth normal magnetic field. */ double magnetic_components2deltaT_sph(const point3dc &Mag_components, double T0_inclina, double T0_declina); /** * @brief Transform magnetic tensor data to delta_T gradient anomalies data. * * @param Mag_tensors Magnetic tensors of the magnetic data * @param deltaTs Output delta_T gradient anomalies data * @param T0_inclina Inclination degrees of the earth normal magnetic field. * @param T0_declina Declination degrees of the earth normal magnetic field. */ void magnetic_tensors2deltaTs_sph(const array &Mag_tensors, const _1d_array &T0_inclina, const _1d_array &T0_declina, array &deltaTs); /** * @brief This function evaluates all of the Schmidt-semi normalized associated Legendre * functions up to degree nMax. The functions are initially scaled by * 10^280 sin^m in order to minimize the effects of underflow at large m * near the poles (see Holmes and Featherstone 2002, J. Geodesy, 76, 279-299). * Note that this function performs the same operation as MAG_PcupLow. * However this function also can be used for high degree (large nMax) models. * * Added into the GCTL by Yi Zhang in March 5th, 2025. * * Adopted from the FORTRAN code written by Mark Wieczorek September 25, 2005. * * Manoj Nair, Nov, 2009 Manoj.C.Nair@Noaa.Gov * * Change from the previous version * The prevous version computes the derivatives as * dP(n,m)(x)/dx, where x = sin(latitude) (or cos(colatitude) ). * However, the WMM Geomagnetic routines requires dP(n,m)(x)/dlatitude. * Hence the derivatives are multiplied by sin(latitude). * Removed the options for CS phase and normalizations. * * Note: In geomagnetism, the derivatives of ALF are usually found with * respect to the colatitudes. Here the derivatives are found with respect * to the latitude. The difference is a sign reversal for the derivative of * the Associated Legendre Functions. * * The derivatives can't be computed for latitude = |90| degrees. * * @param co_lati Co-latitude in degrees. * @param n_max Maximum spherical harmonic degree to compute. * @param P A vector of all associated Legendgre polynomials evaluated at co_lati up to nMax. * The lenght must by greater or equal to (n_max+1)*(n_max+2)/2. The coefficients are stored * in a degree major order. The index of the coefficient at degree n and order m is n*(n+1)/2 + m. * @param dP Derivative of P with respect to latitude. The coefficients are stored in the same * order as in P. */ void schmidt_legendre(double co_lati, int n_max, array &P, array &dP); /** * @brief Calculate the power spectrum of given spherical harmonic coefficients. * * @param P Returned power spectrum * @param C Cosine coefficients strats from 0 order 0 degree up to at least N order N degree * @param S Sine coefficients strats from 0 order 0 degree up to at least N order N degree * @param R Reference radius (usually is 6371.2 km) * @param r Observation radius * @param n Strat order * @param N End order */ void power_spectrum(array &P, const array &C, const array &S, double R, double r, int n, int N); }; #endif // _GCTL_GM_DATA_H