mom_interpolate module reference

This module provides added functionality to the FMS temporal and spatial interpolation routines.

More…

Data Types

forcing_timeseries_dataset

Data type used to store information about forcing datasets that are time series E.g.

Functions/Subroutines

time_interp_external_0d()

Read a scalar field based on model time.

time_interp_external_2d()

Read a 2d field from an external based on model time, potentially including horizontal interpolation and rotation of the data.

time_interp_external_3d()

Read a 3d field based on model time, and rotate to the model grid.

forcing_timeseries_set_time_type_vars()

Set time_type variables in forcing_timeseries_dataset() type based on integer input TODO: make this part of type based on integer input TODO: make this part of forcing_timeseries_dataset() class if OO is okay in MOM6? class if OO is okay in MOM6?

map_model_time_to_forcing_time()

If necessary, apply an offset to convert from model time to forcing time and then ensure result is within acceptable bounds.

Detailed Description

This module provides added functionality to the FMS temporal and spatial interpolation routines.

Type Documentation

type mom_interpolate/forcing_timeseries_dataset

Data type used to store information about forcing datasets that are time series E.g. how do we align the data in the model with the time axis in the file?

Type fields:
  • % file_name [character (len=200)] :: name of file containing river flux forcing

  • % l_time_varying [logical] :: .true. => forcing is dependent on model time, .false. => static forcing

  • % data_forcing [type(time_type)] :: convert data_forcing_year to time type

  • % data_start [type(time_type)] :: convert data_start_year to time type

  • % data_end [type(time_type)] :: convert data_end_year to time type

  • % m2d_offset [type(time_type)] :: add to model time to get data time

Function/Subroutine Documentation

subroutine mom_interpolate/time_interp_external_0d(field, time, data_in, verbose, scale)

Read a scalar field based on model time.

Parameters:
  • field :: [in] Handle for time interpolated field

  • time :: [in] The target time for the data

  • data_in :: [inout] The interpolated value in arbitrary units [A ~> a]

  • verbose :: [in] If true, write verbose output for debugging

  • scale :: [in] A scaling factor that new values of data_in are multiplied by before it is returned [A a-1 ~> 1]

subroutine mom_interpolate/time_interp_external_2d(field, time, data_in, interp, verbose, horz_interp, mask_out, turns, scale)

Read a 2d field from an external based on model time, potentially including horizontal interpolation and rotation of the data.

Parameters:
  • field :: [in] Handle for time interpolated field

  • time :: [in] The target time for the data

  • data_in :: [inout] The array in which to store the interpolated values in arbitrary units [A ~> a]

  • interp :: [in] A flag indicating the temporal interpolation method

  • verbose :: [in] If true, write verbose output for debugging

  • horz_interp :: [in] A structure to control horizontal interpolation

  • mask_out :: [out] An array that is true where there is valid data

  • turns :: [in] Number of quarter turns to rotate the data

  • scale :: [in] A scaling factor that new values of data_in are multiplied by before it is returned [A a-1 ~> 1]

subroutine mom_interpolate/time_interp_external_3d(field, time, data_in, interp, verbose, horz_interp, mask_out, turns, scale)

Read a 3d field based on model time, and rotate to the model grid.

Parameters:
  • field :: [in] Handle for time interpolated field

  • time :: [in] The target time for the data

  • data_in :: [inout] The array in which to store the interpolated values in arbitrary units [A ~> a]

  • interp :: [in] A flag indicating the temporal interpolation method

  • verbose :: [in] If true, write verbose output for debugging

  • horz_interp :: [in] A structure to control horizontal interpolation

  • mask_out :: [out] An array that is true where there is valid data

  • turns :: [in] Number of quarter turns to rotate the data

  • scale :: [in] A scaling factor that new values of data_in are multiplied by before it is returned [A a-1 ~> 1]

subroutine mom_interpolate/forcing_timeseries_set_time_type_vars(data_start_year, data_end_year, data_ref_year, model_ref_year, data_forcing_year, forcing_dataset)

Set time_type variables in forcing_timeseries_dataset() type based on integer input TODO: make this part of type based on integer input TODO: make this part of forcing_timeseries_dataset() class if OO is okay in MOM6? class if OO is okay in MOM6?

Parameters:
  • data_start_year :: [in] first year of data to read (this is ignored for static forcing)

  • data_end_year :: [in] last year of data to read (this is ignored for static forcing)

  • data_ref_year :: [in] for time-varying forcing, align data_ref_year in file with model_ref_year in model

  • model_ref_year :: [in] for time-varying forcing, align data_ref_year in file with model_ref_year in model

  • data_forcing_year :: [in] for static forcing, read file at this date (this is ignored for time-varying forcing)

  • forcing_dataset :: [inout] information about forcing file

function mom_interpolate/map_model_time_to_forcing_time(Time, forcing_dataset) [type(time_type)]

If necessary, apply an offset to convert from model time to forcing time and then ensure result is within acceptable bounds.

Parameters:
  • time :: [in] Model time

  • forcing_dataset :: [in] information about forcing file

Return:

undefined :: time to read forcing file

Called from:

marbl_tracers::marbl_tracers_set_forcing