mom_harmonic_analysis module reference

Inline harmonic analysis (conventional)

More…

Data Types

ha_type

The private control structure for storing the HA info of a particular field.

ha_node

A linked list of control structures that store the HA info of different fields.

harmonic_analysis_cs

The public control structure of the MOM_harmonic_analysis module.

Functions/Subroutines

ha_init()

This subroutine sets static variables used by this module and initializes CSlist.

ha_register()

This subroutine registers each of the fields on which HA is to be performed.

ha_accum_ftf()

This subroutine accumulates the temporal basis functions in FtF.

ha_accum_ftssh()

This subroutine accumulates the temporal basis functions in FtSSH and then calls HA_write to compute harmonic constants and write results.

ha_write()

This subroutine computes the harmonic constants and write output for the current field.

ha_solver()

This subroutine computes the harmonic constants (stored in x) using the dot products of the temporal basis functions accumulated in FtF, and the dot products of the SSH (or other fields) with the temporal basis functions accumulated in FtSSH.

Detailed Description

Inline harmonic analysis (conventional)

Type Documentation

type mom_harmonic_analysis/ha_type

The private control structure for storing the HA info of a particular field.

Type fields:
  • % is [integer,private] :: Lower and upper bounds of input data.

  • % ie [integer,private] :: Lower and upper bounds of input data.

  • % js [integer,private] :: Lower and upper bounds of input data.

  • % je [integer,private] :: Lower and upper bounds of input data.

  • % key [character (len=16),private] :: Name of the field of which harmonic analysis is to be performed.

  • % grid [character (len=1),private] :: The grid on which the field is defined (‘h’, ‘q’, ‘u’, or ‘v’)

  • % old_time [real,private] :: The time of the previous accumulating step [T ~> s].

  • % ref [real(:,:),allocatable, private] :: The initial field in arbitrary units [A].

  • % ftssh [real(:,:,:),allocatable, private] :: Accumulator of (F’ * SSH_in) in arbitrary units [A].

type mom_harmonic_analysis/ha_node

A linked list of control structures that store the HA info of different fields.

Type fields:
  • % this [type( ha_type ),private] :: Control structure of the current field in the list.

  • % next [type( ha_node ),pointer, private] :: The list of other fields.

type mom_harmonic_analysis/harmonic_analysis_cs

The public control structure of the MOM_harmonic_analysis module.

Type fields:
  • % haready [logical] :: If true, perform harmonic analysis.

  • % time_start [type(time_type)] :: Start time of harmonic analysis.

  • % time_end [type(time_type)] :: End time of harmonic analysis.

  • % time_ref [type(time_type)] :: Reference time (t = 0) used to calculate tidal forcing.

  • % freq [real( max_constituents )] :: The frequency of a tidal constituent [T-1 ~> s-1].

  • % phase0 [real( max_constituents )] :: The phase of a tidal constituent at time 0 [rad].

  • % tide_fn [real( max_constituents )] :: Amplitude modulation of tides by nodal cycle [nondim].

  • % tide_un [real( max_constituents )] :: Phase modulation of tides by nodal cycle [rad].

  • % ftf [real(:,:),allocatable] :: Accumulator of (F’ * F) for all fields [nondim].

  • % nc [integer] :: The number of tidal constituents in use.

  • % length [integer] :: Number of fields of which harmonic analysis is to be performed.

  • % const_name [character (len=16)( max_constituents )] :: The name of each constituent.

  • % path [character (len=255)] :: Path to directory where output will be written.

  • % us [type( unit_scale_type )] :: A dimensional unit scaling type.

  • % list [type( ha_node ),pointer] :: A linked list for storing the HA info of different fields.

Function/Subroutine Documentation

subroutine mom_harmonic_analysis/ha_init(Time, US, param_file, time_ref, nc, freq, phase0, const_name, tide_fn, tide_un, CS)

This subroutine sets static variables used by this module and initializes CSlist. THIS MUST BE CALLED AT THE END OF tidal_forcing_init.

Parameters:
  • time :: [in] The current model time

  • time_ref :: [in] Reference time (t = 0) used to calculate tidal forcing

  • us :: [in] A dimensional unit scaling type

  • param_file :: [in] A structure to parse for run-time parameters

  • freq :: [in] The frequency of a tidal constituent [T-1 ~> s-1]

  • phase0 :: [in] The phase of a tidal constituent at time 0 [rad]

  • tide_fn :: [in] Amplitude modulation of tides by nodal cycle [nondim].

  • tide_un :: [in] Phase modulation of tides by nodal cycle [rad].

  • nc :: [in] The number of tidal constituents in use

  • const_name :: [in] The name of each constituent

  • cs :: [out] Control structure of the MOM_harmonic_analysis module

Call to:

mom_error_handler::mom_error mom_error_handler::mom_mesg

Called from:

mom_tidal_forcing::tidal_forcing_init

subroutine mom_harmonic_analysis/ha_register(key, grid, CS)

This subroutine registers each of the fields on which HA is to be performed.

Parameters:
  • key :: [in] Name of the current field

  • grid :: [in] The grid on which the key field is defined

  • cs :: [inout] Control structure of the MOM_harmonic_analysis module

Called from:

mom_tidal_forcing::tidal_forcing_init

subroutine mom_harmonic_analysis/ha_accum_ftf(Time, CS)

This subroutine accumulates the temporal basis functions in FtF. The tidal constituents are those used in MOM_tidal_forcing, plus the mean (of zero frequency). Only the main diagonal and entries below it are calculated, which are needed for Cholesky decomposition.

Parameters:
  • time :: [in] The current model time

  • cs :: [inout] Control structure of the MOM_harmonic_analysis module

Called from:

mom::step_mom

subroutine mom_harmonic_analysis/ha_accum_ftssh(key, data, Time, G, CS)

This subroutine accumulates the temporal basis functions in FtSSH and then calls HA_write to compute harmonic constants and write results. The tidal constituents are those used in MOM_tidal_forcing, plus the mean (of zero frequency).

Parameters:
  • key :: [in] Name of the current field

  • data :: [in] Input data of which harmonic analysis is to be performed [A]

  • time :: [in] The current model time

  • g :: [in] The ocean’s grid structure

  • cs :: [inout] Control structure of the MOM_harmonic_analysis module

Call to:

ha_write mom_error_handler::mom_error

Called from:

mom_barotropic::btstep mom::step_mom

subroutine mom_harmonic_analysis/ha_write(ha1, Time, G, CS)

This subroutine computes the harmonic constants and write output for the current field.

Parameters:
  • ha1 :: [in] Control structure for the current field

  • time :: [in] The current model time

  • g :: [in] The ocean’s grid structure

  • cs :: [in] Control structure of the MOM_harmonic_analysis module

Call to:

mom_io::create_mom_file ha_solver mom_io::var_desc

Called from:

ha_accum_ftssh

subroutine mom_harmonic_analysis/ha_solver(ha1, nc, FtF, x)

This subroutine computes the harmonic constants (stored in x) using the dot products of the temporal basis functions accumulated in FtF, and the dot products of the SSH (or other fields) with the temporal basis functions accumulated in FtSSH. The system is solved by Cholesky decomposition,.

FtF * x = FtSSH,    =>    L * (L' * x) = FtSSH,    =>    L * y = FtSSH,

where L is the lower triangular matrix, y = L’ * x, and x is the solution vector.

Parameters:
  • ha1 :: [in] Control structure for the current field

  • nc :: [in] Number of harmonic constituents

  • ftf :: [in] Accumulator of (F’ * F) for all fields [nondim]

  • x :: [out] Solution vector of harmonic constants [A]

Called from:

ha_write