ex3_1_SubseasonalECVHindcast.md 8.56 KB
Newer Older
nperez's avatar
nperez committed
# Weekly ECV Subseasonal Hindcast Verification
nperez's avatar
nperez committed
#-------

nperez's avatar
nperez committed
This is a practical case to compute skill scores in subseasonal forecast.

  - We will use  ECMWF/S2S-ENSForhc and ERA5.
  - The ECV is air temperature at surface level (tas).


Here, we aim to calculate the skill scores of the corresponding hindcast of 2016.

Note that to do this, we will detect the Mondays and Thursdays, which are the days of the week in which this model is initialized) during the year 2016. However, we will work with the previous year of initialization.

Follow the figure to see details on the subeasonal forecasts storage:

nperez's avatar
nperez committed
<img src="inst/doc/figures/subseasonal_1.png" width="600" />
nperez's avatar
nperez committed

First, load startR package and define the paths to your forecast and reference including labels for $var$ and others:

```{r}
library(startR)

ecmwf_path <- paste0('/esarchive/exp/ecmwf/s2s-monthly_ensforhc/',
                     'weekly_mean/$var$_f24h/$sdate$/$var$_$syear$.nc')
obs_path <- paste0("/esarchive/recon/ecmwf/era5/weekly_mean/",
                   "$var$_f1h-240x121/$var$_$file_date$.nc")
```

Now, create the sequence of start dates in 2016 that is explaind in the figure too:

nperez's avatar
nperez committed
<img src="inst/doc/figures/subseasonal_2.png" width="600" />
nperez's avatar
nperez committed

```{r}
forecast.year <- 2016
# Mondays
sdates.seq.mon <- format(seq(as.Date(paste(forecast.year, 01, 04, sep = '-')),
                             as.Date(paste(forecast.year, 12, 31, sep='-')), 
                             by = 'weeks'), format = '%Y%m%d')
# Thursdays (Monday+3days)
sdates.seq.thu <- format(seq(as.Date(paste(forecast.year, 01, 04, sep = '-')) + 3,
                         as.Date(paste(forecast.year, 12, 31, sep = '-')),
                         by = 'weeks'), format = '%Y%m%d')
# Joint dates in order
sdates.seq <- c(sdates.seq.mon, sdates.seq.thu)
ind <- order(as.Date(sdates.seq, format = '%Y%m%d'))  # dates in order
sdates.seq <- sdates.seq[ind]

# Leap years, remove 29th of February:
pos.bis <- which(sdates.seq == paste0(forecast.year,"0229")) 
if(length(pos.bis) != 0) sdates.seq <- sdates.seq[-pos.bis] 

exp <- Start(dat = ecmwf_path, 
              var = 'tas',
              sdate = sdates.seq,
              syear = 'all',     # before hdate
              time = 'all',
              ensemble = "all",
              latitude = indices(1:121),
              longitude = indices(1:240),
              syear_depends = 'sdate',
              return_vars = list(latitude = NULL,
                                 longitude = NULL,
                                 time = c('sdate', 'syear')), 
              retrieve = FALSE)

# Corresponding ERA5 
dates <- attr(exp, 'Variables')$common$time   #20*4*103=8240  middle day of weekly averages
file_date <- sapply(dates, format, '%Y%m%d')
dim(file_date) <- c(length(sdates.seq), 20, 4)
names(dim(file_date)) <- c('sdate', 'syear', 'time')

obs <- Start(dat = obs_path,
             var = 'tas',
             file_date = file_date,  #  sdate syear  time 
                                     #    103    20     4 
             latitude = indices(1:121),
             longitude = indices(1:240),
             split_multiselected_dims = TRUE,
             retrieve = FALSE)
```

The next step is to define our function to calculate scores:

nperez's avatar
nperez committed
<img src="inst/doc/figures/subseasonal_3.png" width="600" />
nperez's avatar
nperez committed

```{r} 
score_calc <- function(forecast, reference, sdates.seq,
                       scores = 'all') {
  library(easyVerification)
  anomaly_crossval <- function(data) {
    avg <- NA * dim(data)['syear'] # one avg per year
      for (t in 1:dim(data)['syear']) {
        avg[t]<- mean(data[,-t,])   
      }                        
    data <- Apply(data, c('syear'), function(x) x - avg)[[1]]
    return(data)
  }
  reference <- s2dv::InsertDim(reference, pos = 3, len = 1, name = 'ensemble')
  forecast <- anomaly_crossval(forecast)
  reference <- anomaly_crossval(reference)
      
  # create objects to store the outputs
  Scores <- NULL
  Skill_Scores <- NULL

  for(month in 1:12) {
    # take only data of 1 month (indices of each month in start_dates)
    startdates <- which(as.integer(substr(sdates.seq, 5,6)) == month)
    forecast_month <- s2dverification::Subset(forecast, 'sdate', list(startdates))
    reference_month <- s2dverification::Subset(reference, 'sdate', list(startdates))
    forecast_month_reshaped <- array(forecast_month, 
	          c(dim(forecast_month)['sdate'] * dim(forecast_month)['syear'],
                    dim(forecast_month)['ensemble']))
    reference_month_reshaped <- array(reference_month, 
	          c(dim(reference_month)['sdate'] * dim(reference_month)['syear'],
	            dim(reference_month)['ensemble']))
    if (any(c('fairrpss', 'frpss', 'all') %in% tolower(scores))) {
      Scores <- c(Scores, unlist(easyVerification::veriApply("FairRpss",
                  fcst = forecast_month_reshaped, obs = reference_month_reshaped,
                  prob = c(1/3, 2/3), tdim = 1, ensdim = 2)))
    } 
    if (any(c('faircrpss', 'fcrpss', 'all') %in% tolower(scores))) {
      Scores <- c(Scores, unlist(easyVerification::veriApply("FairCrpss",
                  fcst = forecast_month_reshaped, obs = reference_month_reshaped,
                  tdim = 1, ensdim = 2)))
    } 
    if (any(c('enscorr', 'corr', 'all') %in% tolower(scores))) {
      fcst_mean <- Apply(forecast_month_reshaped, 'ensemble', mean)[[1]]
      Scores <- c(Scores,
                  cor = cor(fcst_mean, reference_month_reshaped,
                            use = "complete.obs"),
                  cor.pv = cor.test(fcst_mean, reference_month_reshaped,
                           use = "complete.obs")$p.value)
    }
    if (any(c('bss10', 'fbss10', 'all') %in% scores)) {
      Scores <- c(Scores, unlist(easyVerification::veriApply("FairRpss",
                  fcst = forecast_month_reshaped, obs = reference_month_reshaped, 
                  prob = (1/10), tdim = 1, ensdim = 2)))
    }
    if (any(c('bss90', 'fbss90', 'all') %in% scores)) {
      Scores <- c(Scores, unlist(easyVerification::veriApply("FairRpss",
                  fcst = forecast_month_reshaped, obs = reference_month_reshaped,
                  prob = (9/10), tdim = 1, ensdim = 2)))
    }
    Skill_Scores <- cbind(Skill_Scores, Scores)
  }
  return(Skill_Scores)
}
```

The workflow is created with Step() and AddStep()

nperez's avatar
nperez committed
<img src="inst/doc/figures/subseasonal_4.png" width="600" />
nperez's avatar
nperez committed


```{r}
step <- Step(fun = score_calc, 
             target_dims = list(forecast = c('sdate','syear','ensemble'),
                                reference = c('sdate', 'syear')),
             output_dims = list(c('scores', 'month')),
             use_libraries = c('easyVerification',
                               'SpecsVerification',
                               's2dverification'))
# Workflow:
wf <- AddStep(list(exp,obs), step, sdates.seq = sdates.seq, scores = 'all')
```

Finally, execute the analysis, for instance, in Nord3:


```{r}
#-----------modify according to your personal info---------
  queue_host = 'nord3'   #your own host name for power9
  temp_dir =  '/gpfs/scratch/bsc32/bsc32339/startR_hpc/'
  ecflow_suite_dir = '/home/Earth/nperez/startR_local/'  #your own local directory
#------------------------------------------------------------
res <- Compute(wf,
               chunks = list(time = 4, longitude = 2, latitude = 2),
               threads_load = 1,
               threads_compute = 12,
               cluster = list(queue_host = queue_host,
                              queue_type = 'lsf',
                              extra_queue_params = list('#BSUB -q bsc_es'),
                              cores_per_job = 12,
                              temp_dir = temp_dir,
                              polling_period = 10,
                              job_wallclock = '03:00',
                              max_jobs = 16,
                              bidirectional = FALSE),
               ecflow_suite_dir = ecflow_suite_dir,
               wait = TRUE)
```

If 'all' scores are requested the order in the 'scores' dimension of the output will be:
  - FairRpss_skillscore, FairRpss_sd, 
  - FairCrpss_skillscore, FairCrpss_sd, 
  - EnsCorr, EnsCorr_pv
  - FairBSS10, FairBSS10_pv
  - FairBSS90, FairBSS90_pv

If you choose a subset of 'scores', it will follow the same order omitting the non-declared scores.

```
library(s2dv)
PlotEquiMap(res$output1[1,1,1,1,1,,], 
            lon = attributes(exp)$Variables$common$longitude,
            lat = attributes(exp)$Variables$common$latitude,
            filled.con = FALSE,
            toptitle = "FairRPSS January 2016 Hindcast leadtime 1")
```

nperez's avatar
nperez committed
<img src="inst/doc/figures/subseasonal_5.png" width="600" />
nperez's avatar
nperez committed


Credits to
Original code: Andrea Manrique
Adaptation: Núria Pérez-Zanón