Contact: B. Wyman Reviewers: Change history: WebCVS Log for bgrid_core_driver.f90
Provides high-level interfaces to the FMS B-grid dynamical core that allow easy initialization, integration, and diagnostics.
This module performs the following: - reads namelist input, - calls routines to set up B-grid core constants, - read and write restart files, and - output diagnostic fields.
bgrid_core_mod bgrid_horiz_mod bgrid_vert_mod bgrid_prog_var_mod bgrid_diagnostics_mod bgrid_integrals_mod bgrid_conserve_energy_mod field_manager_mod tracer_manager_mod time_manager_mod fms_mod
use bgrid_core_driver_mod [ ,only: bgrid_dynam_type, bgrid_core_driver_init, bgrid_core_driver, bgrid_core_driver_end, bgrid_core_time_diff, get_bottom_data, put_bottom_data ] bgrid_dynam_type A data structure (derived-type variable) that contains various constants used by the B-grid dynamical core. For more details see the documentation for module bgrid_core_mod. bgrid_core_driver_init Initializes the B-grid dynamical core. This interface returns values for data structures of type "bgrid_dynam_type" and "prog_var_type". Also internally initialized are other B-grid derived-type variables for the horizontal and vertical grid constants. bgrid_core_driver A wrapper for integrating the dynamical core one (atmospheric) time step. Note that only the tendencies of prognostic variables are updated. bgrid_core_time_diff Performs the time differencing of the prognostics variables and then outputs diagnostics for the B-grid dynamical core. bgrid_core_driver_end A wrapper for terminating the dynamical core. get_bottom_data, put_bottom_data Routines for "extracting" or "inserting" a 2-D field of data at the lowest model level "from" or "into" a 3-dimensional data field. NOTES * A namelist interface called bgrid_core_driver_nml is read from file input.nml.
type (bgrid_dynam_type) This data structure is inherited from and initialize in module bgrid_core_mod. See bgrid_core_mod for details.
call bgrid_core_driver_init ( Time_init, Time, Time_step, Var, Var_dt, Dynam, phys_axes ) DESCRIPTION Returns initialized/allocated values for the "bgrid_dynam_type" and "prog_var_type" data structures. Also internally initialized are other B-grid data structures for the horizontal and vertical grid constants. INPUT Time_init The initial (or base) time. [time_type] Time The current time. [time_type] Time_step The atmospheric model/physics time step. [time_type] INPUT/OUTPUT Var A derived-type variable that contains the prognostic variables for the B-grid dynamical core. The returned values will have been initialized by prog_var_mod (most likely read from a restart file). [type(prog_var_type)] Var_dt A derived-type variable that contains the prognostic variable time tendencies. The returned value is zero. [type(prog_var_type)] Dynam A derived-type variable that contains almost everything needed by the dynamical core. [type(bgrid_dynam_type)] OUTPUT phys_axes Axis identifiers as returned by the diagnostics manager and needed for subsequent calls to the diagnostics manager. [integer, dimension(4)]
call bgrid_core_driver ( Time_diag, Var, Var_dt, Dynam, omega ) DESCRIPTION Updates the prognostic variable tendencies with the dynamical core tendencies for the current atmospheric time step. Also calls diagnostics routines for outputting the dynamical core tendencies. It is assumed that the initial prognostic variable tendencies are zero. INPUT Time_diag The diagnostics time, usually the current time + time step. [type(time_type)] Var A derived-type variable that contains the B-grid's prognostic variables. [type(prog_var_type)] INPUT/OUTPUT Var_dt A derived-type variable that contains the TENDENCIES for the B-grid's prognostic variables. [type(prog_var_type)] Dynam The derived-type variable returned by a previous call to bgrid_core_driver_init (see above). [type(bgrid_dynam_type)] OUTPUT omega The omega diagnostic (from the thermodynamic equation) in pascals per second. The array should have horizontal dimensions that are consistent with the data domain of the B-grid dynamical core. [real, dimension(ilb:,jlb:,:)]
call bgrid_core_time_diff ( omega, Time_diag, Dynam, Var, Var_dt ) DESCRIPTION Performs the time differencing of the prognostics variables and outputs diagnostics for the B-grid dynamical core. INPUT omega The pressure vertical velocity in Pascals/second. This is only needed for diagnostic purposes. [real, dimension(:,:,:)] Time_diag The diagnostics time, usually the current time + time step. [type(time_type)] Dynam The derived-type variable returned by a previous call to bgrid_core_driver_init (see above). [type(bgrid_dynam_type)] INPUT/OUTPUT Var The prognostic variables. On input quantities are at the current time step and on output they are at the next time step. [type(prog_var_type)] Var_dt The time tendencies for the prognostic variables. The output tendencies will have been set to zero. [type(prog_var_type)]
call bgrid_core_driver_end (Var, Dynam) DESCRIPTION Termination routine for the B-grid dynamical core. INPUT Var The prognostic variables. [type(prog_var_type)] INPUT/OUTPUT Dynam The derived-type variable returned by a previous call to bgrid_core_driver_init (see above). [type(bgrid_dynam_type)]
call get_bottom_data ( a, b, a_bot, b_bot, [, k_bot] ) DESCRIPTION Given a pair of 3-dimensional model fields this interface returns the 2-dimensional fields at the model level closest to the ground. If optional argument "kbot" is NOT present the returned field will be the 2-d field at k = size(a,3). INPUT a, b Three-dimension fields on the model grid. The last dimension varies from the top of the atmosphere towards the surface. [real, dimension(:,:,:)] OUTPUT a_bot, b_bot Data located at the model level closest to the ground. Must have the same size as the first two dimensions of a and b. [real, dimension(:,:)] OPTIONAL INPUT k_bot The vertical index for the model level closest to the ground. Must have the same size as a_bot and b_bot. [integer, dimension(:,:)]
call put_bottom_data ( a_bot, b_bot, a, b [, k_bot] ) DESCRIPTION Inserts two-dimensional data at the lowest model level into 3-dimensional model fields. INPUT a_bot, b_bot Data located at the model level closest to the ground. This data will be inserted into arrays a and b. [real, dimension(:,:)] INPUT/OUTPUT a, b Three-dimension fields on the model grid. [real, dimension(:,:,:)] OPTIONAL INPUT k_bot The vertical index for the model level closest to the ground. Must have the same size as a_bot and b_bot. [integer, dimension(:,:)]
&bgrid_core_driver_nml num_adjust_dt The number of adjustment time steps for each advection time step, where num_adjust_dt >= 1. [integer, default: num_adjust_dt = 3] num_advec_dt The number of advection/dynamics time steps for each atmospheric/physics time step, where num_advec_dt >= 1. [integer, default: num_advec_dt = 1] decomp The domain decomposition, where decomp(1) = x-axis decomposition, decomp(2) = y-axis decomposition. If decomp(1)*decomp(2) does not equal the number of processors the model will fail. If decomp(1)=decomp(2)=0, then decomposition is determined by the MPP_DOMAINS module. [integer, dimension(2), default: decomp = 0,0] filter_option Determines how polar filtering is performed. filter_option = 0, no polar filtering (decrease time step) = 1, obsolete scheme (NO NOT USE) = 2, default scheme (refer to the technical documentation) [integer, default: filter_option = 2] filter_weight Weight applied to the polar filter that will increase (or decrease) the strength of the standard Arakawa and Lamb filter response function. SS(new) = SS(std)**filter_weight, where SS(std) is the standard response function. [integer, default: filter_weight = 1 ] ref_lat_filter The reference latitude for polar filtering. Poleward of this latitude (in both hemispheres) the polar filter is applied. Setting this argument >= 90. will turn off polar filtering. [real, default: ref_lat_filter = 60.] do_conserve_energy If TRUE the temperature tendency will be updated to guarantee that the dynamical core conserves total energy. The correction is applied as a uniform global value. [logical, default: do_conserve_energy=.false.] pgf_scheme The scheme used to compute the pressure gradient. Specify one of the following: 'DEFAULT', 'FINITE_VOLUME' (case-insensitive). The default scheme is that of Simmons and Burridge. [character, default: pgf_scheme = 'default'] restart_output_format Format used for the output restart file. only possible values are: 'native' or 'netcdf'. [character, default: restart_output_format='netcdf'] do_average_omega If TRUE the omega diagostic returned by the dynamical core is averaged over all adjustment time steps. If FALSE then omega for the last adjustment step is returned. [logical, default: do_average_omega=.false.] ddamp_coeff Damping coefficient for divergence damping. The coefficent has been normalized by the maximum value so that ddamp_coeff should be in the range [0,1]. Use ddamp_coeff > 0 for second-order, and ddamp_coeff < 0 for fourth-order. If ddamp_coeff = 0 then no divergence damping is done. [real, default: ddamp_coeff = 0.] verbose Flag that control additional printed output. Currently, this option is not being used. [integer, default: verbose = 0] NOTES
None.
None.
None.