mom_harmonic_analysis module reference
Inline harmonic analysis (conventional)
Data Types
The private control structure for storing the HA info of a particular field. |
|
A linked list of control structures that store the HA info of different fields. |
|
The public control structure of the MOM_harmonic_analysis module. |
Functions/Subroutines
This subroutine sets static variables used by this module and initializes CSlist. |
|
This subroutine registers each of the fields on which HA is to be performed. |
|
This subroutine accumulates the temporal basis functions in FtF and FtSSH and then calls HA_write to compute harmonic constants and write results. |
|
This subroutine computes the harmonic constants and write output for the current field. |
|
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, privateLower and upper bounds of input data.% ie ::
integer, privateLower and upper bounds of input data.% js ::
integer, privateLower and upper bounds of input data.% je ::
integer, privateLower and upper bounds of input data.% key ::
character(len=16), privateName of the field of which harmonic analysis is to be performed.% grid ::
character(len=1), privateThe grid on which the field is defined (‘h’, ‘q’, ‘u’, or ‘v’)% old_time ::
real, privateThe time of the previous accumulating step [T ~> s].% ref ::
real, dimension(:,:), allocatable, privateThe initial field in arbitrary units [A].% ftf ::
real, dimension(:,:), allocatable, privateAccumulator of (F’ * F) [nondim].% ftssh ::
real, dimension(:,:,:), allocatable, privateAccumulator 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), privateControl structure of the current field in the list.% next ::
type(ha_node), pointer, privateThe list of other fields.
- type mom_harmonic_analysis/harmonic_analysis_cs
The public control structure of the MOM_harmonic_analysis module.
- Type fields:
% haready ::
logicalIf 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, dimension(:), allocatableThe frequency of a tidal constituent [T-1 ~> s-1].% phase0 ::
real, dimension(:), allocatableThe phase of a tidal constituent at time 0 [rad].% tide_fn ::
real, dimension(:), allocatableAmplitude modulation of tides by nodal cycle [nondim].% tide_un ::
real, dimension(:), allocatablePhase modulation of tides by nodal cycle [rad].% nc ::
integerThe number of tidal constituents in use.% length ::
integerNumber of fields of which harmonic analysis is to be performed.% const_name ::
character(len=4), dimension(:), allocatableThe 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), pointerA linked list for storing the HA info of different fields.
Function/Subroutine Documentation
- subroutine mom_harmonic_analysis/ha_init(Time, US, param_file, nc, 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
us :: [in] A dimensional unit scaling type
param_file ::
param_file[in] A structure to parse for run-time parametersnc ::
nc[in] The number of tidal constituents in usecs :: [out] Control structure of the MOM_harmonic_analysis module
- Call to:
mom_tidal_forcing::astro_longitudes_initmom_tidal_forcing::eq_phaseha_registermom_error_handler::mom_errormom_error_handler::mom_mesgmom_tidal_forcing::nodal_fumom_tidal_forcing::tidal_frequency- Called from:
mom_dynamics_split_rk2::initialize_dyn_split_rk2mom_dynamics_split_rk2b::initialize_dyn_split_rk2b
- 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 ::
key[in] Name of the current fieldgrid ::
grid[in] The grid on which the key field is definedcs :: [inout] Control structure of the MOM_harmonic_analysis module
- Called from:
- subroutine mom_harmonic_analysis/ha_accum(key, data, Time, G, CS)
This subroutine accumulates the temporal basis functions in FtF and 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). For FtF, only the main diagonal and entries below it are calculated, which are needed for Cholesky decomposition.
- Parameters:
key ::
key[in] Name of the current fielddata ::
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:
- Called from:
- 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 ::
ha1[in] Control structure for the current fieldtime :: [in] The current model time
g :: [in] The ocean’s grid structure
cs :: [in] Control structure of the MOM_harmonic_analysis module
- Call to:
- Called from:
- 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 ::
ha1[in] Control structure for the current fieldnc ::
nc[in] Number of harmonic constituentsftf :: [in] Accumulator of (F’ * F) for all fields [nondim]
x ::
x[out] Solution vector of harmonic constants [A]
- Called from: