Vgrid/lowlevel
|
The lower-level interface of the Vertical Grid Descriptors module is provided for users who wish to have a greater level of control over the actions of the Vertical Grid Descriptor module than is available using the (very thin) high-level interface. Using the methods described here, the user can query and set all components of the vertical descriptors, and employ user-defined reference fields for calculations.
All methods (functions) in the Vertical Grid Descriptor module are prefixed with vgd_, as are all public class parameters (VGD_). All functions return a status code as VGD_OK for successful completion or VGD_ERROR if an error has occurred.
Contents[hide] |
1 Vertical Descriptor Type
This derived type contains information about the vertical grid descriptors of a grid. It can be instatiated (constructed, or created) using the vgd_new() function described below. The declaration of a set of vertical grid descriptors (an instance of the Vertical Grid Descriptor class) looks like,type(vgrid_descriptor) :: grid_desc
2 Type Constructor
In order to begin using a Vertical Grid Descriptor object, it must be instatiated (constructed). The constructor in the Vertical Grid Descriptor Module is called vgd_new since it creates a new instance of class. The standard invocation of the constructor consists of providing the ip1 and ip2 values of the desired vertical descriptors, along with the file unit number.status = vgd_new(grid_desc,unit=10,format='fst',ip1=100,ip2=200)
An pair of alternative constructors are provided for users that wish to create a vertical grid descriptor without RPN Standard file input.
3 Available Methods
The following methods are available to all procedures that use-associate the Vertical Grid Descriptor module. Note that the constructor should be called in order to fill the vgrid_descriptor derived type (instatiate the Vertical Grid Descriptor object) before any other methods are called. All methods return the a status integer variable that contains either VGD_OK or VGD_ERROR (see variable descriptions) depending on the completion status of the method. The name self is used for the Vertical Grid Descriptor object in all of the prototypes listed below, but the variable can be of any name in the calling sequence.
3.1 Constructors
The primary interface used to initialize a Vertical Grid Descriptor object involves using existing vertical descriptors in a file. This is the standard constructor that will be used by most applications.
status = vgd_new(self,unit,format,ip1,ip2)
type(vgrid_descriptor), intent(out) :: self !Vertical descriptor instance
integer, intent(in) :: unit !File unit to read descriptor information from
character(len=*), intent(in) :: format !File format ('fst' or 'bin')
integer, intent(in) :: ip1,ip2 !ip1,2 values of the desired descriptors
A secondary interface is provided to allow for the generation of a Vertical Grid Descriptor object given the vertical coordinate information contained in the model's settings file (i.e. gem_settings.nml). This constructor will be used by applications that wish to create a specific type of vertical coordinate.
status = vgd_new(self,kind,version,hyb,rcoef1,rcoef2,ptop_8,pref_8,ip1,ip2,ip3,stdout_unit) result(status)
type(vgrid_descriptor), intent(out) :: self !Vertical descriptor instance
integer, intent(in) :: kind,version !Kind,version to create
real, dimension(:), intent(in) :: hyb !List of hybrid levels
real, intent(in) :: rcoef1,rcoef2 !R-coefficient values for rectification
real*8, intent(in) :: ptop_8 !Top-level pressure (Pa)
real*8, intent(in) :: pref_8 !Reference-level pressure (Pa)
integer, optional, intent(in) :: ip1,ip2,ip3 !IP1,2,3 values for FST file record [0,0,0]
integer, intent(in), optional :: stdout_unit !Unit number for verbose output [STDERR]
The following table gives which parameter must be set for a given coordinate. An other convenient way to get this list is to call the function with only its mandatory parameters (self,kind,version and hyb), the function will return the list of parameter needed for the given kind and version.
|
library | hyb | rcoef1 | rcoef2 | ptop_8 | pref_8 | ip1 | ip2 | ip3 |
---|---|---|---|---|---|---|---|---|---|
kind=1, Version=1 sigma | 2.0.4 and above | yes | no | no | no | no | optional | optional | optional |
kind=1, Version=2 eta | 2.0.3 and above | yes | no | no | yes | no | optional | optional | optional |
kind=1, Version=3 hybrid normalized | 2.0.2 and above | yes | yes | no | yes | yes | optional | optional | optional |
kind=2, Version=1 pressure | 1.0.0 and above | yes | no | no | no | no | optional | optional | optional |
kind=5, Version=1 hybrid un-staggered | 2.0.2 and above | yes | yes | no | yes | yes | optional | optional | optional |
kind=5, Version=2 hybrid staggered, first level is a thermo level | 2.0.2 and avobe | yes | yes | yes | yes | yes | optional | optional | optional |
kind=5, Version=3 hybrid staggered, first level is a thermo level, unstaggered last Thermo level | 2.0.5 and avobe | yes | yes | yes | yes | yes | optional | optional | optional |
kind=5, Version=4 hybrid staggered, first level is a momentum level, same number of thermo and momentum levels | 5.1.0 and avobe | yes | yes | yes | yes, possible values -1, -2 or any pressure value. Recommended value is -2 which gives a flat first (top) momentum level [B(1)=0], this is what GEM 4.7 uses. | yes | optional | optional | optional |
Alternatively, a complete constructor is provided to initialize a pressure-type vertical coordinate component of a Vertical Grid Descriptor directly. The use of this constructor requires an intimate knowledge of the coordinate structure. Specifically, different sets of arguments are required for different kind and version combinations. In general, only the model will use this form of the Vertical Grid Descriptor constructor. The package will alert users if required arguments are missing, so a call to vgd_new with a minimal (non-optional) argument list will usually be the basis for development.
status = vgd_new(self,kind,version,nk,ip1,ip2,ptop_8,pref_8,rcoef1,rcoef2, &
a_m_8,b_m_8,a_t_8,b_t_8,ip1_m,ip1_t) result(status)
type(vgrid_descriptor), intent(out) :: self !Vertical descriptor instance
integer, intent(in) :: kind,version !Kind,version to create
integer, intent(in) :: nk !Number of levels
integer, optional, intent(in) :: ip1,ip2 !IP1,2 values for FST file record [0,0]
real, optional, intent(in) :: rcoef1,rcoef2 !R-coefficient values for rectification
real*8, optional, intent(in) :: ptop_8 !Top-level pressure (Pa)
real*8, optional, intent(in) :: pref_8 !Reference-level pressure (Pa)
real*8, optional, dimension(:) :: a_m_8,a_t_8 !A-coefficients for momentum(m),thermo(t) levels
real*8, optional, dimension(:) :: b_m_8,b_t_8 !B-coefficients for momentum(m),thermo(t) levels
integer, optional, dimension(:) :: ip1_m,ip1_t !Level ID (IP1) for momentum(m),thermo(t) levels
In some context, like in a MPI environment, the real*8 table containing the vgrid_descriptor information may become available a may need to be reintroduce in the vgrid_descriptor structure. For such a context the following may be used :
status = vgd_new(self,table)
type(vgrid_descriptor), intent(inout) :: self !Vertical descriptor instance
real(kind=8), dimension(:,:,:), pointer :: table !Raw table of vgrid records
Since no FST-specific information is contained in the table (etiket record stamp string, ig1-4 values, et cetera), all such record-specific information is set to default values when this constructor is used. If this information is known, then the definition methods can be used to set the required values correctly after construction. For this reason, the use of this constructor should be restricted to situations in which it is essential.
3.2 Garbage Collecting Method
Once a Vertical Grid Descriptor object had been constucted, the allocated memory can be released by use the vgd_free method like so :
stat = vgd_free(vgd)
3.3 Query Methods
The vgd_get function allows the user to query the Vertical Grid Descriptor object for information.status = vgd_get(self,key,value) type(vgrid_descriptor), intent(in) :: self !Vertical descriptor instance character(len=*), intent(in) :: key !Descriptor key to retrieve TYPE, intent(out) :: value !Retrieved value
Where TYPE can be any of the valid data types in the object. For a complete list of Vertical Grid Descriptor keys, see the reference list.
3.4 Definition Methods
The vgd_put function allows the user to set values in the Vertical Grid Descriptor object (instance variables) after construction.status = vgd_put(self,key,value) type(vgrid_descriptor), intent(inout) :: self !Vertical descriptor instance character(len=*), intent(in) :: key !Descriptor key to set TYPE, intent(in) :: value !Value to set
Where TYPE can be any of the valid data types in the object. Under extremely special circumstances, it may be necessary to modify information in the vertical descriptor record itself (for example, its etiket in an RPN Standard File). For a complete list of Vertical Grid Descriptor keys, see the reference list
3.5 Write Methods
This method allows the user to write the current contents of the Vertical Grid Descriptor object to a file as a vertical grid descriptor record (!!).status = vgd_write(self,unit,format) type(vgrid_descriptor), intent(in) :: self !Vertical descriptor instance integer, intent(in) :: unit !File unit to write to character(len=*), intent(in) :: format !File format ('fst' or 'bin')
3.6 Position Calculations
These methods allow the user to compute the position of points in space. The outputs fields from these methods should be declared in the calling procedure as pointer arrays to allow for correctly-sized allocation within the module.
3.6.1 Vertical Levels
To compute the hydrostatic pressure (Pa).
The use of a prototype field (actually, the FST key of the prototype field) allows the user to specify which vertical coordinate should be used to compute level information. This is especially useful for staggered grids, where different variables appear at different levels with the same vertical descriptor. Note that in this case, any required reference field is read automatically from the input file.status = vgd_levels(self,unit,fstkeys,levels,in_log) type(vgrid_descriptor), intent(in) :: self !Vertical descriptor instance integer, intent(in) :: unit !File unit associated with the key integer, dimension(:), intent(in) :: fstkeys !Key of prototype field real, dimension(:,:,:), intent(out) :: levels !Physical level values logical, optional, intent(in) :: in_log !Compute levels in ln() [.false.]
If no prototype field is available, there is no input file (i.e. the vertical descriptor has been constructed manually), or a user-defined reference field is to be provided, then a second interface to the level calculation is provided. Waring: data in the sfc_field must by in Pa.
status = vgd_levels(self,ip1_list,levels,sfc_field,in_log) result(status) type(vgrid_descriptor), intent(in) :: self !Vertical descriptor instance integer, dimension(:), intent(in) :: ip1_list !Key of prototype field real, dimension(:,:,:), pointer :: levels !Physical level values real, dimension(:,:), optional, intent(in) :: sfc_field !Surface field reference for coordinate [none] logical, optional, intent(in) :: in_log !Compute levels in ln() [.false.]
status = vgd_levels(self,ip1_list,levels,sfc_field,in_log) result(status) type(vgrid_descriptor), intent(in) :: self !Vertical descriptor instance integer, dimension(:), intent(in) :: ip1_list !Key of prototype field real, dimension(:), pointer :: levels !Physical level values real, optional, intent(in) :: sfc_field !Surface field reference for coordinate [none] logical, optional, intent(in) :: in_log !Compute levels in ln() [.false.]
Note : real*8 interfaces are also provided. To access these interface, make pointer levels and sfc_field of type real*8.
3.7 hydrostatic pressure partial derivative with respect to surface hydrostatic pressure
To compute the hydrostatic pressure partial derivative with respect to surface hydrostatic pressure (more info here) use the following method:
vgd_dpidpis(self,ip1_list,dpidpis,sfc_field) result(status) type(vgrid_descriptor), intent(in) :: self !Vertical descriptor instance integer, dimension(:), intent(in) :: ip1_list !Key of prototype field real, dimension(:,:,:), pointer :: dpidpis !pressure derivative with respect to sfc pressure real, dimension(:,:), optional, intent(in) :: sfc_field !Surface field reference for coordinate [none]
A very similar interface is provided for users who wish to compute the field on profiles rather than the full three-dimensional fields. In this case, a scalar sfc_field value is provided (if required for the coordinate in question) and the returned levels is a one-dimensional real pointer.
vgd_dpidpis(self,ip1_list,dpidpis,sfc_field) result(status) type(vgrid_descriptor), intent(in) :: self !Vertical descriptor instance integer, dimension(:), intent(in) :: ip1_list !Key of prototype field real, dimension(:), pointer :: dpidpis !Derivative values real, optional, intent(in) :: sfc_field !Surface field reference for coordinate [none]
Note : real*8 interfaces are also provided. To access these interface, make pointer levels and sfc_field of type real*8.
3.8 Diagnostic Methods
These methods allow the user to obtain a dump of the information contained in the Vertical Grid Descriptor object. While not normally useful, they can be important methods during development since they give all available information about the vertical grid descriptors in use.
status = vgd_print(self,stdout) type(vgrid_descriptor), intent(in) :: self !Vertical descriptor instance integer, intent(in), optional :: stdout !Output unit to write to [6]
An alternative interface to vgd_print is provided (version 5.1.0 and up) and allows to get basic information about any Vcode
status = vgd_print(vcode) integer, intent(in) :: vcode
4 Valid Keys
A set of keys are available to allow the user to reference the instance variables (elements) of the Grid Desriptor object. These keys can be used in the vgd_get and vgd_put methods to retrieve and set values of components of the object. Each key is identified using the first four characters of the string. After that, the string may contain any text that makes the code as self-documenting as possible. In the first column of the table shown here, suggested long-form names for the keys are provided for reference. All key names are case-independent.
4.1 Standard File Record Keys
The ip1 and ip2 entries in the standard file record are used to determine which fields are associated with which set of vertical descriptors. To modify these values, use the entries in the following table.
Key Name | Type and Attributes | Description |
---|---|---|
DATE - date stamp | integer | Date stamp of descriptor |
ETIK - record stamp | character(len=*) | Stamp (etiket) of descriptor |
IG_1 - record ig1 | integer | IG1 value of descriptor |
IG_2 - record ig2 | integer | IG2 value of descriptor |
IG_3 - record ig3 | integer | IG3 value of descriptor |
IG_4 - record ig4 | integer | IG4 value of descriptor |
IP_1 - record ip1 | integer | IP1 value of vertical descriptor |
IP_2 - record ip2 | integer | IP2 value of vertical descriptor |
IP_3 - record ip3 | integer | IP3 value of vertical descriptor |
NAME - descriptor name | character(len=*) | Name of the descriptor |
status = vgd_put(d,key='IP_1 - record ip1',value=100)
4.2 Vertical Coordinate Keys
Key Name | Type and Attributes | Description |
---|---|---|
KIND - vertical coordinate ip1 kind | integer | Kind of the vertical coordinate ip1 |
VERS - vertical coordinate version | integer | For a given kind there may be many versions, example kind=5 version=2 is hyb staggered GEM4.1 |
NL_M - number of momentum levels | integer | Number of momentum levels (verison 3.2.0 and up) |
NL_T - number of thermodynamic levels | integer | Number of thermodynamic levels (version 3.2.0 and up) |
CA_M - vertical A coefficient (m) | real*8, dimension(:), pointer | Values of coefficient A on momentum levels |
CA_T - vertical A coefficient (t) | real*8, dimension(:), pointer | Values of coefficient A on thermodynamic levels |
CB_M - vertical B coefficient (m) | real*8, dimension(:), pointer | Values of coefficient B on momentum levels |
CB_T - vertical B coefficient (t) | real*8, dimension(:), pointer | Values of coefficient B on thermodynamic levels |
COFA - vertical A coefficient | real*8, dimension(:), pointer | Values of coefficient A in unstaggered levelling |
COFB - vertical B coefficient | real*8, dimension(:), pointer | Values of coefficient B in unstaggered levelling |
DIPM - IP1 of diagnostic level (m) | integer | The IP1 value of the momentum diagnostic level |
DIPT - IP1 of diagnostic level (t) | integer | The IP1 value of the thermodynamic diagnostic level |
DHM - height of the diagnostic level (m)[1] | real, dimension(:), pointer | Height of the momentum diagonstic level (m) |
DHT - height of the diagnostic level (t)[2] | real, dimension(:), pointer | Height of the thermodynamic diagonstic level (m) |
PREF - reference pressure | real*8 | Pressure of the reference level (Pa) |
PTOP - top level pressure | real*8 | Pressure of the top level (Pa) |
RC_1 - first R-coef value | real , dimension(:), pointer | First coordinate recification R-coefficient |
RC_2 - second R-coef value | real , dimension(:), pointer | Second coordinate recification R-coefficient |
RFLD - reference field name | character(len=VGD_LEN_NAME) | Name of the reference field for the vertical coordinate in the FST file |
VCDM - vertical coordinate (m)[3] | real, dimension(:), pointer | List of momentum coordinate values |
VCDT - vertical coordinate (t)[4] | real, dimension(:), pointer | List of thermodynamic coordinate values |
VIPM - level ip1 list (m) | integer, dimension(:), pointer | List of IP1 momentum values associated with this coordinate |
VIPT - level ip1 list (t) | integer, dimension(:), pointer | List of IP1 thermodynamic values associated with this coordinate |
VTBL - vgrid_descriptor table | real*8, dimension(:,:,:), pointer | real*8 Fortran array containing all vgrid_descriptor information. more info here |
LOGP - furmula gives log(p) T/F (version 1.0.3 and greater) | logical | True -> Formula with A and B gives log(p),
False -> Formula with A and B gives p |
5 Class Variables
The following class variables are available to all procedures that use-associate the Vertical Grid Descriptor module.
Class Variable Name | Description |
---|---|
VGD_ERROR | Value of return status for an error method completion. |
VGD_LEN_NAME | Length of field names for vertical grid descriptors. |
VGD_MISSING | Missing value returned from requests to retrieve invalid elements. |
VGD_OK | Value of return status for successful method completion. |
6 Options
6.1 Setting and getting options
Functions vgd_getopt and vgd_putopt allows to query and set option value, for example, option ALLOW_RESHAPE can be queried and set with the following code:
logical :: my_bol_value
...
stat = vgd_getopt("ALLOW_RESHAPE",my_bol_value)
stat = vgd_putopt("ALLOW_RESHAPE",.true.)
stat = vgd_putopt("ALLOW_RESHAPE",.false.)
6.2 ALLOW_RESHAPE
The option ALLOW_RESHAPE controls vectors/tables reshaping by the package (version 4.2.0 and up). The default value is .false., therefore, if a user allocates a vector to the incorrect length, the package will not reshape it and will return the error code VGD_ERROR and will print an error message similar to the following:
ERROR: (vgd) 3D pointer levels already allocated with a different length, will not reallocate in levels_withref
On the other hand, if the ALLOW_RESHAPE is set to .true. (see code below), the package will reshape an incorrect length vector and will print a warning.
stat = vgd_getopt("ALLOW_RESHAPE",my_bol_value)
stat = vgd_putopt("ALLOW_RESHAPE",.true.)
stat = vgd_putopt("ALLOW_RESHAPE",.false.)
- Footnotes
- ↑ This value can be retrieved with using vgd_get(), but not set since it is not a member of the vgrid_descriptor type. The DIPM key should be used to set the diagnostic level IP1 value via vgd_put() instead.
- ↑ This value can be retrieved with using vgd_get(), but not set since it is not a member of the vgrid_descriptor type. The DIPT key should be used to set the diagnostic level IP1 value via vgd_put() instead.
- ↑ Value may be retrieved using vgd_get(), but not set since it is not a member of the vgrid_descriptor type. The VIPM key should be used to set levels in vgd_put() instead.
- ↑ Value may be retrieved using vgd_get(), but not set since it is not a member of the vgrid_descriptor type. The VIPT key should be used to set levels in vgd_put() instead.