grace_input_months
Reads GRACE/GRACE-FO/Swarm files for a specified spherical harmonic degree and order and for a specified date range
Calling Sequence
from gravity_toolkit.grace_input_months import grace_input_months
GRACE_Ylms = grace_input_months(base_dir, PROC, DREL, DSET, LMAX,
start_mon, end_mon, missing, SLR_C20, DEG1, SLR_C30=SLR_C30)
- gravity_toolkit.grace_input_months(base_dir, PROC, DREL, DSET, LMAX, start_mon, end_mon, missing, SLR_C20, DEG1, **kwargs)[source]
Reads GRACE/GRACE-FO files for a spherical harmonic degree and order and a date range
Can include geocenter values for degree 1 coefficients
Can replace C20 with SLR values for all months
Can replace low-degree harmonics with SLR values for months 179+
Can correct for ECMWF atmospheric “jumps” using GAE/GAF/GAG files [18]
Can correct for Pole Tide drift following Wahr et al. [74]
- Parameters:
- base_dir: str
Working data directory for GRACE/GRACE-FO data
- PROC: str
GRACE/GRACE-FO/Swarm data processing center
'CSR': University of Texas Center for Space Research'GFZ': German Research Centre for Geosciences (GeoForschungsZentrum)'JPL': Jet Propulsion Laboratory'CNES': French Centre National D’Etudes Spatiales'GRAZ': Institute of Geodesy from GRAZ University of Technology'COSTG': Combination Service for Time-variable Gravity Fields'Swarm': Time-variable gravity data from Swarm satellites
- DREL: str
GRACE/GRACE-FO/Swarm data release
- DSET: str
GRACE/GRACE-FO/Swarm data product
'GAA': non-tidal atmospheric correction'GAB': non-tidal oceanic correction'GAC': combined non-tidal atmospheric and oceanic correction'GAD': ocean bottom pressure product'GSM': corrected monthly static gravity field product
- LMAX: int
Upper bound of Spherical Harmonic Degrees
- start_mon: int
starting month to consider in analysis
- end_mon: int
ending month to consider in analysis
- missing: list
missing months to not consider in analysis
- SLR_C20: str
Replaces C20 with SLR values
None: use original values'CSR': use values from CSR (TN-07, TN-09, TN-11)'GFZ': use values from GFZ'GSFC': use values from GSFC (TN-14)
- DEG1: str
Use Degree 1 coefficients
- MMAX: int or NoneType, default None
Upper bound of Spherical Harmonic Orders
- SLR_21: str or NoneType, default ‘’
Replace C21 and S21 with SLR values
None: use original values'CSR': use values from CSR'GFZ': use values from GFZ GravIS'GSFC': use values from GSFC
- SLR_22: str or NoneType, default ‘’
Replace C22 and S22 with SLR values
None: use original values'CSR': use values from CSR'GSFC': use values from GSFC
- SLR_C30: str or NoneType, default ‘’
Replace C30 with SLR values
None: use original values'CSR': use values from CSR (5x5 with 6,1)'GFZ': use values from GFZ GravIS'GSFC': use values from GSFC (TN-14)
- SLR_C40: str or NoneType, default ‘’
Replace C40 with SLR values
None: use original values'CSR': use values from CSR (5x5 with 6,1)'GSFC': use values from GSFC
- SLR_C50: str or NoneType, default ‘’
Replace C50 with SLR values
None: use original values'CSR': use values from CSR (5x5 with 6,1)'GSFC': use values from GSFC
- POLE_TIDE: bool, default False
Correct GSM data with pole tides following Wahr et al. [74]
- ATM: bool, default False
Correct data with ECMWF “jump” corrections following Fagiolini et al. [18]
- DEG1_FILE: str or NoneType, default None
full path to degree 1 coefficients file
- MODEL_DEG1: bool, default False
least-squares model missing degree 1 coefficients
- Returns:
- clm: np.ndarray
cosine spherical harmonics to degree/order
LMAXandMMAX- slm: np.ndarray
sine spherical harmonics to degree/order
LMAXandMMAX- eclm: np.ndarray
uncalibrated cosine spherical harmonic errors
- eslm: np.ndarray
uncalibrated sine spherical harmonic errors
- time: np.ndarray
time of each measurement (mid-month)
- month: np.ndarray
GRACE/GRACE-FO months of input datasets
- l: np.ndarray
spherical harmonic degree to
LMAX- m: np.ndarray
spherical harmonic order to
MMAX- title: str
Processing string denoting low degree zonals replacement, geocenter usage and corrections
- directory: str
Directory of exact GRACE/GRACE-FO/Swarm product
- attributes: dict
Attributes of input files and corrections
- gravity_toolkit.grace_input_months.read_ecmwf_corrections(base_dir, LMAX, months, MMAX=None)[source]
Read atmospheric jump corrections from [18]
- Parameters:
- base_dir: str
Working data directory for GRACE/GRACE-FO data
- LMAX: int
Upper bound of Spherical Harmonic Degrees
- months: list
list of GRACE/GRACE-FO months
- MMAX: int or NoneType, default None
Upper bound of Spherical Harmonic orders
- Returns:
- clm: float
atmospheric correction cosine spherical harmonics
- slm: float
atmospheric correction sine spherical harmonics
- files: list
atmospheric correction files
- gravity_toolkit.grace_input_months.regress_model(t_in, d_in, t_out, ORDER=2, CYCLES=None, RELATIVE=0.0)[source]
Calculates a regression model for extrapolating values
- Parameters:
- t_in: float
input time array
- d_in: float
input data array
- t_out: float
time array for output regressed values
- ORDER: int, default 2
maximum polynomial order for regression model
- CYCLES: list or NoneType, default None
list of cyclical terms to include in regression model
- RELATIVE: float, default 0.0
relative time for polynomial coefficients in fit
- Returns:
- d_out: float
output regressed value data array