added README_changes.txt
[wrffire.git] / WPS / metgrid / src / datatype_module.f90
blob066ba9b03967e1c9ffb7711d2e3effb0de019575
1 module datatype_module
3 use bitarray_module
4 use module_debug
6 ! Return values for comparison functions primary_cmp() and secondary_cmp()
7 integer, parameter :: LESS = -1, &
8 EQUAL = 0, &
9 GREATER = 1, &
10 NOT_EQUAL = 2
12 type header_info
13 integer :: version
15 ! YYYY?MM?DD?HH?mm?ss
16 character (len=32) :: date
17 logical :: time_dependent, mask_field
19 ! Set = 0 if this is an analysis.
20 real :: forecast_hour
22 ! AVN, GFS, ETA???, ARW, NMM, AGRMET, NAM, RUC, SST
23 character (len=32) :: fg_source
25 character (len=128) :: field
26 character (len=128) :: units
27 character (len=128) :: description
29 ! PRESSURE, SIGMA, NATIVE, HYBRID
30 character (len=32) :: vertical_coord
31 integer :: vertical_level
33 ! XY, YX - ENOUGH INFO?
34 character (len=32) :: array_order
35 integer, dimension(2) :: dim1, dim2
37 logical :: is_wind_grid_rel
38 logical :: array_has_missing_values
40 integer :: sr_x, sr_y
41 end type header_info
43 type map_info
44 ! Mercator, Polar Stereographic, Lambert, Gaussian, Lat Lon
45 character (len=32) :: projection
47 integer :: projection_flag
49 ! For ARW: M, U, or V; for NMM: H or V
50 integer :: stagger
52 real :: knownlat, knownlon, deltalat, deltalon
53 real :: deltax, deltay, xlonc, truelat1, truelat2
54 real :: lat1, lon1
55 end type map_info
57 ! This is the datatype that is understood by data_storage module
58 type fg_input
59 ! BEGIN any types we want to keep and use for sorting in storage module
60 type (header_info) :: header
61 type (map_info) :: map
62 ! END any types we want to keep and use for sorting in storage module
64 real, dimension(:,:), pointer :: r_arr !!!!! REQUIRED !!!!!
65 type (bitarray), pointer :: valid_mask, modified_mask
66 end type fg_input
68 ! This type is used for the nodes of the secondary linked lists, the ones that
69 ! actually store data
70 type data_node
71 type (fg_input) :: fg_data
73 type (data_node), pointer :: next, prev
74 integer, dimension(2) :: field_shape
76 ! If non-zero, the array is actually stored in a Fortran unit
77 integer :: filenumber
79 ! The following two are used by heaps
80 integer :: last_used
81 integer :: heap_index
82 end type data_node
84 ! This type is used for the nodes in the primary linked lists, and thus has head
85 ! and tail pointers for secondary linked lists
86 type head_node
87 type (fg_input) :: fg_data
88 type (head_node), pointer :: next, prev
89 type (data_node), pointer :: fieldlist_head, fieldlist_tail
90 end type head_node
93 contains
96 ! Compares two fg_input types; returns EQUAL if the two should
97 ! belong to the same secondary linked list, and NOT_EQUAL otherwise
98 function primary_cmp(a, b)
100 implicit none
102 ! Arguments
103 type (fg_input), intent(in) :: a, b
105 ! Return value
106 integer :: primary_cmp
108 ! if ((a%header%date == b%header%date) .and. &
109 ! (a%header%forecast_hour == b%header%forecast_hour) .and. &
110 ! (a%header%fg_source == b%header%fg_source) .and. &
111 ! (a%header%field == b%header%field)) then
112 if (a%header%field == b%header%field) then
113 primary_cmp = EQUAL
114 else
115 primary_cmp = NOT_EQUAL
116 end if
118 end function primary_cmp
121 ! Compares two fg_input types; returns EQUAL if the two belong
122 ! at the same position in a secondary linked list, LESS if "a" belongs
123 ! after "b", and GREATER if "a" belongs before "b"
124 function secondary_cmp(a, b)
126 implicit none
128 ! Arguments
129 type (fg_input), intent(in) :: a, b
131 ! Return value
132 integer :: secondary_cmp
134 ! BUG: Eventually, we only want to sort pressure-level data this way, and
135 ! all others the opposite way, as in the else case below.
136 if (a%header%time_dependent) then
137 if (a%header%vertical_level > b%header%vertical_level) then
138 secondary_cmp = LESS
139 else if (a%header%vertical_level == b%header%vertical_level) then
140 secondary_cmp = EQUAL
141 else
142 secondary_cmp = GREATER
143 end if
145 else
146 if (a%header%vertical_level < b%header%vertical_level) then
147 secondary_cmp = LESS
148 else if (a%header%vertical_level == b%header%vertical_level) then
149 secondary_cmp = EQUAL
150 else
151 secondary_cmp = GREATER
152 end if
153 end if
155 end function secondary_cmp
158 ! Duplicates an fg_input type
159 subroutine dup(src, dst)
161 implicit none
163 ! Arguments
164 type (fg_input), intent(in) :: src
165 type (fg_input), intent(out) :: dst
167 dst%header = src%header
168 dst%map = src%map
169 dst%r_arr => src%r_arr
170 dst%valid_mask => src%valid_mask
171 dst%modified_mask => src%modified_mask
173 end subroutine dup
176 function is_time_dependent(a)
178 implicit none
180 ! Arguments
181 type (fg_input), intent(in) :: a
183 ! Return value
184 logical :: is_time_dependent
186 is_time_dependent = a%header%time_dependent
188 end function is_time_dependent
191 function is_mask_field(a)
193 implicit none
195 ! Arguments
196 type (fg_input), intent(in) :: a
198 ! Return value
199 logical :: is_mask_field
201 is_mask_field = a%header%mask_field
203 end function is_mask_field
206 ! Returns the vertical level of an fg_input type
207 function get_level(a)
209 implicit none
211 ! Arguments
212 type (fg_input), intent(in) :: a
214 ! Return value
215 integer :: get_level
217 get_level = a%header%vertical_level
219 end function get_level
222 ! Returns the description string of an fg_input type
223 function get_description(a)
225 implicit none
227 ! Arguments
228 type (fg_input), intent(in) :: a
230 ! Return value
231 character (len=128) :: get_description
233 get_description = a%header%description
235 end function get_description
238 ! Returns the units string of an fg_input type
239 function get_units(a)
241 implicit none
243 ! Arguments
244 type (fg_input), intent(in) :: a
246 ! Return value
247 character (len=128) :: get_units
249 get_units = a%header%units
251 end function get_units
254 ! Returns the field staggering an fg_input type
255 function get_staggering(a)
257 implicit none
259 ! Arguments
260 type (fg_input), intent(in) :: a
262 ! Return value
263 integer :: get_staggering
265 get_staggering = a%map%stagger
267 end function get_staggering
270 ! Returns the fieldname string of an fg_input type
271 function get_fieldname(a)
273 implicit none
275 ! Arguments
276 type (fg_input), intent(in) :: a
278 ! Return value
279 character (len=128) :: get_fieldname
281 get_fieldname = a%header%field
283 end function get_fieldname
286 ! Gives starting and ending indices for a field
287 subroutine get_dims(a, start_mem_1, end_mem_1, start_mem_2, end_mem_2)
289 implicit none
291 ! Arguments
292 type (fg_input), intent(in) :: a
293 integer, intent(out) :: start_mem_1, end_mem_1, start_mem_2, end_mem_2
295 start_mem_1 = a%header%dim1(1)
296 end_mem_1 = a%header%dim1(2)
297 start_mem_2 = a%header%dim2(1)
298 end_mem_2 = a%header%dim2(2)
300 end subroutine get_dims
303 ! Prints relevant information from the headers of an fg_input type; mainly
304 ! used for debugging
305 subroutine print_header(a)
307 implicit none
309 ! Arguments
310 type (fg_input), intent(in) :: a
312 call mprintf(.true.,DEBUG,'FIELD : %s',s1=trim(a%header%field))
313 call mprintf(.true.,DEBUG,'DATE : %s',s1=trim(a%header%date))
314 call mprintf(.true.,DEBUG,'SOURCE : %s',s1=trim(a%header%fg_source))
315 call mprintf(.true.,DEBUG,'FCST HR: %f',f1=a%header%forecast_hour)
317 end subroutine print_header
319 end module datatype_module