knitr::opts_chunk$set(echo = TRUE)

Codecov test coverage CRAN Downloads DOI R build status

Raytracing Documentation

The identification of the atmospheric Rossby wave ray paths is of paramount importance for atmospheric scientists, climatologists, meteorologists, and students in order to research, assign, monitor and forecast weather and climate patterns. raytracing is designed to detect Rossby waves trigger regions, characteristics, and paths, using a zonally symmetric basic state.

It relies on the theory described in:

For a brief review, see Rehbein et al. (2020).

Including Code

The raytracing include the below functions:

library(knitr)
df <- data.frame(
  Functions = c("betaks", "ray", "ray_path", "ray_source", "trin", "ypos"),
  Arguments = c("u, lat, lon, uname, ofile, a, plots, show.warnings", 
                "betam, u, lat, x0, y0, K, dt, itime, direction, interpolation, tl, a, verbose, ofile",
                "x, y",
                "betam, u, lat, x0, y0, K, dt, itime, direction, interpolation, tl, a, verbose, ofile",
                "y, yk, mercator",
                "y, lat, yk, mercator"), 
  Description = c("Calculates Beta and Ks",
                  "Calculates the Rossby waves ray paths",
                  "Calculate the ray paths / great circles",
                  "Calculate the Rossby waves ray paths over a source region",
                  "Performs trigonometric interpolation", 
                  "Interpolation selecting the nearest neighbor")
)
knitr::kable(x = df, caption = "**Table** **1** **-** ```raytracing``` functions, arguments, and its description.")

Instalation

To install the CRAN version:

install.packages("raytracing")

To install the developing version:

remotes::install_github("salvatirehbein/raytracing")

Example

Simple usage example, reproduced from Coelho et al. (2015).

library(raytracing)
input <- system.file("extdata",
                     "uwnd.mon.mean_200hPa_2014JFM.nc",
                     package = "raytracing")
b <- betaks(u = input)
rt <- ray_source(betam = b$betam,
                 u = b$u,
                 lat = b$lat,
                 K = 3,
                 itime = 10*4,
                 x0 = -c(130, 135),
                 y0 = -30,
                 dt = 6,
                 direction = -1,
                 interpolation = "trin")

The ray or the ray_source functions return a sf data.table as below:

head(rt)
tail(rt)
# Large stationary wave numbers (Ks) or short wavelength are set to be Ks = 10
b$sfpoly$ksm <- ifelse(b$sfpoly$ksm > 10 & b$sfpoly$ksm <= 20, 10, 
                ifelse(b$sfpoly$ksm >=30 | b$sfpoly$ksm < 0, NA, 
                b$sfpoly$ksm))
# Plot Ks
plot(b$sfpoly["ksm"],
     lty = 0,
     axes = TRUE, 
     reset = FALSE,
     breaks = seq(0, 10, 2),
     pal= c("#f0ff00", "#ffce00", "#ff9a00", "#ff5a00", "#ff0000"),
     main = "Stationary Wave Number and Ray Tracings: JFM/2014")

# Plot contour maps
data(coastlines)
plot(coastlines, add = TRUE)
# Select linestrings for plotting only the lines
li <- rt[st_is(rt$geometry, "LINESTRING"), ]

# Plot ray traces
plot(li["lon_ini"],
    add = TRUE,
    lwd = 2,
    pal = colorRampPalette(c("black", "blue")))

This figure is a reproduction of the Figure 9 from Coelho et al. (2015).

citation

To cite raytracing in publications use this:

Rehbein, A., Ambrizzi, T., Ibarra-Espinosa, S., Dutra, L. M. M.: Rossby Wave Ray Tracing v0.1.0. https://github.com/salvatirehbein/raytracing, 2020.

A BibTeX entry for LaTeX users is

@Manual{ray,
    title = {raytracing: An R package for identification and tracking the atmospheric Rossby waves},
    author = {Amanda Rehbein and Tercio Ambrizzi and Sergio Ibarra-Espinosa and Livia M. M. Dutra},
    year = {2020},
    url = {https://github.com/salvatirehbein/raytracing},
  }

Thanks



salvatirehbein/raytracing documentation built on June 11, 2022, 11:36 a.m.