Title: | Parallel Implementation of 'ETERNA 3.40' for Prediction and Analysis of Earth Tides |
---|---|
Description: | This is a port of 'Fortran ETERNA 3.4' <http://igets.u-strasbg.fr/soft_and_tool.php> by H.G. Wenzel for calculating synthetic Earth tides using the Hartmann and Wenzel (1994) <doi:10.1029/95GL03324> or Kudryavtsev (2004) <doi:10.1007/s00190-003-0361-2> tidal catalogs. |
Authors: | Jonathan Kennel [aut, cre, trl], Beth Parker [ths], Wenzel Hans-Georg [ctb] |
Maintainer: | Jonathan Kennel <[email protected]> |
License: | GPL-3 |
Version: | 0.0.15 |
Built: | 2024-11-12 09:32:55 UTC |
Source: | https://github.com/jkennel/earthtide |
The goal of this package is to generate synthetic earth tides for use in the R programming language and in particular environmental models. Code was parallized and refactored to minimize duplication, and to allow for future improvements.
You can learn about the earthtide package in the vignettes:
browseVignettes(package = "earthtide")
Maintainer: Jonathan Kennel [email protected] [translator]
Other contributors:
Beth Parker [email protected] [thesis advisor]
Wenzel Hans-Georg [contributor]
Hartmann, T., Wenzel, H.-G., 1995. The HW95 tidal potential catalogue. Geophys. Res. Lett. 22, 3553-3556. doi:10.1029/95GL03324
Kudryavtsev, S.M., 2004. Improved harmonic development of the Earth tide-generating potential. J. Geod. 77, 829-838. doi:10.1007/s00190-003-0361-2
Wenzel, H.G., 1996. The nanogal software: Earth tide data processing package ETERNA 3.30. Bull. Inf. Marées Terrestres, 124, pp.9425-9439. https://www.eas.slu.edu/GGP/ETERNA34/MANUAL/ETERNA33.HTM
Useful links:
This is a wrapper to the Earthtide R6 class for the prediction of Earth tides. This function is provided for users who would prefer a more typical R function.
calc_earthtide( utc, do_predict = TRUE, method = "gravity", astro_update = 1, latitude = 0, longitude = 0, elevation = 0, azimuth = 0, gravity = 0, earth_radius = 6378136.3, earth_eccen = 0.0066943979514, cutoff = 1e-06, wave_groups = NULL, catalog = "ksm04", eop = NULL, return_matrix = FALSE, scale = TRUE, ... )
calc_earthtide( utc, do_predict = TRUE, method = "gravity", astro_update = 1, latitude = 0, longitude = 0, elevation = 0, azimuth = 0, gravity = 0, earth_radius = 6378136.3, earth_eccen = 0.0066943979514, cutoff = 1e-06, wave_groups = NULL, catalog = "ksm04", eop = NULL, return_matrix = FALSE, scale = TRUE, ... )
utc |
The date-time in UTC (POSIXct vector). |
do_predict |
run in predict or analyze mode |
method |
One or more of "gravity", "tidal_potential", "tidal_tilt", "vertical_displacement", "horizontal_displacement", "n_s_displacement", "e_w_displacement", "vertical_strain", "areal_strain", "volume_strain", "horizontal_strain", or "ocean_tides", "pole_tide", "lod_tide". The pole tide and lod_tide are used in predict mode even if do_predict is FALSE. More than one value can only be used if do_predict == TRUE. |
astro_update |
Integer that determines how often to phases are updated in number of samples. Defaults to 1 (every sample), but speed gains are realized with larger values. Typically updating every hour will have speed gains and keep precision (ie 3600 for one second data, 60 for minute data, 1 for hourly data). |
latitude |
The station latitude (numeric) defaults to 0. |
longitude |
The station longitude (numeric) defaults to 0. |
elevation |
The station elevation (m) (numeric) defaults to 0. |
azimuth |
Earth azimuth (numeric) defaults to 0. |
gravity |
Gravity at the station (m/s^2) (numeric) 0 to estimate gravity from elevation and latitude. |
earth_radius |
Radius of earth (m) (numeric) defaults to 6378136.3 |
earth_eccen |
Eccentricity of earth (numeric) defaults to 6.69439795140e-3 |
cutoff |
Cutoff amplitude for constituents (numeric) defaults to 1e-6. |
wave_groups |
Two column data.frame having start and end of frequency groups (data.frame). This data.frame must have two columns with the names 'start', and 'end' signifying the start and end of the wave groupings. An optional third column 'multiplier' can be provided to scale the particular wave group. If column names do no match, the inferred column positions are start, end, multiplier. |
catalog |
Use the "hw95s" catalog or "ksm04" catalog (character). |
eop |
User defined Earth Orientation Parameter (EOP) data.frame with the following columns: datetime, ddt, ut1_utc, lod, x, y, dx, dy |
return_matrix |
Return a matrix of tidal values instead of data.frame. The datetime column will not be present in this case (logical). |
scale |
Scale results when do_predict is FALSE |
... |
Currently not used. |
data.frame of tidal results
tms <- as.POSIXct('1990-01-01', tz = 'UTC') + c(0, 3600) wave_groups = data.frame(start = 0, end = 8, multiplier = 1.5) et <- calc_earthtide(utc = tms, do_predict = TRUE, method = c('tidal_potential', 'lod_tide', 'pole_tide'), astro_update = 1, latitude = 52.3868, longitude = 9.7144, elevation = 110, gravity = 9.8127, cutoff = 1.0e-5, catalog = 'ksm04', wave_groups = wave_groups)
tms <- as.POSIXct('1990-01-01', tz = 'UTC') + c(0, 3600) wave_groups = data.frame(start = 0, end = 8, multiplier = 1.5) et <- calc_earthtide(utc = tms, do_predict = TRUE, method = c('tidal_potential', 'lod_tide', 'pole_tide'), astro_update = 1, latitude = 52.3868, longitude = 9.7144, elevation = 110, gravity = 9.8127, cutoff = 1.0e-5, catalog = 'ksm04', wave_groups = wave_groups)
Class to generate synthetic earthtide signals.
An R6Class
generator object
et <- Earthtide$new( utc = as.POSIXct("2017-01-01", tz = "UTC") + 0:(24 * 7) * 3600, latitude = 52.3868, longitude = 9.7144, catalog = "ksm04", wave_groups = data.frame(start = 0.0, end = 6.0)) et$predict(method = "gravity", astro_update = 1) et$analyze(method = "gravity", astro_update = 1) et$lod_tide() et$pole_tide() et$tide() et$print()
Earthtide$new
et: An Earthtide
object.
utc: The date-time in UTC (POSIXct vector).
latitude: The station latitude (WGS84) (degree) (numeric) defaults to 0.
longitude: The station longitude (WGS84) (degree) (numeric) defaults to 0.
elevation: The station ellipsoidal height (WGS84) (m) (numeric) defaults to 0.
azimuth: Earth azimuth (numeric) defaults to 0 (degrees)
gravity: Gravity at the station (m/s^2) (numeric) 0 to estimate gravity from elevation and latitude.
earth_radius: Radius of earth (m) (numeric) defaults to 6378136.3
earth_eccen: Eccentricity of earth (numeric) defaults to 6.69439795140e-3
cutoff: Cutoff amplitude for constituents (numeric) defaults to 1e-6
wave_groups: Two column data.frame having start and end of frequency groups (data.frame). This data.frame must have two columns with the names 'start', and 'end' signifying the start and end of the wave groupings. An optional third column 'multiplier' can be provided to scale the particular wave group. If column names do no match, the inferred column positions are start, end, multiplier.
catalog: Use the "hw95s" catalog or "ksm04" catalog (character).
eop: User defined Earth Orientation Parameter (EOP) data.frame with the following columns: datetime, ddt, ut1_utc, lod, x, y, dx, dy
...: Currently not used.
Earthtide$predict, Earthtide$analyze
method: For predict
and analyze
. One of "gravity",
"tidal_potential", "tidal_tilt", "vertical_displacement",
"horizontal_displacement", "n_s_displacement", "e_w_displacement",
"vertical_strain", "areal_strain", "volume_strain", "horizontal_strain"
or "ocean_tides".
astro_update: For predict
and analyze
. Integer that
determines how often to phases are updated in number of samples. Defaults
to 1 (every sample), but speed gains are realized with larger values.
Typically updating every hour will have speed gains and keep precision
(ie 3600 for one second data, 60 for minute data, 1 for hourly data).
return_matrix: For predict
and analyze
. Return a
matrix of tidal values instead of data.frame. The datetime column will
not be present in this case (logical).
$new(utc, latitude, longitude, elevation, azimuth, gravity,
earth_radius, earth_eccen, cutoff, wave_groups, catalog, ...)
create a new Earthtide
object and initialize catalog, station and times.
$predict(method, astro_argument, return_matrix)
generate a combined
synthetic Earth tide.
$analyze(method, astro_argument, return_matrix, scale)
generate
components of the Earth tide for analysis.
$lod_tide()
generate components of the LOD (Length Of Day) tide.
$pole_tide()
generate components of the pole tide.
$tide()
get the tide data.frame
.
$print()
print the Earthtide
object.
Hartmann, T., Wenzel, H.-G., 1995. The HW95 tidal potential catalogue. Geophys. Res. Lett. 22, 3553-3556. doi:10.1029/95GL03324
Kudryavtsev, S.M., 2004. Improved harmonic development of the Earth tide-generating potential. J. Geod. 77, 829-838. doi:10.1007/s00190-003-0361-2
Wenzel, H.G., 1996. The nanogal software: Earth tide data processing package ETERNA 3.30. Bull. Inf. Marées Terrestres, 124, pp.9425-9439. https://www.eas.slu.edu/GGP/ETERNA34/MANUAL/ETERNA33.HTM
et <- Earthtide$new( utc = as.POSIXct("2017-01-01", tz = "UTC") + 0:(24 * 7) * 3600, latitude = 52.3868, longitude = 9.7144, catalog = "ksm04", wave_groups = data.frame(start = 0.0, end = 6.0)) et$predict(method = "gravity", astro_update = 1) plot(gravity~datetime, et$tide(), type='l')
et <- Earthtide$new( utc = as.POSIXct("2017-01-01", tz = "UTC") + 0:(24 * 7) * 3600, latitude = 52.3868, longitude = 9.7144, catalog = "ksm04", wave_groups = data.frame(start = 0.0, end = 6.0)) et$predict(method = "gravity", astro_update = 1) plot(gravity~datetime, et$tide(), type='l')
This data.frame contains wavegroups for different data time spans. The wavegroups should be subset prior to use and the 'time' column provides guidelines based on your input time span.
eterna_wavegroups
eterna_wavegroups
A data.frame
The columns are:
name
wave group name
start
lowest frequency of the wave group
end
highest frequency of the wave group
time
applicable to data of what length
utils::data(eterna_wavegroups)
utils::data(eterna_wavegroups)
get_iers
returns a data.frame
of earth orientation
parameters from (1962-present). This function requires an active internet connection.
Bulletins A and B are combined giving precedence to B.
Approximately (~ 7 MB) of data are downloaded. This function is brittle and
may fail when data sources change.
get_iers(a_path = NULL, b_path = NULL, daily_path = NULL, tai_utc_path = NULL)
get_iers(a_path = NULL, b_path = NULL, daily_path = NULL, tai_utc_path = NULL)
a_path |
ftp or http path to download IERS bullitin A |
b_path |
ftp or http path to download IERS bullitin B |
daily_path |
ftp or http path to download IERS daily data |
tai_utc_path |
ftp or http path to tai-utc data |
data.frame
of earth orientation parameters with the following
columns: datetime, ddt, ut1_utc, lod, x, y, dx, dy.
## Not run: eop <- get_iers() ## End(Not run)
## Not run: eop <- get_iers() ## End(Not run)
Get the frequency of the wave with the maximum amplitude in a range.
get_main_frequency(start, end)
get_main_frequency(start, end)
start |
the starting frequency in cycles per day (numeric) |
end |
the ending frequency in cycles per day (numeric) |
the main frequency between start and end