1 ! Module for handling all associated scale_factor, add_offset, and
2 ! attributes for individual files across various possible
3 ! National Water Model output files. In the future, this will move to a
4 ! table that the user will be able to switch on/off variables for
5 ! outputting. For now, attributes, etc will be stored here.
8 ! National Center for Atmospheric Research
9 ! Research Applications Laboratory
10 ! karsten at ucar dot edu
12 module module_NWM_io_dict
14 use module_version, only: get_model_version
15 use config_base, only: nlst
16 use module_hydro_stop, only: HYDRO_stop
20 ! Declare parameter values for module.
21 integer, parameter :: numChVars = 11
22 integer, parameter :: numLdasVars = 116
23 ! Note: if more ldas variables are added the logic will need to be changed in
24 ! module_NWM_io.F:output_NoahMP_NWM for when to close the restart file
25 integer, parameter :: numLdasVars_crocus_off = 98
27 integer, parameter :: numRtDomainVars = 5
28 integer, parameter :: numLakeVars = 2
29 integer, parameter :: numChGrdVars = 1
30 integer, parameter :: numLsmVars = 14
31 integer, parameter :: numChObsVars = 1
32 integer, parameter :: numGwVars = 4
34 !integer :: nsoil = nlst_rt(1)%nsoil
37 ! Declare public types that will hold metadata
38 public :: chrtMeta ! Public CHRTOUT metadata for NWM output.
39 public :: ldasMeta ! Public LDASOUT metadata for NWM output.
40 public :: rtDomainMeta ! Public RT_DOMAIN metadata for NWM output.
41 public :: lakeMeta ! Public lake metadata for NWM output
42 public :: chrtGrdMeta ! Public CHRTOUT_GRID metadata for NWM output.
43 public :: lsmMeta ! Public metadata for LSMOUT output.
44 public :: chObsMeta ! Public CHANOBS metadata for NWM output.
45 public :: gwMeta ! Public groundwater metadata.
47 real*8 :: one_dbl = 1.0d0
49 ! Establish types for each output type
52 character (len=64), dimension(numChVars) :: varNames
53 integer :: numVars = numChVars
54 character (len=512) :: modelOutputType = "channel_rt"
55 character (len=64) :: modelConfigType
56 ! Output variable attributes
57 real, dimension(numChVars) :: scaleFactor ! scale_factor values used for each
58 ! variable to converte from real to
60 real, dimension(numChVars) :: addOffset ! add_offset values for each variable.
61 character (len=64), dimension(numChVars) :: longName ! Long names for each variable.
62 character (len=64), dimension(numChVars) :: units ! Units for each variable.
63 character (len=64), dimension(numChVars) :: coordNames ! Coordinate names for each variable.
64 integer(kind=4), dimension(numChVars) :: validMinComp ! Valid min (after conversion to integer)
65 integer(kind=4), dimension(numChVars) :: validMaxComp ! Valid max (after conversion to integer)
66 real*8, dimension(numChVars) :: validMinDbl ! Valid minimum (before conversion to integer)
67 real*8, dimension(numChVars) :: validMaxDbl ! Valid maximum (before converstion to integer)
68 integer(kind=4), dimension(numChVars) :: missingComp ! Missing value attribute (after conversion to integer)
69 real, dimension(numChVars) :: missingReal ! Missing value attribute (before conversion to integer)
70 real, dimension(numChVars) :: fillReal ! Fill value (before conversion to integer)
71 integer(kind=4), dimension(numChVars) :: fillComp ! Fill value (after conversion to integer)
72 integer, dimension(numChVars) :: outFlag ! 0/1 flag to turn outputting off/on
73 integer, dimension(numChVars) :: timeZeroFlag ! 0/1 flag to either set variable to all NDV values or output
74 ! the actual data. This was done because time 0
75 ! output does not all contain valid data.
76 real :: modelNdv ! NDV value represented within the model code
77 ! Time variable attribues
78 character (len=64) :: timeLName ! long_name - usually valid output time
79 character (len=64) :: timeUnits ! Usually seconds since 1/1/1970
80 character (len=64) :: timeStName ! standard_name - usually time
81 integer :: timeValidMin ! the minimum time each configuration can have, time of the first output file
82 integer :: timeValidMax ! the maximum time each configuration can have, time of the last output file
83 ! Reference time attributes
84 character (len=64) :: rTimeLName ! long_name - usually model initialization time
85 character (len=64) :: rTimeStName ! standard_name - usually forecast_reference_time
86 character (len=64) :: rTimeUnits ! Usually seconds since 1/1/1970
87 ! feature_id attributes
88 character (len=256) :: featureIdLName ! long_name - usually Reach ID
89 character (len=256) :: featureIdComment ! Comment attribute
90 character (len=64) :: cfRole ! cf_role attribute
91 ! latitude variable attributes
92 character (len=64) :: latLName ! long_name
93 character (len=64) :: latUnits ! units
94 character (len=64) :: latStName ! Standard Name
95 ! longitude variable attributes
96 character (len=64) :: lonLName ! long_name
97 character (len=64) :: lonUnits ! units
98 character (len=64) :: lonStName ! Standard Name
99 ! Elevation variable attributes
100 character (len=64) :: elevLName ! long_name
101 character (len=64) :: elevUnits ! units
102 character (len=64) :: elevStName ! Standard Name
104 ! Order variable attributes
105 character (len=64) :: orderLName ! long_name
106 character (len=64) :: orderStName ! Standard Name
108 character (len=128) :: title ! file title
109 character (len=128) :: fType ! featureType attribute
110 character (len=128) :: proj4 ! proj4 attribute
111 character (len=128) :: initTime ! model_initialization_time attribute
112 character (len=128) :: validTime ! model_output_valid_time attribute
113 character (len=128) :: stDim ! station_dimension attribute
114 integer :: stOrder ! stream_order_output attribute
115 character (len=128) :: cdm ! cdm_datatype attribute
116 character (len=1024) :: esri ! esri_pe_string attribute
117 character (len=128) :: conventions ! Conventions string
118 integer :: totalValidTime ! # number of valid time (#of output files)
125 character (len=64), dimension(numLdasVars) :: varNames
126 integer :: numVars = numLdasVars
127 integer :: numSnowLayers = 3
128 integer :: numSoilLayers ! Fill from Namelist
130 integer :: numSpectrumBands = 2
131 character (len=512) :: modelOutputType = "land"
132 character (len=64) :: modelConfigType
134 ! Output variable attributes
135 real, dimension(numLdasVars) :: scaleFactor ! scale_factor values used for each
136 ! variable to converte from
139 real, dimension(numLdasVars) :: addOffset ! add_offset values for each variable.
140 character (len=128), dimension(numLdasVars) :: longName ! Long names for each variable.
141 character (len=64), dimension(numLdasVars) :: units ! Units for each variable.
142 integer, dimension(numLdasVars) :: numLev ! Number of levels for each variable.
143 integer(kind=4), dimension(numLdasVars) :: validMinComp ! Valid min (after conversion to integer)
144 integer(kind=4), dimension(numLdasVars) :: validMaxComp ! Valid max (after conversion to integer)
145 real*8, dimension(numLdasVars) :: validMinDbl ! Valid minimum (before conversion to integer)
146 real*8, dimension(numLdasVars) :: validMaxDbl ! Valid maximum (before converstion to integer)
147 integer(kind=4), dimension(numLdasVars) :: missingComp ! Missing value attribute (after conversion to integer)
148 real, dimension(numLdasVars) :: missingReal ! Missing value attribute (before conversion to integer)
149 real, dimension(numLdasVars) :: fillReal ! Fill value (before conversion to integer)
150 integer(kind=4), dimension(numLdasVars) :: fillComp ! Fill value (after conversion to integer)
151 integer, dimension(numLdasVars) :: outFlag ! 0/1 flag to turn outputting off/on
152 integer, dimension(numLdasVars) :: timeZeroFlag ! 0/1 flag to either set variable to all NDV values or output
153 ! the actual data. This was done because time 0
154 ! output does not all contain valid data.
155 real :: modelNdv ! NDV value represented within the model code
156 real :: modelNdv2 ! Alternative NDV value in NoahMP
157 real :: modelNdv3 ! Alternative NDV value in NoahMP
158 integer :: modelNdvInt ! NDV value represented in model as integer
159 ! Time variable attribues
160 character (len=64) :: timeLName ! long_name - usually valid output time
161 character (len=64) :: timeUnits ! Usually seconds since 1/1/1970
162 character (len=64) :: timeStName ! standard_name - usually time
163 integer :: timeValidMin ! the minimum time each configuration can have, time of the first output file
164 integer :: timeValidMax ! the maximum time each configuration can have, time of the last output file
166 ! Reference time attributes
167 character (len=64) :: rTimeLName ! long_name - usually model initialization time
168 character (len=64) :: rTimeStName ! standard_name - usually forecast_reference_time
169 character (len=64) :: rTimeUnits ! Usually seconds since 1/1/1970
170 ! Establish a proj4 string
171 character (len=1024) :: proj4
172 ! Establish crs-related variables.
173 character(len=1024), dimension(20) :: crsCharAttNames,crsFloatAttNames
174 character(len=1024), dimension(20) :: crsCharAttVals
175 real, dimension(20,20) :: crsRealAttVals
176 integer, dimension(20) :: crsRealAttLen
177 integer :: nCrsRealAtts, nCrsCharAtts
178 ! Establish x/y related variables.
179 character(len=1024), dimension(20) :: xCharAttNames,xFloatAttNames
180 character(len=1024), dimension(20) :: xCharAttVals
181 real, dimension(20,20) :: xRealAttVals
182 integer, dimension(20) :: xRealAttLen
183 integer :: nxRealAtts, nxCharAtts
184 character(len=1024), dimension(20) :: yCharAttNames,yFloatAttNames
185 character(len=1024), dimension(20) :: yCharAttVals
186 real, dimension(20,20) :: yRealAttVals
187 integer, dimension(20) :: yRealAttLen
188 integer :: nyRealAtts, nyCharAtts
191 character (len=128) :: title ! file title
192 character (len=128) :: initTime ! model_initialization_time attribute
193 character (len=128) :: validTime ! model_output_valid_time attribute
194 character (len=128) :: conventions ! Conventions string
195 integer :: totalValidTime ! # number of valid time (#of output files)
201 character (len=64), dimension(numRtDomainVars) :: varNames
202 integer :: numVars = numRtDomainVars
203 integer :: numSoilLayers ! Fill from Namelist
204 character (len=512) :: modelOutputType = "terrain_rt"
205 character (len=64) :: modelConfigType
207 ! Output variable attributes
208 real, dimension(numRtDomainVars) :: scaleFactor ! scale_factor values used for each
209 ! variable to converte from
212 real, dimension(numRtDomainVars) :: addOffset ! add_offset values for each variable.
213 character (len=64), dimension(numRtDomainVars) :: longName ! Long names for each variable.
214 character (len=64), dimension(numRtDomainVars) :: units ! Units for each variable.
215 character (len=64), dimension(numRtDomainVars) :: coordNames ! Coordinate names for each variable.
216 integer(kind=4), dimension(numRtDomainVars) :: validMinComp ! Valid min (after conversion to integer)
217 integer(kind=4), dimension(numRtDomainVars) :: validMaxComp ! Valid max (after conversion to integer)
218 real*8, dimension(numRtDomainVars) :: validMinDbl ! Valid minimum (before conversion to integer)
219 real*8, dimension(numRtDomainVars) :: validMaxDbl ! Valid maximum (before converstion to integer)
220 integer(kind=4), dimension(numRtDomainVars) :: missingComp ! Missing value attribute (after conversion to integer)
221 real, dimension(numRtDomainVars) :: missingReal ! Missing value attribute (before conversion to integer)
222 real, dimension(numRtDomainVars) :: fillReal ! Fill value (before conversion to integer)
223 integer(kind=4), dimension(numRtDomainVars) :: fillComp ! Fill value (after conversion to integer)
224 integer, dimension(numRtDomainVars) :: outFlag ! 0/1 flag to turn outputting off/on
225 integer, dimension(numRtDomainVars) :: timeZeroFlag ! 0/1 flag to either set variable to all NDV values or output
226 ! the actual data. This was done because time 0
227 ! output does not all contain valid data.
228 real :: modelNdv ! NDV value represented within the model code
229 ! Time variable attribues
230 character (len=64) :: timeLName ! long_name - usually valid output time
231 character (len=64) :: timeUnits ! Usually seconds since 1/1/1970
232 character (len=64) :: timeStName ! standard_name - usually time
233 integer :: timeValidMin ! the minimum time each configuration can have, time of the first output file
234 integer :: timeValidMax ! the maximum time each configuration can have, time of the last output file
235 ! Reference time attributes
236 character (len=64) :: rTimeLName ! long_name - usually model initialization time
237 character (len=64) :: rTimeStName ! standard_name - usually forecast_reference_time
238 character (len=64) :: rTimeUnits ! Usually seconds since 1/1/1970
239 real :: yRes, xRes ! Resoltution in meters
241 ! Establish a proj4 string
242 character (len=1024) :: proj4
243 ! Establish crs-related variables.
244 character(len=1024), dimension(20) :: crsCharAttNames,crsFloatAttNames
245 character(len=1024), dimension(20) :: crsCharAttVals
246 real, dimension(20,20) :: crsRealAttVals
247 integer, dimension(20) :: crsRealAttLen
248 integer :: nCrsRealAtts, nCrsCharAtts
249 ! Establish x/y related variables.
250 character(len=1024), dimension(20) :: xCharAttNames,xFloatAttNames
251 character(len=1024), dimension(20) :: xCharAttVals
252 real, dimension(20,20) :: xRealAttVals
253 integer, dimension(20) :: xRealAttLen
254 integer :: nxRealAtts, nxCharAtts
255 character(len=1024), dimension(20) :: yCharAttNames,yFloatAttNames
256 character(len=1024), dimension(20) :: yCharAttVals
257 real, dimension(20,20) :: yRealAttVals
258 integer, dimension(20) :: yRealAttLen
259 integer :: nyRealAtts, nyCharAtts
262 integer :: decimation ! Decimation factor
263 character (len=128) :: title ! file title
264 character (len=128) :: initTime ! model_initialization_time attribute
265 character (len=128) :: validTime ! model_output_valid_time attribute
266 character (len=128) :: conventions ! Conventions string
267 integer :: totalValidTime ! # number of valid time (#of output files)
269 end type rtDomainMeta
273 character (len=64), dimension(numLakeVars) :: varNames
274 integer :: numVars = numLakeVars
275 character (len=512) :: modelOutputType = "reservoir"
276 character (len=64) :: modelConfigType
278 ! Output variable attributes
279 real, dimension(numLakeVars) :: scaleFactor ! scale_factor values used for each
280 ! variable to converte from
283 real, dimension(numLakeVars) :: addOffset ! add_offset values for each variable.
284 character (len=64), dimension(numLakeVars) :: longName ! Long names for each variable.
285 character (len=64), dimension(numLakeVars) :: units ! Units for each variable.
286 character (len=64), dimension(numLakeVars) :: coordNames ! Coordinate names for each variable.
287 integer(kind=4), dimension(numLakeVars) :: validMinComp ! Valid min (after conversion to integer)
288 integer(kind=4), dimension(numLakeVars) :: validMaxComp ! Valid max (after conversion to integer)
289 real*8, dimension(numLakeVars) :: validMinDbl ! Valid minimum (before conversion to integer)
290 real*8, dimension(numLakeVars) :: validMaxDbl ! Valid maximum (before converstion to integer)
291 integer(kind=4), dimension(numLakeVars) :: missingComp ! Missing value attribute (after conversion to integer)
292 real, dimension(numLakeVars) :: missingReal ! Missing value attribute (before conversion to integer)
293 real, dimension(numLakeVars) :: fillReal ! Fill value (before conversion to integer)
294 integer(kind=4), dimension(numLakeVars) :: fillComp ! Fill value (after conversion to integer)
295 integer, dimension(numLakeVars) :: outFlag ! 0/1 flag to turn outputting off/on
296 integer, dimension(numLakeVars) :: timeZeroFlag ! 0/1 flag to either set variable to all NDV values or output
297 ! the actual data. This was done because time 0
298 ! output does not all contain valid data.
299 real :: modelNdv ! NDV value represented within the model code
300 ! Time variable attribues
301 character (len=64) :: timeLName ! long_name - usually valid output time
302 character (len=64) :: timeUnits ! Usually seconds since 1/1/1970
303 character (len=64) :: timeStName ! standard_name - usually time
304 integer :: timeValidMin ! the minimum time each configuration can have, time of the first output file
305 integer :: timeValidMax ! the maximum time each configuration can have, time of the last output file
306 ! Reference time attributes
307 character (len=64) :: rTimeLName ! long_name - usually model initialization time
308 character (len=64) :: rTimeStName ! standard_name - usually forecast_reference_time
309 character (len=64) :: rTimeUnits ! Usually seconds since 1/1/1970
311 character (len=64) :: lakeIdLName ! long_name - usually Lake COMMON ID
312 character (len=256) :: LakeIdComment ! Comment attribute
313 ! feature_id attributes
314 character (len=256) :: featureIdLName ! long_name - usually lake COMMON ID
315 character (len=256) :: featureIdComment ! Comment attribute
316 character (len=64) :: cfRole ! cf_role attribute
317 ! reservoir_type attributes
318 character (len=64) :: reservoirTypeLName ! long_name - usually "reservoir_type"
319 integer, dimension(4) :: reservoirTypeFlagValues ! valid flags attribute
320 character (len=64) :: reservoirTypeFlagMeanings ! flag meanings attribute
321 ! reservoir_assimilated_value attributes
322 character (len=64) :: reservoirAssimilatedValueLName ! long_name - usually "reservoir_assimilated_value, m3 s-1"
323 character (len=64) :: reservoirAssimilatedValueUnits ! units
324 ! reservoir_assimilated_source_file attributes
325 character (len=64) :: reservoirAssimilatedSourceFileLName ! long_name - usually "reservoir_assimilated_source_file"
326 ! latitude variable attributes
327 character (len=64) :: latLName ! long_name
328 character (len=64) :: latUnits ! units
329 character (len=64) :: latStName ! Standard Name
330 ! longitude variable attributes
331 character (len=64) :: lonLName ! long_name
332 character (len=64) :: lonUnits ! units
333 character (len=64) :: lonStName ! Standard Name
334 ! elevation variable attributes
335 character (len=64) :: elevLName ! long_name
336 character (len=64) :: elevUnits ! units
337 character (len=64) :: elevStName ! Standard Name
338 character (len=256) :: elevComment ! Comment
340 character (len=128) :: title ! file title
341 character (len=128) :: fType ! featureType attribute
342 character (len=128) :: proj4 ! proj4 attribute
343 character (len=128) :: initTime ! model_initialization_time attribute
344 character (len=128) :: validTime ! model_output_valid_time attribute
345 character (len=128) :: lakeDim ! lake_dimension attribute
346 character (len=128) :: cdm ! cdm_datatype attribute
347 character (len=1024) :: esri ! esri_pe_string attribute
348 character (len=128) :: conventions ! Conventions string
349 integer :: totalValidTime ! # number of valid time (#of output files)
355 character (len=64), dimension(numChGrdVars) :: varNames
356 integer :: numVars = numChGrdVars
357 character (len=512) :: modelOutputType = "channel_rt"
358 character (len=64) :: modelConfigType
360 ! Output variable attributes
361 real, dimension(numChGrdVars) :: scaleFactor ! scale_factor values for each
362 ! variable to convert from
364 real, dimension(numChGrdVars) :: addOffset ! add_offset values for each variable.
365 character (len=64), dimension(numChGrdVars) :: longName ! longname for each variable
366 character (len=64), dimension(numChGrdVars) :: units ! Units for each variable.
367 character (len=64), dimension(numChGrdVars) :: coordNames ! Coordinate names for each variable.
368 integer(kind=4), dimension(numChGrdVars) :: validMinComp ! Valid min (after conversion to integer)
369 integer(kind=4), dimension(numChGrdVars) :: validMaxComp ! Valid max (after conversion to integer)
370 real*8, dimension(numChGrdVars) :: validMinDbl ! Valid minimum (before conversion to integer)
371 real*8, dimension(numChGrdVars) :: validMaxDbl ! Valid maximum (before converstion to integer)
372 integer(kind=4), dimension(numChGrdVars) :: missingComp ! Missing value attribute (after conversion to integer)
373 real, dimension(numChGrdVars) :: missingReal ! Missing value attribute (before conversion to integer)
374 real, dimension(numChGrdVars) :: fillReal ! Fill value (before conversion to integer)
375 integer(kind=4), dimension(numChGrdVars) :: fillComp ! Fill value (after conversion to integer)
376 integer, dimension(numChGrdVars) :: outFlag ! 0/1 flag to turn outputting off/on
377 integer, dimension(numChGrdVars) :: timeZeroFlag ! 0/1 flag to either set variable to all NDV values or output
378 ! the actual data. This
379 ! was done because time 0
380 ! output does not all
381 ! contain valid data.
382 real :: modelNdv ! NDV value represented within the model code
383 ! Time variable attribues
384 character (len=64) :: timeLName ! long_name - usually valid output time
385 character (len=64) :: timeUnits ! Usually seconds since 1/1/1970
386 character (len=64) :: timeStName ! standard_name - usually time
387 integer :: timeValidMin ! the minimum time each configuration can have, time of the first output file
388 integer :: timeValidMax ! the maximum time each configuration can have, time of the last output file
389 ! Reference time attributes
390 character (len=64) :: rTimeLName ! long_name - usually model initialization time
391 character (len=64) :: rTimeStName ! standard_name - usually forecast_reference_time
392 character (len=64) :: rTimeUnits ! Usually seconds since 1/1/1970
393 real :: yRes, xRes ! Resoltution in meters
395 ! Establish a proj4 string
396 character (len=1024) :: proj4
397 ! Establish crs-related variables.
398 character(len=1024), dimension(20) :: crsCharAttNames,crsFloatAttNames
399 character(len=1024), dimension(20) :: crsCharAttVals
400 real, dimension(20,20) :: crsRealAttVals
401 integer, dimension(20) :: crsRealAttLen
402 integer :: nCrsRealAtts, nCrsCharAtts
403 ! Establish x/y related variables.
404 character(len=1024), dimension(20) :: xCharAttNames,xFloatAttNames
405 character(len=1024), dimension(20) :: xCharAttVals
406 real, dimension(20,20) :: xRealAttVals
407 integer, dimension(20) :: xRealAttLen
408 integer :: nxRealAtts, nxCharAtts
409 character(len=1024), dimension(20) :: yCharAttNames,yFloatAttNames
410 character(len=1024), dimension(20) :: yCharAttVals
411 real, dimension(20,20) :: yRealAttVals
412 integer, dimension(20) :: yRealAttLen
413 integer :: nyRealAtts, nyCharAtts
416 integer :: decimation ! Decimation factor
417 character (len=128) :: title ! file title
418 character (len=128) :: initTime ! model_initialization_time attribute
419 character (len=128) :: validTime ! model_output_valid_time attribute
420 character (len=128) :: conventions ! Conventions string
421 integer :: totalValidTime ! # number of valid time (#of output files)
427 character (len=64), dimension(numLsmVars) :: varNames
428 integer :: numVars = numLsmVars
429 integer :: numSnowLayers = 3
430 integer :: numSoilLayers ! Fill from Namelist
432 character (len=512) :: modelOutputType = "land"
433 character (len=64) :: modelConfigType
435 ! Output variable attributes
436 real, dimension(numLsmVars) :: scaleFactor ! scale_factor values used for each
437 ! variable to converte from
440 real, dimension(numLsmVars) :: addOffset ! add_offset values for each variable.
441 character (len=64), dimension(numLsmVars) :: longName ! Long names for each variable.
442 character (len=64), dimension(numLsmVars) :: units ! Units for each variable.
443 integer, dimension(numLsmVars) :: numLev ! Number of levels for each variable.
444 integer(kind=4), dimension(numLsmVars) :: validMinComp ! Valid min (after conversion to integer)
445 integer(kind=4), dimension(numLsmVars) :: validMaxComp ! Valid max (after conversion to integer)
446 real*8, dimension(numLsmVars) :: validMinDbl ! Valid minimum (before conversion to integer)
447 real*8, dimension(numLsmVars) :: validMaxDbl ! Valid maximum (before converstion to integer)
448 integer(kind=4), dimension(numLsmVars) :: missingComp ! Missing value attribute (after conversion to integer)
449 real, dimension(numLsmVars) :: missingReal ! Missing value attribute (before conversion to integer)
450 real, dimension(numLsmVars) :: fillReal ! Fill value (before conversion to integer)
451 integer(kind=4), dimension(numLsmVars) :: fillComp ! Fill value (after conversion to integer)
452 integer, dimension(numLsmVars) :: outFlag ! 0/1 flag to turn outputting off/on
453 integer, dimension(numLsmVars) :: timeZeroFlag ! 0/1 flag to either set variable to all NDV values or output
454 ! the actual data. This was
455 ! done because time 0
456 ! output does not all contain
458 real :: modelNdv ! NDV value represented within the model code
459 integer :: modelNdvInt ! NDV value represented in model as integer
460 ! Time variable attribues
461 character (len=64) :: timeLName ! long_name - usually valid output time
462 character (len=64) :: timeUnits ! Usually seconds since 1/1/1970
463 character (len=64) :: timeStName ! standard_name - usually time
464 integer :: timeValidMin ! the minimum time each configuration can have, time of the first output file
465 integer :: timeValidMax ! the maximum time each configuration can have, time of the last output file
466 ! Reference time attributes
467 character (len=64) :: rTimeLName ! long_name - usually model initialization time
468 character (len=64) :: rTimeStName ! standard_name - usually forecast_reference_time
469 character (len=64) :: rTimeUnits ! Usually seconds since 1/1/1970
471 ! Establish a proj4 string
472 character (len=1024) :: proj4
473 ! Establish crs-related variables.
474 character(len=1024), dimension(20) :: crsCharAttNames,crsFloatAttNames
475 character(len=1024), dimension(20) :: crsCharAttVals
476 real, dimension(20,20) :: crsRealAttVals
477 integer, dimension(20) :: crsRealAttLen
478 integer :: nCrsRealAtts, nCrsCharAtts
479 ! Establish x/y related variables.
480 character(len=1024), dimension(20) :: xCharAttNames,xFloatAttNames
481 character(len=1024), dimension(20) :: xCharAttVals
482 real, dimension(20,20) :: xRealAttVals
483 integer, dimension(20) :: xRealAttLen
484 integer :: nxRealAtts, nxCharAtts
485 character(len=1024), dimension(20) :: yCharAttNames,yFloatAttNames
486 character(len=1024), dimension(20) :: yCharAttVals
487 real, dimension(20,20) :: yRealAttVals
488 integer, dimension(20) :: yRealAttLen
489 integer :: nyRealAtts, nyCharAtts
492 character (len=128) :: title ! file title
493 character (len=128) :: initTime ! model_initialization_time attribute
494 character (len=128) :: validTime ! model_output_valid_time attribute
495 character (len=128) :: conventions ! Conventions string
496 integer :: totalValidTime ! # number of valid time (#of output files)
502 character (len=64), dimension(numChObsVars) :: varNames
503 integer :: numVars = numChObsVars
504 character (len=512) :: modelOutputType = "channel_rt"
505 character (len=64) :: modelConfigType
507 ! Output variable attributes
508 real, dimension(numChObsVars) :: scaleFactor ! scale_factor values used for
509 ! each variable to convert
510 ! from real to integer
511 real, dimension(numChObsVars) :: addOffset ! add_offset values for each variable.
512 character (len=64), dimension(numChObsVars) :: longName ! Long names for each variable.
513 character (len=64), dimension(numChObsVars) :: units ! Units for each variable.
514 character (len=64), dimension(numChObsVars) :: coordNames ! Coordinate names for each variable.
515 integer(kind=4), dimension(numChObsVars) :: validMinComp ! Valid min (after conversion to integer)
516 integer(kind=4), dimension(numChObsVars) :: validMaxComp ! Valid max (after conversion to integer)
517 real*8, dimension(numChObsVars) :: validMinDbl ! Valid minimum (before conversion to integer)
518 real*8, dimension(numChObsVars) :: validMaxDbl ! Valid maximum (before converstion to integer)
519 integer(kind=4), dimension(numChObsVars) :: missingComp ! Missing value attribute (after conversion to integer)
520 real, dimension(numChObsVars) :: missingReal ! Missing value attribute (before conversion to integer)
521 real, dimension(numChObsVars) :: fillReal ! Fill value (before conversion to integer)
522 integer(kind=4), dimension(numChObsVars) :: fillComp ! Fill value (after conversion to integer)
523 integer, dimension(numChObsVars) :: outFlag ! 0/1 flag to turn outputting off/on
524 integer, dimension(numChObsVars) :: timeZeroFlag ! 0/1 flag to either set variable to all NDV values or output
525 ! the actual data. This was
526 ! done because time 0
527 ! output does not all contain valid data.
528 real :: modelNdv ! NDV value represented within the model code
529 ! Time variable attribues
530 character (len=64) :: timeLName ! long_name - usually valid output time
531 character (len=64) :: timeUnits ! Usually seconds since 1/1/1970
532 character (len=64) :: timeStName ! standard_name - usually time
533 integer :: timeValidMin ! the minimum time each configuration can have, time of the first output file
534 integer :: timeValidMax ! the maximum time each configuration can have, time of the last output file
535 ! Reference time attributes
536 character (len=64) :: rTimeLName ! long_name - usually model initialization time
537 character (len=64) :: rTimeStName ! standard_name - usually forecast_reference_time
538 character (len=64) :: rTimeUnits ! Usually seconds since 1/1/1970
539 ! feature_id attributes
540 character (len=256) :: featureIdLName ! long_name - usually Reach ID
541 character (len=256) :: featureIdComment ! Comment attribute
542 character (len=64) :: cfRole ! cf_role attribute
543 ! latitude variable attributes
544 character (len=64) :: latLName ! long_name
545 character (len=64) :: latUnits ! units
546 character (len=64) :: latStName ! Standard Name
547 ! longitude variable attributes
548 character (len=64) :: lonLName ! long_name
549 character (len=64) :: lonUnits ! units
550 character (len=64) :: lonStName ! Standard Name
551 ! Elevation variable attributes
552 character (len=64) :: elevLName ! long_name
553 character (len=64) :: elevUnits ! units
554 character (len=64) :: elevStName ! Standard Name
555 ! Order variable attributes
556 character (len=64) :: orderLName ! long_name
557 character (len=64) :: orderStName ! Standard Name
559 character (len=128) :: title ! file title
560 character (len=128) :: fType ! featureType attribute
561 character (len=128) :: proj4 ! proj4 attribute
562 character (len=128) :: initTime ! model_initialization_time attribute
563 character (len=128) :: validTime ! model_output_valid_time attribute
564 character (len=128) :: stDim ! station_dimension attribute
565 integer :: stOrder ! stream_order_output attribute
566 character (len=128) :: cdm ! cdm_datatype attribute
567 character (len=1024) :: esri ! esri_pe_string attribute
568 character (len=128) :: conventions ! Conventions string
569 integer :: totalValidTime ! # number of valid time (#of output files)
575 character (len=64), dimension(numGwVars) :: varNames
576 integer :: numVars = numGwVars
577 character (len=512) :: modelOutputType = "groundwater_rt"
578 character (len=64) :: modelConfigType
580 ! Output variable attributes
581 real, dimension(numGwVars) :: scaleFactor ! scale_factor values used for each
582 ! variable to converte from
585 real, dimension(numGwVars) :: addOffset ! add_offset values for each variable.
586 character (len=64), dimension(numGwVars) :: longName ! Long names for each variable.
587 character (len=64), dimension(numGwVars) :: units ! Units for each variable.
588 !character (len=64), dimension(numGwVars) :: coordNames ! Coordinate names for each variable.
589 integer(kind=4), dimension(numGwVars) :: validMinComp ! Valid min (after conversion to integer)
590 integer(kind=4), dimension(numGwVars) :: validMaxComp ! Valid max (after conversion to integer)
591 real*8, dimension(numGwVars) :: validMinDbl ! Valid minimum (before conversion to integer)
592 real*8, dimension(numGwVars) :: validMaxDbl ! Valid maximum (before converstion to integer)
593 integer(kind=4), dimension(numGwVars) :: missingComp ! Missing value attribute (after conversion to integer)
594 real, dimension(numGwVars) :: missingReal ! Missing value attribute (before conversion to integer)
595 real, dimension(numGwVars) :: fillReal ! Fill value (before conversion to integer)
596 integer(kind=4), dimension(numGwVars) :: fillComp ! Fill value (after conversion to integer)
597 integer, dimension(numGwVars) :: outFlag ! 0/1 flag to turn outputting off/on
598 integer, dimension(numGwVars) :: timeZeroFlag ! 0/1 flag to either set variable to all NDV values or output
599 ! the actual data. This was
600 ! done because time 0
601 ! output does not all contain
603 real :: modelNdv ! NDV value represented within the model code
604 ! Time variable attribues
605 character (len=64) :: timeLName ! long_name - usually valid output time
606 character (len=64) :: timeUnits ! Usually seconds since 1/1/1970
607 character (len=64) :: timeStName ! standard_name - usually time
608 integer :: timeValidMin ! the minimum time each configuration can have, time of the first output file
609 integer :: timeValidMax ! the maximum time each configuration can have, time of the last output file
610 ! Reference time attributes
611 character (len=64) :: rTimeLName ! long_name - usually model initialization time
612 character (len=64) :: rTimeStName ! standard_name - usually forecast_reference_time
613 character (len=64) :: rTimeUnits ! Usually seconds since 1/1/1970
615 character (len=64) :: gwIdLName ! long_name - usually gw bucket ID
616 character (len=256) :: gwIdComment ! Comment attribute
617 ! feature_id attributes
618 character (len=256) :: featureIdLName ! long_name - usually gw bucket ID
619 character (len=256) :: featureIdComment ! Comment attribute
620 character (len=64) :: cfRole ! cf_role attribute
621 ! latitude variable attributes
622 !character (len=64) :: latLName ! long_name
623 !character (len=64) :: latUnits ! units
624 !character (len=64) :: latStName ! Standard Name
625 ! longitude variable attributes
626 !character (len=64) :: lonLName ! long_name
627 !character (len=64) :: lonUnits ! units
628 !character (len=64) :: lonStName ! Standard Name
629 ! elevation variable attributes
630 !character (len=64) :: elevLName ! long_name
631 !character (len=64) :: elevUnits ! units
632 !character (len=64) :: elevStName ! Standard Name
634 character (len=128) :: title ! file title
635 character (len=128) :: fType ! featureType attribute
636 !character (len=128) :: proj4 ! proj4 attribute
637 character (len=128) :: initTime ! model_initialization_time attribute
638 character (len=128) :: validTime ! model_output_valid_time attribute
639 character (len=128) :: gwDim ! gw_dimension attribute
640 !character (len=128) :: cdm ! cdm_datatype attribute
641 !character (len=1024) :: esri ! esri_pe_string attribute
642 character (len=128) :: conventions ! Conventions string
643 integer :: totalValidTime ! # number of valid time (#of output files)
649 subroutine initChrtDict(chrtOutDict,diagFlag,procId)
650 use config_base, only: nlst
654 type(chrtMeta), intent(inout) :: chrtOutDict
655 integer, intent(inout) :: diagFlag
656 integer, intent(inout) :: procId
659 integer :: ftnMeta,iret
662 chrtOutDict%modelNdv = -9.E15
665 ! Pull spatial metadata information about the modeling domain from the land
666 ! spatial metadata file.
667 if(procId .eq. 0) then
668 iret = nf90_open(trim(nlst(1)%land_spatial_meta_flnm),NF90_NOWRITE,ncid=ftnMeta)
670 ! Spatial metadata file not found for land grid.
671 call postDiagMsg(diagFlag,'WARNING: Unable to open LAND spatial metadata file.')
672 chrtOutDict%proj4 = ''
673 chrtOutDict%esri = ''
675 ! First pull metadata on coordinate system.
676 iret = nf90_inq_varid(ftnMeta,'crs',projVarId)
678 call postDiagMsg(diagFlag,'WARNING: Unable to find crs in LAND spatial metadata file')
679 chrtOutDict%proj4 = ''
680 chrtOutDict%esri = ''
682 iret = nf90_get_att(ftnMeta,projVarId,'esri_pe_string',chrtOutDict%esri)
684 call postDiagMsg(diagFlag,'WARNING: Unable to find esri_pe_string in LAND spatial metadata file.')
685 chrtOutDict%esri = ''
687 iret = nf90_get_att(ftnMeta,NF90_GLOBAL,'proj4',chrtOutDict%proj4)
688 ! We are going to put a relaxed constraint on the proj4 string.
690 call postDiagMsg(diagFlag,'WARNING: proj4 string not found. Defaulting to blank string.')
691 chrtOutDict%proj4 = ''
695 iret = nf90_close(ftnMeta)
696 call nwmCheck(diagFlag,iret,'ERROR: Unable to close LAND spatial metadata file.')
700 ! NOTE !!!!! If you see PLC, this means OWP has no desire to output these,
701 ! which means meta-data standards have yet to be determined
702 ! for these variables. Fill in if it's desired to output....
703 ! First establish global attributes for the channel output files
704 chrtOutDict%title = "OUTPUT FROM " // trim(get_model_version())
705 chrtOutDict%fType = 'timeSeries'
706 chrtOutDict%initTime = '1970-01-01_00:00:00' ! This will be calculated in I/O code
707 chrtOutDict%validTime = '1970-01-01_00:00:00' ! This will be calculated in I/O code
708 chrtOutDict%stDim = 'feature_id'
709 chrtOutDict%stOrder = 1
710 chrtOutDict%cdm = 'Station'
711 chrtOutDict%conventions = 'CF-1.6'
713 ! Next establish time attribues
714 chrtOutDict%timeLName = 'valid output time'
715 chrtOutDict%timeUnits = 'minutes since 1970-01-01 00:00:00 UTC'
716 chrtOutDict%timeStName = 'time'
717 chrtOutDict%rTimeLName = 'model initialization time'
718 chrtOutDict%rTimeStName = 'forecast_reference_time'
719 chrtOutDict%rTimeUnits = 'minutes since 1970-01-01 00:00:00 UTC'
721 ! Esatablish lat/lon attributes
722 chrtOutDict%latLName = "Feature latitude"
723 chrtOutDict%latUnits = "degrees_north"
724 chrtOutDict%latStName = "latitude"
725 chrtOutDict%lonLName = "Feature longitude"
726 chrtOutDict%lonUnits = "degrees_east"
727 chrtOutDict%lonStName = "longitude"
729 ! Establish streamflw order attributes
730 chrtOutDict%orderLName = "Streamflow Order"
731 chrtOutDict%orderStName = "order"
733 ! Establish point elevation attributes
734 chrtOutDict%elevLName = "Feature Elevation"
735 chrtOutDict%elevUnits = "meters"
736 chrtOutDict%elevStName = "Elevation"
738 ! Next establish feature_id attributes. Given this is merged community, we
739 chrtOutDict%featureIdLName = 'Reach ID'
740 chrtOutDict%featureIdComment = 'NHDPlusv2 ComIDs within CONUS, arbitrary Reach IDs outside of CONUS'
741 chrtOutDict%cfRole = 'timeseries_id'
743 ! Now establish attributes for output variables.
744 !chrtOutDict%varNames(:) = (/"streamflow","nudge","q_lateral","velocity",&
745 ! "Head","qSfcLatRunoff","qBucket",&
746 ! "qBtmVertRunoff","AccSfcLatRunoff","accBucket"/)
747 chrtOutDict%varNames(:) = [character(len=64) :: "streamflow","nudge","q_lateral","velocity",&
748 "Head","qSfcLatRunoff","qBucket",&
749 "qBtmVertRunoff","AccSfcLatRunoff","accBucket","qloss"]
750 chrtOutDict%longName(:) = [character(len=64) :: "River Flow","Amount of stream flow alteration",&
751 "Runoff into channel reach","River Velocity",&
752 "River Stage","Runoff from terrain routing",&
753 "Flux from gw bucket",&
754 "Runoff from bottom of soil to bucket",&
755 "Accumulated runoff from terrain routing",&
756 "Accumulated runoff from gw bucket",&
757 "Channel Infiltration"]
758 chrtOutDict%units(:) = [character(len=64) :: "m3 s-1","m3 s-1","m3 s-1","m s-1",&
759 "meter","m3 s-1","m3 s-1",&
760 "m3","m3","m3","m3 s-1"]
761 chrtOutDict%coordNames(:) = [character(len=64) :: "latitude longitude","latitude longitude",&
762 "latitude longitude","latitude longitude",&
763 "latitude longitude","latitude longitude",&
764 "latitude longitude","latitude longitude",&
765 "latitude longitude","latitude longitude",&
766 "latitude longitude"]
767 chrtOutDict%scaleFactor(:) = [0.01,0.01,0.1,0.01,0.01,0.00001,0.00001,0.001,&
769 chrtOutDict%addOffset(:) = [0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,&
771 ! Initialize all output flags to 0. Modify (if absolutely necessary) in the
773 chrtOutDict%outFlag(:) = [0,0,0,0,0,0,0,0,0,0,0]
774 chrtOutDict%timeZeroFlag(:) = [1,1,1,1,1,1,1,1,1,1,1]
775 chrtOutDict%fillReal(:) = [-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,&
776 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0]
777 chrtOutDict%missingReal(:) = [-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,&
778 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0,&
780 chrtOutDict%validMinDbl(:) = [0.0d0, -50000.0d0, 0.0d0, 0.0d0, 0.0d0, 0.0d0, 0.0d0, 0.0d0, &
782 chrtOutDict%validMaxDbl(:) = [50000.0d0, 50000.0d0, 50000.0d0, 50000.0d0, 50000.0d0, &
783 20000.0d0, 20000.0d0, 20000.0d0, 50000.0d0, &
784 50000.0d0, 50000.0d0]
785 ! Loop through and calculate missing/fill/min/max values that will be placed
786 ! into the NetCDF attributes after scale_factor/add_offset are applied.
788 chrtOutDict%fillComp(i) = &
789 nint((chrtOutDict%fillReal(i) + chrtOutDict%addOffset(i)), 8) * &
790 nint(one_dbl / chrtOutDict%scaleFactor(i), 8)
791 chrtOutDict%missingComp(i) = &
792 nint((chrtOutDict%missingReal(i) + chrtOutDict%addOffset(i)), 8) * &
793 nint(one_dbl / chrtOutDict%scaleFactor(i), 8)
794 chrtOutDict%validMinComp(i) = &
795 nint((chrtOutDict%validMinDbl(i) + chrtOutDict%addOffset(i)) * &
796 nint(one_dbl / chrtOutDict%scaleFactor(i), 8), 8)
797 chrtOutDict%validMaxComp(i) = &
798 nint((chrtOutDict%validMaxDbl(i) + chrtOutDict%addOffset(i)) * &
799 nint(one_dbl / chrtOutDict%scaleFactor(i), 8), 8)
801 end subroutine initChrtDict
803 subroutine initLdasDict(ldasOutDict,procId,diagFlag)
804 use config_base, only: nlst
808 type(ldasMeta), intent(inout) :: ldasOutDict
809 integer, intent(inout) :: procId
810 integer, intent(inout) :: diagFlag
811 integer :: ftnMeta,projVarId,xVarId,yVarId
813 integer :: crsRealAttCnt,xRealAttCnt,yRealAttCnt
814 integer :: crsCharAttCnt,xCharAttCnt,yCharAttCnt
815 integer :: i, nCrsAtts,nxAtts,nyAtts
818 character(len=512) :: tmpAttName
824 ldasOutDict%numSoilLayers = nlst(1)%nsoil
825 ldasOutDict%act_lev = nlst(1)%act_lev
828 ldasOutDict%modelNdv = 9.9692099683868690E36
829 ldasOutDict%modelNdv2 = -1.E33
830 ldasOutDict%modelNdv3 = -1.E36
831 ldasOutDict%modelNdvInt = -2147483647
833 ! First establish global attributes.
834 ldasOutDict%title = "OUTPUT FROM " // trim(get_model_version())
835 ldasOutDict%initTime = "1970-01-01_00:00:00" ! This will be calculated in I/O code.
836 ldasOutDict%validTime = "1970-01-01_00:00:00" ! This will be calculated in I/O code.
837 ldasOutDict%conventions = "CF-1.6"
839 ! Next establish time attributes
840 ldasOutDict%timeLName = "valid output time"
841 ldasOutDict%timeUnits = "minutes since 1970-01-01 00:00:00 UTC"
842 ldasOutDict%timeStName = "time"
843 ldasOutDict%rTimeLName = "model initialization time"
844 ldasOutDict%rTimeUnits = "minutes since 1970-01-01 00:00:00 UTC"
845 ldasOutDict%rTimeStName = "forecast_reference_time"
854 ! Pull spatial metadata information about the modeling domain from the land
855 ! spatial metadata file.
856 if(procId .eq. 0) then
857 iret = nf90_open(trim(nlst(1)%land_spatial_meta_flnm),NF90_NOWRITE,ncid=ftnMeta)
859 ! Spatial metadata file not found for land grid.
860 call postDiagMsg(diagFlag,'WARNING: Unable to open LAND spatial metadata file. No crs variable or attributes will be created.')
861 ldasOutDict%nCrsCharAtts = 0
862 ldasOutDict%nCrsRealAtts = 0
863 ldasOutDict%nxCharAtts = 0
864 ldasOutDict%nxRealAtts = 0
865 ldasOutDict%nyCharAtts = 0
866 ldasOutDict%nyRealAtts = 0
867 ldasOutDict%proj4 = ''
869 iret = nf90_get_att(ftnMeta,NF90_GLOBAL,'proj4',ldasOutDict%proj4)
871 call postDiagMsg(diagFlag,'WARNING: Unable to find proj4 global attribute. Defaulting to blank string.')
872 ldasOutDict%proj4 = ''
876 ! Find the crs variable and pull out the attributes, their names, and
877 ! their values. This will be translated to output files.
878 iret = nf90_inq_varid(ftnMeta,'crs',projVarId)
879 ! For now we are going to allow the code to move forward without
880 ! finding this variable. In the future, we will probably restrict the
881 ! code to ensure things are more seamless.
883 call postDiagMsg(diagFlag,'WARNING: Unable to locate the crs variable. No crs variable or attributes will be created.')
884 ldasOutDict%nCrsCharAtts = 0
885 ldasOutDict%nCrsRealAtts = 0
887 iret = nf90_inquire_variable(ftnMeta,projVarId,nAtts=nCrsAtts)
888 call nwmCheck(diagFlag,iret,'ERROR: Unable to find crs number of attributes')
890 iret = nf90_inq_attname(ftnMeta,projVarId,i,name=tmpAttName)
891 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from crs variable.')
892 iret = nf90_inquire_attribute(ftnMeta,projVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
893 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from crs variable.')
894 select case (xtypeTmp)
916 if(charFlag .eq. 1) then
917 ldasOutDict%crsCharAttNames(crsCharAttCnt+1) = trim(tmpAttName)
918 iret = nf90_get_att(ftnMeta,projVarId,trim(tmpAttName),ldasOutDict%crsCharAttVals(crsCharAttCnt+1))
919 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull crs attributes')
920 crsCharAttCnt = crsCharAttCnt + 1
922 ldasOutDict%crsFloatAttNames(crsRealAttCnt+1) = trim(tmpAttName)
923 ldasOutDict%crsRealAttLen(crsRealAttCnt+1) = tmpLen
924 iret = nf90_get_att(ftnMeta,projVarId,trim(tmpAttName),ldasOutDict%crsRealAttVals(crsRealAttCnt+1,1:tmpLen))
925 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull crs attributes')
926 crsRealAttCnt = crsRealAttCnt + 1
931 ldasOutDict%nCrsRealAtts = crsRealAttCnt
932 ldasOutDict%nCrsCharAtts = crsCharAttCnt
936 ! Next pull the attributes from the x/y dimensions
939 iret = nf90_inq_varid(ftnMeta,'x',xVarId)
941 call postDiagMsg(diagFlag,'WARNING: Unable to locate the x variable. No x variable or attributes will be created.')
942 ldasOutDict%nxCharAtts = 0
943 ldasOutDict%nxRealAtts = 0
945 iret = nf90_inquire_variable(ftnMeta,xVarId,nAtts=nxAtts)
946 call nwmCheck(diagFlag,iret,'ERROR: Unable to find x number of attributes')
948 iret = nf90_inq_attname(ftnMeta,xVarId,i,name=tmpAttName)
949 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from x variable.')
950 iret = nf90_inquire_attribute(ftnMeta,xVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
951 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from x variable.')
952 select case (xtypeTmp)
974 if(charFlag .eq. 1) then
975 ldasOutDict%xCharAttNames(xCharAttCnt+1) = trim(tmpAttName)
976 iret = nf90_get_att(ftnMeta,xVarId,trim(tmpAttName),ldasOutDict%xCharAttVals(xCharAttCnt+1))
977 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull x attributes')
978 xCharAttCnt = xCharAttCnt + 1
980 ldasOutDict%xFloatAttNames(xRealAttCnt+1) = trim(tmpAttName)
981 ldasOutDict%xRealAttLen(xRealAttCnt+1) = tmpLen
982 iret = nf90_get_att(ftnMeta,xVarId,trim(tmpAttName),ldasOutDict%xRealAttVals(xRealAttCnt+1,1:tmpLen))
983 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull x attributes')
984 xRealAttCnt = xRealAttCnt + 1
989 ldasOutDict%nxRealAtts = xRealAttCnt
990 ldasOutDict%nxCharAtts = xCharAttCnt
996 iret = nf90_inq_varid(ftnMeta,'y',yVarId)
998 call postDiagMsg(diagFlag,'WARNING: Unable to locate the y variable. No y variable or attributes will be created.')
999 ldasOutDict%nyCharAtts = 0
1000 ldasOutDict%nyRealAtts = 0
1002 iret = nf90_inquire_variable(ftnMeta,yVarId,nAtts=nyAtts)
1003 call nwmCheck(diagFlag,iret,'ERROR: Unable to find y number of attributes')
1005 iret = nf90_inq_attname(ftnMeta,yVarId,i,name=tmpAttName)
1006 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from y variable.')
1007 iret = nf90_inquire_attribute(ftnMeta,yVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
1008 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from y variable.')
1009 select case (xtypeTmp)
1031 if(charFlag .eq. 1) then
1032 ldasOutDict%yCharAttNames(yCharAttCnt+1) = trim(tmpAttName)
1033 iret = nf90_get_att(ftnMeta,yVarId,trim(tmpAttName),ldasOutDict%yCharAttVals(yCharAttCnt+1))
1034 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull y attributes')
1035 yCharAttCnt = yCharAttCnt + 1
1037 ldasOutDict%yFloatAttNames(yRealAttCnt+1) = trim(tmpAttName)
1038 ldasOutDict%yRealAttLen(yRealAttCnt+1) = tmpLen
1039 iret = nf90_get_att(ftnMeta,yVarId,trim(tmpAttName),ldasOutDict%yRealAttVals(yRealAttCnt+1,1:tmpLen))
1040 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull y attributes')
1041 yRealAttCnt = yRealAttCnt + 1
1046 ldasOutDict%nyRealAtts = yRealAttCnt
1047 ldasOutDict%nyCharAtts = yCharAttCnt
1051 iret = nf90_close(ftnMeta)
1052 call nwmCheck(diagFlag,iret,'ERROR: Unable to close LAND spatial metadata file.')
1056 ! Now establish metadata attributes for variables. ESRI string is defined
1057 ! above and applies to all gridded variables for LDASOUT.
1058 ldasOutDict%varNames(:) = [character(len=64) :: &
1059 "IVGTYP","ISLTYP","FVEG","LAI","SAI", & !1-5
1060 "SWFORC","COSZ","LWFORC","RAINRATE","EMISS", & !6-10
1061 "FSA","FIRA","GRDFLX","HFX","LH", & !11-15
1062 "ECAN","EDIR","ALBEDO","ETRAN","UGDRNOFF", & !16-20
1063 "SFCRNOFF","CANLIQ","CANICE","ZWT","WA", & !21-25
1064 "WT","ACCPRCP","ACCECAN","ACCEDIR","ACCETRAN", & !26-30
1065 "SAV","TR","EVC","IRC","SHC", & !31-35
1066 "IRG","SHG","EVG","GHV","SAG", & !36-40
1067 "IRB","SHB","EVB","GHB","TRAD", & !41-45
1068 "TG","TV","TAH","TGV","TGB", & !46-50
1069 "T2MV","T2MB","Q2MV","Q2MB","EAH", & !51-55
1070 "FWET","ZSNSO_SN","SNICE","SNLIQ","SOIL_T", & !56-60
1071 "SOIL_W","SNOW_T","SOIL_M","SNOWH","SNEQV", & !60-65
1072 "QSNOW","ISNOW","FSNO","ACSNOW","ACSNOM", & !66-70
1073 "CM","CH","CHV","CHB","CHLEAF", & !71-75
1074 "CHUC","CHV2","CHB2","LFMASS","RTMASS", & !76-80
1075 "STMASS","WOOD","STBLCP","FASTCP","NEE", & !81-85
1076 "GPP","NPP","PSN","APAR","ACCET", & !86-90
1077 "CANWAT","SOILICE","SOILSAT_TOP","SOILSAT","SNOWT_AVG", & !91-95
1078 "ALBSND","ALBSNI","QRAIN",& !96-98
1079 "glacier", "glacier_thickness" ,"PSNOWALB",& !99-101
1080 "PSNOWTHRUFAL" ,"PSNOWHEIGHT","PSNOWTOTSWE" ,& !102-104
1081 "PSNOWGRAN1","PSNOWGRAN2","PSNOWAGE",& !105-107
1082 "PSNOWTEMP","PSNOWDZ","PSNOWHIST",& !108-110
1083 "PSNOWLIQ","PSNOWHEAT","PSNOWRHO",& !111-113
1084 "PSNOWSWE", "FLOW_ICE", "FLOW_SNOW"] !114-116
1086 ldasOutDict%longName(:) = [character(len=128) :: &
1087 "Dominant vegetation category",& !1
1088 "Dominant soil category",& !2
1089 "Green Vegetation Fraction",& !3
1090 "Leaf area index",& !4
1091 "Stem area index",& !5
1092 "Shortwave forcing",& !6
1093 "Cosine of zenith angle",& !7
1094 "Longwave forcing",& !8
1095 "Precipitation rate",& !9
1096 "Grid emissivity",& !10
1097 "Total absorbed SW radiation",& !11
1098 "Total net LW radiation to atmosphere",& !12
1099 "Heat flux into the soil",& !13
1100 "Total sensible heat to the atmosphere",& !14
1101 "Total latent heat to the atmosphere",& !15
1102 "Canopy water evaporation rate",& !16
1103 "Direct from soil evaporation rate",& !17
1104 "Surface albedo",& !18
1105 "Transpiration rate",& !19
1106 "Accumulated underground runoff",& !20
1107 "Accumulated surface runoff",& !21
1108 "Canopy liquid water content",& !22
1109 "Canopy ice water content",& !23
1110 "Depth to water table",& !24
1111 "Water in aquifer",& !25
1112 "Water in aquifer and saturated soil",& !26
1113 "Accumulated precip",& !27
1114 "Accumulated canopy evap",& !28
1115 "Accumulated direct soil evap",& !29
1116 "Accumulated transpiration",& !30
1117 "Solar radiative heat flux absorbed by vegetation",& !31
1118 "Transpiration heat",& !32
1119 "Canopy evap heat",& !33
1120 "Canopy net LW rad",& !34
1121 "Canopy sensible heat",& !35
1122 "Ground net LW rad",& !36
1123 "Ground sensible heat",& !37
1124 "Ground evap heat",& !38
1125 "Ground heat flux + to soil vegetated",& !39
1126 "Solar radiative heat flux absorbed by ground",& !40
1127 "Net LW rad to atm bare",& !41
1128 "Sensible heat atm bare",& !42
1129 "Evaporation heat to atm bare",& !43
1130 "Ground heat flux + to soil bare",& !44
1131 "Surface radiative temperature",& !45
1132 "Ground temperature",& !46
1133 "Vegetation temperature",& !47
1134 "Canopy air temperature",& !48
1135 "Ground surface Temp vegetated",& !49
1136 "Ground surface Temp bare",& !50
1137 "2m Air Temp vegetated",& !51
1138 "2m Air Temp bare",& !52
1139 "2m mixing ratio vegetated",& !53
1140 "2m mixing ratio bare",& !54
1141 "Canopy air vapor pressure",& !55
1142 "Wetted or snowed fraction of canopy",& !56
1143 "Snow layer depths from snow surface",& !57
1144 "Snow layer ice",& !58
1145 "Snow layer liquid water",& !59
1146 "soil temperature",& !60
1147 "liquid volumetric soil moisture",& !61
1148 "snow temperature",& !62
1149 "volumetric soil moisture, the dimensionless ratio of water volume (m3) to soil volume (m3)",& !63
1151 "Snow water equivalent",& !65
1152 "Snowfall rate on the ground",& !66
1153 "Number of snow layers",& !67
1154 "Snow-cover fraction on the ground",& !68
1155 "accumulated snow fall",& !69
1156 "accumulated melting water out of snow bottom",& !70
1157 "Momentum drag coefficient",& !71
1158 "Sensible heat exchange coefficient",& !72
1159 "Exchange coefficient vegetated",& !73
1160 "Exchange coefficient bare",& !74
1161 "Exchange coefficient leaf",& !75
1162 "Exchange coefficient bare",& !76
1163 "Exchange coefficient 2-meter vegetated",&!77
1164 "Exchange coefficient 2-meter bare",& !78
1166 "Mass of fine roots",& !80
1168 "Mass of wood and woody roots",& !82
1169 "Stable carbon in deep soil",& !83
1170 "Short-lived carbon in shallow soil",& !84
1171 "Net ecosystem exchange",& !85
1172 "Net instantaneous assimilation",& !86
1173 "Net primary productivity",& !87
1174 "Total photosynthesis",& !88
1175 "Photosynthesis active energy by canopy",&!89
1176 "Accumulated total ET",& !90
1177 "Total canopy water (liquid + ice)",& !91
1178 "fraction of soil moisture that is ice",& !92
1179 "fraction of soil saturation, top 2 layers",& !93
1180 "fraction of soil saturation, column integrated",& !94
1181 "average snow temperature (by layer mass)",& !95
1182 "snowpack albedo, direct",& !96
1183 "snowpack albedo, diffuse",& !97
1184 "Rainfall rate on the ground",& !98
1185 "Glacier grid point",& !99
1186 "Glacier height",& !100
1187 "Snow albedo",& !101
1188 "Runoff from glacier",& !102
1189 "Total snow height",& !103
1190 "Total snow swe", & !104
1191 "Snow gran 1, optiacal diameter",& !105
1192 "Snow gran 2, sphericity",& !106
1194 "Snow temperature",& !108
1195 "Snow thickness",& !109
1196 "Snow history",& !110
1197 "Liquid in snow",& !111
1199 "Snow density",& !113
1200 "Snow water equivalent", & !114
1201 "Accumulated glacier melt from ice", & !115
1202 "Accumulated glacier melt from snow"] !116
1204 ldasOutDict%units(:) = [character(len=64) :: &
1205 "category","category","-","-","-", & !1-5
1206 "W m-2","-","W m-2","kg m-2 s-1","-", & !6-10
1207 "W m-2","W m-2","W m-2","W m-2","W m-2", & !11-15
1208 "kg m-2 s-1","kg m-2 s-1","-","kg m-2 s-1","mm", & !16-20
1209 "mm","mm","mm","m","kg m-2", & !21-25
1210 "kg m-2","mm","mm","mm","mm", & !26-30
1211 "W m-2","W m-2","W m-2","W m-2","W m-2", & !31-35
1212 "W m-2","W m-2","W m-2","W m-2","W m-2", & !36-40
1213 "W m-2","W m-2","W m-2","W m-2","K", & !41-45
1214 "K","K","K","K","K", & !46-50
1215 "K","K","kg/kg","kg/kg","Pa", & !51-55
1216 "fraction","m","mm","mm","K", & !56-60
1217 "m3 m-3","K","m3 m-3","m","kg m-2", & !61-65
1218 "mm s-1","count","1","mm","mm", & !66-70
1219 "-","-","m s-1","m s-1","m s-1", & !71-75
1220 "m s-1","m s-1","m s-1","g m-2","g m-2", & !76-80
1221 "g m-2","g m-2","g m-2","g m-2","g m-2s-1 CO2", & !81-85
1222 "g m-2s-1 C","g m-2s-1 C","umol CO2 m-2 s-1","W m-2","mm", & !86-90
1223 "mm","1","1","1","K", & !91-95
1224 "-","-","mm s-1",& !96-98
1225 "-","m","-", & !99-101
1226 "kg/(m2 s)","m","kg m-2",& !102-104
1227 "m","-","days since snowfall",& !105-107
1228 "K","m","-",& !108-110
1229 "kg/m3","J/m2","m",& !111-113
1230 "kg m-2","kg/m2","kg/m2"] !114-116
1232 ldasOutDict%scaleFactor(:) = [1.0, 1.0, 0.01, 0.1, 0.1, & !1-5
1233 0.1, 0.01, 0.1, 0.00001, 0.01, & !6-10
1234 0.1, 0.1, 0.1, 0.1, 0.1, & !11-15
1235 0.00001, 0.00001, 0.01, 0.00001, 0.01, & !16-20
1236 0.001, 0.01, 0.01, 0.00001, 0.01, & !21-25
1237 0.01, 0.01, 0.01, 0.01, 0.01, & !26-30
1238 0.1, 0.1, 0.1, 0.1, 0.1, & !31-35
1239 0.1, 0.1, 0.1, 0.1, 0.1, & !36-40
1240 0.1, 0.1, 0.1, 0.1, 0.1, & !41-45
1241 0.1, 0.1, 0.1, 0.1, 0.1, & !46-50
1242 0.1, 0.1, 0.0001, 0.0001, 0.1, & !51-55
1243 0.01, 0.00001, 0.01, 0.01, 0.1, & !56-60
1244 0.01, 0.1, 0.01, 0.0001, 0.1, & !61-65
1245 0.00001, 1.0, 0.001, 0.01, 0.01, & !66-70
1246 0.00001, 0.00001, 0.00001, 0.00001, 0.00001, &!71-75
1247 0.00001, 0.00001, 0.00001, 0.01, 0.01, & !76-80
1248 0.01, 0.01, 0.01, 0.01, 0.01, & !81-85
1249 0.01, 0.01, 0.01, 0.01, 0.01, & !86-90
1250 0.01, 0.01, 0.001, 0.001, 0.1, & !91-95
1251 0.01, 0.01, 0.00001, & !96-98
1252 1.0, 0.1, 0.01, & !99-101
1253 0.0001, 0.0001, 0.1, & !102-104
1254 0.01, 0.01, 0.01, & !105-107
1255 0.1, 0.0001, 0.01, & !108-110
1256 0.001, 1000.0, 0.1, & !111-113
1257 0.1, 0.001, 0.001 ] !114-116
1259 ldasOutDict%addOffset(:) = [0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0, & !1-10
1260 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0, & !11-20
1261 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0, & !21-30
1262 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0, & !31-40
1263 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0, & !41-50
1264 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0, & !51-60
1265 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0, & !61-70
1266 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0, & !71-80
1267 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0, & !81-90
1268 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0, & !91-98
1269 0.0,0.0,0.0, & !99-101
1270 0.0,0.0,0.0, & !102-104
1271 0.0,0.0,0.0, & !105-107
1272 0.0,0.0,0.0, & !108-110
1273 0.0,0.0,0.0, & !111-113
1274 0.0,0.0,0.0] !114-116
1276 ! Note that output flags will be set in the the output routine, and will vary
1277 ! by the IOC flag specified in hydro.namelist.
1278 ldasOutDict%outFlag(:) = [0,0,0,0,0,0,0,0,0,0, & !1-10
1279 1,0,0,0,0,0,0,0,0,0, & !11-20
1280 0,0,0,0,0,0,0,0,0,0, & !21-30
1281 0,0,0,0,0,0,0,0,0,0, & !31-40
1282 0,0,0,0,0,0,0,0,0,0, & !41-50
1283 0,0,0,0,0,0,0,0,0,0, & !51-60
1284 0,0,0,0,0,0,0,0,0,0, & !61-70
1285 0,0,0,0,0,0,0,0,0,0, & !71-80
1286 0,0,0,0,0,0,0,0,0,0, & !81-90
1287 0,0,0,0,0,0,0,0, & !91-98
1295 ldasOutDict%timeZeroFlag(:) = [1,1,1,1,1,1,1,1,1,1, & !1-10
1296 1,1,1,1,1,1,1,1,1,1, & !11-20
1297 1,1,1,1,1,1,1,1,1,1, & !21-30
1298 1,1,1,1,1,1,1,1,1,1, & !31-40
1299 1,1,1,1,1,1,1,1,1,1, & !41-50
1300 1,1,1,1,1,1,1,1,1,1, & !51-60
1301 1,1,1,1,1,1,1,1,1,1, & !61-70
1302 1,1,1,1,1,1,1,1,1,1, & !71-80
1303 1,1,1,1,1,1,1,1,1,1, & !81-90
1304 1,1,1,1,1,1,1,1, & !91-98
1312 ldasOutDict%numLev(:) = [1,1,1,1,1,1,1,1,1,1, & !1-10
1313 1,1,1,1,1,1,1,1,1,1, & !11-20
1314 1,1,1,1,1,1,1,1,1,1, & !21-30
1315 1,1,1,1,1,1,1,1,1,1, & !31-40
1316 1,1,1,1,1,1,1,1,1,1, & !41-50
1317 1,1,1,1,1,1,3,3,3,4, & !51-60
1318 4,3,4,1,1,1,1,1,1,1, & !61-70
1319 1,1,1,1,1,1,1,1,1,1, & !71-80
1320 1,1,1,1,1,1,1,1,1,1, & !81-90
1321 1,1,1,1,1,2,2,1, & !91-98
1324 40,40,40, & !105-107
1325 40,40,40, & !108-110
1326 40,40,40, & !111-113
1328 ldasOutDict%numLev(105:114) = ldasOutDict%act_lev ! Set crocus levels to number from namelist
1330 ldasOutDict%missingReal(:) = [-9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !1-5
1331 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !6-10
1332 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !11-15
1333 -9999.0, -999.0,-9999.0,-9999.0,-9999.0, & !16-20
1334 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !21-25
1335 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !26-30
1336 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !31-35
1337 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !36-40
1338 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !41-45
1339 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !46-50
1340 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !51-55
1341 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !56-60
1342 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !61-65
1343 -999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !66-70
1344 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !71-75
1345 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !76-80
1346 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !81-85
1347 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !86-90
1348 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !91-95
1349 -9999.0,-9999.0,-999.0, & !96-98
1350 -9999.0,-9999.0,-9999.0, & !99-101
1351 -9999.0,-9999.0,-9999.0, & !102-104
1352 -9999.0,-9999.0,-9999.0, & !105-107
1353 -9999.0,-9999.0,-9999.0, & !108-110
1354 -9999.0, 9999000.0,-9999.0, & !111-113
1355 -9999.0,-9999.0,-9999.0 ] !114-116
1357 ldasOutDict%fillReal(:) = [-9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !1-5
1358 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !6-10
1359 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !11-15
1360 -9999.0, -999.0,-9999.0,-9999.0,-9999.0, & !16-20
1361 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !21-25
1362 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !26-30
1363 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !31-35
1364 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !36-40
1365 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !41-45
1366 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !46-50
1367 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !51-55
1368 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !56-60
1369 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !61-65
1370 -999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !66-70
1371 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !71-75
1372 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !76-80
1373 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !81-85
1374 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !86-90
1375 -9999.0,-9999.0,-9999.0,-9999.0,-9999.0, & !91-95
1376 -9999.0,-9999.0,-999.0, & !96-98
1377 -9999.0,-9999.0,-9999.0, & !99-101
1378 -9999.0,-9999.0,-9999.0, & !102-104
1379 -9999.0,-9999.0,-9999.0, & !105-107
1380 -9999.0,-9999.0,-9999.0, & !108-110
1381 -9999.0, 9999000.0,-9999.0, & !111-113
1382 -9999.0,-9999.0,-9999.0 ] !114-116
1384 ldasOutDict%validMinDbl(:) = [0.0d0, 0.0d0, 0.0d0, 0.0d0, 0.0d0, & !1-5
1385 -1000.0d0, -1.0d0, -1500.0d0, 0.0d0, 0.0d0, & !6-10
1386 -1500.0d0, -1500.0d0, -1500.0d0, -1500.0d0, -1500.0d0, & !11-15
1387 -100.0d0, -100.0d0, 0.0d0, -100.0d0, -100.0d0, & !16-20
1388 0.0d0, -5.0d0, -5.0d0, 0.0d0, 0.0d0, & !21-25
1389 0.0d0, 0.0d0, -100.0d0, -100.0d0, -100.0d0, & !26-30
1390 -1500.0d0, -1500.0d0, -1500.0d0, -1500.0d0, -1500.0d0, & !31-35
1391 -1500.0d0, -1500.0d0, -1500.0d0, -1500.0d0, -1500.0d0, & !36-40
1392 -1500.0d0, -1500.0d0, -1500.0d0, -1500.0d0, 0.0d0, & !41-45
1393 0.0d0, 0.0d0, 0.0d0, 0.0d0, 0.0d0, & !46-50
1394 0.0d0, 0.0d0, 0.0d0, 0.0d0, -1000.0d0, & !51-55
1395 0.0d0, -100.0d0, 0.0d0, 0.0d0, 0.0d0, & !56-60
1396 0.0d0, 0.0d0, 0.0d0, 0.0d0, 0.0d0, & !61-65
1397 0.0d0, -10.d0, 0.0d0, 0.0d0, 0.0d0, & !66-70
1398 -5.0d0, -5.0d0, -5.0d0, -5.0d0, -5.0d0, & !71-75
1399 -5.0d0, -5.0d0, -5.0d0, 0.0d0, 0.0d0, & !76-80
1400 0.0d0, 0.0d0, 0.0d0, 0.0d0, 0.0d0, & !81-85
1401 0.0d0, 0.0d0, 0.0d0, 0.0d0, -1000.0d0, & !86-90
1402 -5.0d0, 0.0d0, 0.0d0, 0.0d0, 0.0d0, & !91-95
1403 0.0d0, 0.0d0, 0.0d0,& !96-98
1404 !! NBNB Check these values
1405 0.0d0, 0.0d0, 0.0d0, & !99-101
1406 0.0d0, 0.0d0, 0.0d0, & !102-104
1407 -1.0d2, 0.0d0, 0.0d0, & !105-107
1408 0.0d0, 0.0d0, 0.0d0, & !108-110
1409 0.0d0, -2.0d12, 0.0d0, & !111-113
1410 0.0d0, 0.0d0, 0.0d0 ] !114-116
1412 ldasOutDict%validMaxDbl(:) = [100.0d0, 100.0d0, 1.0d0, 20.0d0, 20.0d0, & !1-5
1413 3000.0d0, 1.0d0, 1500.0d0, 100.0d0, 1.0d0, & !6-10
1414 1500.0d0, 1500.0d0, 1500.0d0, 1500.0d0, 1500.0d0, & !11-15
1415 100.0d0, 100.0d0, 1.0d0, 100.0d0, 100000.0d0, & !16-20
1416 100000.0d0, 30000.0d0, 30000.0d0, 10.0d0, 10000.0d0, & !21-25
1417 10000.0d0, 1.0D+6, 1.0D+6, 1.0D+6, 1.0D+6, & !26-30
1418 1500.0d0, 1500.0d0, 1500.0d0, 1500.0d0, 1500.0d0, & !31-35
1419 1500.0d0, 1500.0d0, 1500.0d0, 1500.0d0, 1500.0d0, & !36-40
1420 1500.0d0, 1500.0d0, 1500.0d0, 1500.0d0, 400.0d0, & !41-45
1421 400.0d0, 400.0d0, 400.0d0, 400.0d0, 400.0d0, & !46-50
1422 400.0d0, 400.0d0, 1.0d0, 1.0d0, 100000.0d0, & !51-55
1423 1.0d0, 100.0d0, 100000.0d0, 100000.0d0, 400.0d0, & !56-60
1424 1.0d0, 400.0d0, 1.0d0, 10000.0d0, 100000000.0d0, & !61-65
1425 100.0d0, 10.0d0, 1.0d0, 1.0D+6, 100000.0d0, & !66-70
1426 5.0d0, 5.0d0, 5.0d0, 5.0d0, 5.0d0, & !71-76
1427 5.0d0, 5.0d0, 5.0d0, 1000.0d0, 1000.0d0, & !76-80
1428 1000.0d0, 1000.0d0, 5000.0d0, 5000.0d0, 1000.0d0, & !81-85
1429 1000.0d0, 1000.0d0, 1000.0d0, 1000.0d0, 1.0D+6, & !86-90
1430 30000.0d0, 1.0d0, 1.0d0, 1.0d0, 400.0d0, & !91-95
1431 1.0d0, 1.0d0, 100.0d0,& !96-98
1432 ! NBNB Check these values
1433 1.0d0, 1.0d4, 1.0d0, & !99-101
1434 1.0d5, 1.0d4, 1.0d8, & !102-104
1435 1.0D2, 1.0D2, 1.0d5, & !105-107
1436 300.0d0, 1.0d4, 1.0d2, & !108-110
1437 1.0d5, 0.0d0, 1.0d3, & !111-113
1438 1.0d6, 1.d5, 1.d5] !114-116
1440 ! Loop through and calculate missing/fill/min/max values that will be placed
1441 ! into the NetCDF attributes after scale_factor/add_offset are applied.
1443 if (ldasOutDict%scaleFactor(i) .le. 1. .and. ldasOutDict%scaleFactor(i) .gt. 0.) then
1444 ldasOutDict%fillComp(i) = &
1445 nint((ldasOutDict%fillReal(i) + ldasOutDict%addOffset(i)), 4) * &
1446 nint(one_dbl / ldasOutDict%scaleFactor(i), 4)
1447 ldasOutDict%missingComp(i) = &
1448 nint((ldasOutDict%missingReal(i) + ldasOutDict%addOffset(i)), 4) * &
1449 nint(one_dbl / ldasOutDict%scaleFactor(i), 4)
1450 ldasOutDict%validMinComp(i) = &
1451 nint((ldasOutDict%validMinDbl(i) + ldasOutDict%addOffset(i)) * &
1452 nint(one_dbl / ldasOutDict%scaleFactor(i), 4), 4)
1453 ldasOutDict%validMaxComp(i) = &
1454 nint((ldasOutDict%validMaxDbl(i) + ldasOutDict%addOffset(i)) * &
1455 nint(one_dbl / ldasOutDict%scaleFactor(i), 4), 4)
1456 else if (ldasOutDict%scaleFactor(i) .gt. 1.) then
1457 ldasOutDict%fillComp(i) = &
1458 nint( (ldasOutDict%fillReal(i) + ldasOutDict%addOffset(i)) * &
1459 (one_dbl / ldasOutDict%scaleFactor(i)), 4 )
1460 ldasOutDict%missingComp(i) = &
1461 nint( (ldasOutDict%missingReal(i) + ldasOutDict%addOffset(i)) * &
1462 (one_dbl / ldasOutDict%scaleFactor(i)), 4 )
1463 ldasOutDict%validMinComp(i) = &
1464 nint( (ldasOutDict%validMinDbl(i) + ldasOutDict%addOffset(i)) * &
1465 (one_dbl / ldasOutDict%scaleFactor(i)), 4 )
1466 ldasOutDict%validMaxComp(i) = &
1467 nint( (ldasOutDict%validMaxDbl(i) + ldasOutDict%addOffset(i)) * &
1468 (one_dbl / ldasOutDict%scaleFactor(i)), 4 )
1472 end subroutine initLdasDict
1474 subroutine initRtDomainDict(rtDomainDict,procId,diagFlag)
1475 use config_base, only: nlst
1479 type(rtDomainMeta), intent(inout) :: rtDomainDict
1480 integer, intent(inout) :: procId,diagFlag
1481 integer :: ftnMeta,projVarId,xVarId,yVarId,ftnGeo
1482 integer :: xDimId,numColLand,numColHydro
1483 real :: resLand,resHydro,aggFactor
1485 integer :: crsRealAttCnt,xRealAttCnt,yRealAttCnt
1486 integer :: crsCharAttCnt,xCharAttCnt,yCharAttCnt
1487 integer :: i, nCrsAtts,nxAtts,nyAtts
1489 integer :: floatFlag
1490 character(len=512) :: tmpAttName
1495 rtDomainDict%numSoilLayers = nlst(1)%nsoil
1497 rtDomainDict%modelNdv = -9.E15
1499 ! First establish global attributes.
1500 rtDomainDict%title = "OUTPUT FROM " // trim(get_model_version())
1501 rtDomainDict%initTime = "1970-01-01_00:00:00" ! This will be calculated in I/O code.
1502 rtDomainDict%validTime = "1970-01-01_00:00:00" ! This will be calculated in I/O code.
1503 rtDomainDict%decimation = 1
1504 rtDomainDict%conventions = "CF-1.6"
1506 ! Next establish time attributes
1507 rtDomainDict%timeLName = "valid output time"
1508 rtDomainDict%timeUnits = "minutes since 1970-01-01 00:00:00 UTC"
1509 rtDomainDict%timeStName = "time"
1510 rtDomainDict%rTimeLName = "model initialization time"
1511 rtDomainDict%rTimeUnits = "minutes since 1970-01-01 00:00:00 UTC"
1512 rtDomainDict%rTimeStName = "forecast_reference_time"
1521 ! Pull spatial metadata information about the modeling domain from the
1523 if(procId .eq. 0) then
1524 iret = nf90_open(trim(nlst(1)%geo_finegrid_flnm),NF90_NOWRITE,ncid=ftnMeta)
1525 if(iret .ne. 0) then
1526 ! Spatial metadata file not found for hydro grid.
1527 call postDiagMsg(diagFlag,'WARNING: Unable to open Fulldom metadata file. No crs variable or attributes will be created.')
1528 rtDomainDict%nCrsCharAtts = 0
1529 rtDomainDict%nCrsRealAtts = 0
1530 rtDomainDict%nxCharAtts = 0
1531 rtDomainDict%nxRealAtts = 0
1532 rtDomainDict%nyCharAtts = 0
1533 rtDomainDict%nyRealAtts = 0
1534 rtDomainDict%proj4 = ''
1536 iret = nf90_get_att(ftnMeta,NF90_GLOBAL,'proj4',rtDomainDict%proj4)
1537 if(iret .ne. 0) then
1538 call postDiagMsg(diagFlag,'WARNING: Unable to find proj4 global attribute. Defaulting to blank string.')
1539 rtDomainDict%proj4 = ''
1543 ! Find the crs variable and pull out the attributes, their names, and
1544 ! their values. This will be translated to output files.
1545 iret = nf90_inq_varid(ftnMeta,'crs',projVarId)
1546 ! For now we are going to allow the code to move forward without
1547 ! finding this variable. In the future, we will probably restrict the
1548 ! code to ensure things are more seamless.
1550 if(iret .ne. 0) then
1551 call postDiagMsg(diagFlag,'WARNING: Unable to locate the crs variable. No crs variable or attributes will be created.')
1552 rtDomainDict%nCrsCharAtts = 0
1553 rtDomainDict%nCrsRealAtts = 0
1555 iret = nf90_inquire_variable(ftnMeta,projVarId,nAtts=nCrsAtts)
1556 call nwmCheck(diagFlag,iret,'ERROR: Unable to find crs number of attributes')
1558 iret = nf90_inq_attname(ftnMeta,projVarId,i,name=tmpAttName)
1559 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from crs variable.')
1560 iret = nf90_inquire_attribute(ftnMeta,projVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
1561 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from crs variable.')
1562 select case (xtypeTmp)
1584 if(charFlag .eq. 1) then
1585 rtDomainDict%crsCharAttNames(crsCharAttCnt+1) = trim(tmpAttName)
1586 iret = nf90_get_att(ftnMeta,projVarId,trim(tmpAttName),rtDomainDict%crsCharAttVals(crsCharAttCnt+1))
1587 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull crs attributes')
1588 crsCharAttCnt = crsCharAttCnt + 1
1590 rtDomainDict%crsFloatAttNames(crsRealAttCnt+1) = trim(tmpAttName)
1591 rtDomainDict%crsRealAttLen(crsRealAttCnt+1) = tmpLen
1592 iret = nf90_get_att(ftnMeta,projVarId,trim(tmpAttName),rtDomainDict%crsRealAttVals(crsRealAttCnt+1,1:tmpLen))
1593 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull crs attributes')
1594 crsRealAttCnt = crsRealAttCnt + 1
1599 rtDomainDict%nCrsRealAtts = crsRealAttCnt
1600 rtDomainDict%nCrsCharAtts = crsCharAttCnt
1604 ! Next pull the attributes from the x/y dimensions
1607 iret = nf90_inq_varid(ftnMeta,'x',xVarId)
1608 if(iret .ne. 0) then
1609 call postDiagMsg(diagFlag,'WARNING: Unable to locate the x variable. No x variable or attributes will be created.')
1610 rtDomainDict%nxCharAtts = 0
1611 rtDomainDict%nxRealAtts = 0
1613 iret = nf90_inquire_variable(ftnMeta,xVarId,nAtts=nxAtts)
1614 call nwmCheck(diagFlag,iret,'ERROR: Unable to find x number of attributes')
1616 iret = nf90_inq_attname(ftnMeta,xVarId,i,name=tmpAttName)
1617 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from x variable.')
1618 iret = nf90_inquire_attribute(ftnMeta,xVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
1619 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from x variable.')
1620 select case (xtypeTmp)
1642 if(charFlag .eq. 1) then
1643 rtDomainDict%xCharAttNames(xCharAttCnt+1) = trim(tmpAttName)
1644 iret = nf90_get_att(ftnMeta,xVarId,trim(tmpAttName),rtDomainDict%xCharAttVals(xCharAttCnt+1))
1645 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull x attributes')
1646 xCharAttCnt = xCharAttCnt + 1
1648 rtDomainDict%xFloatAttNames(xRealAttCnt+1) = trim(tmpAttName)
1649 rtDomainDict%xRealAttLen(xRealAttCnt+1) = tmpLen
1650 iret = nf90_get_att(ftnMeta,xVarId,trim(tmpAttName),rtDomainDict%xRealAttVals(xRealAttCnt+1,1:tmpLen))
1651 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull x attributes')
1652 xRealAttCnt = xRealAttCnt + 1
1657 rtDomainDict%nxRealAtts = xRealAttCnt
1658 rtDomainDict%nxCharAtts = xCharAttCnt
1664 iret = nf90_inq_varid(ftnMeta,'y',yVarId)
1665 if(iret .ne. 0) then
1666 call postDiagMsg(diagFlag,'WARNING: Unable to locate the y variable. No y variable or attributes will be created.')
1667 rtDomainDict%nyCharAtts = 0
1668 rtDomainDict%nyRealAtts = 0
1670 iret = nf90_inquire_variable(ftnMeta,yVarId,nAtts=nyAtts)
1671 call nwmCheck(diagFlag,iret,'ERROR: Unable to find y number of attributes')
1673 iret = nf90_inq_attname(ftnMeta,yVarId,i,name=tmpAttName)
1674 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from y variable.')
1675 iret = nf90_inquire_attribute(ftnMeta,yVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
1676 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from y variable.')
1677 select case (xtypeTmp)
1699 if(charFlag .eq. 1) then
1700 rtDomainDict%yCharAttNames(yCharAttCnt+1) = trim(tmpAttName)
1701 iret = nf90_get_att(ftnMeta,yVarId,trim(tmpAttName),rtDomainDict%yCharAttVals(yCharAttCnt+1))
1702 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull y attributes')
1703 yCharAttCnt = yCharAttCnt + 1
1705 rtDomainDict%yFloatAttNames(yRealAttCnt+1) = trim(tmpAttName)
1706 rtDomainDict%yRealAttLen(yRealAttCnt+1) = tmpLen
1707 iret = nf90_get_att(ftnMeta,yVarId,trim(tmpAttName),rtDomainDict%yRealAttVals(yRealAttCnt+1,1:tmpLen))
1708 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull y attributes')
1709 yRealAttCnt = yRealAttCnt + 1
1714 rtDomainDict%nyRealAtts = yRealAttCnt
1715 rtDomainDict%nyCharAtts = yCharAttCnt
1718 iret = nf90_close(ftnMeta)
1719 call nwmCheck(diagFlag,iret,'ERROR: Unable to close Fulldom file.')
1722 ! Next get the number of columns on the land grid. This will be used to
1723 ! calculate the resolution of the routing grid in meters and aggfactor.
1724 iret = nf90_open(trim(nlst(1)%land_spatial_meta_flnm),NF90_NOWRITE,ncid=ftnGeo)
1725 if(iret .ne. 0) then
1726 call postDiagMsg(diagFlag,'WARNING: Unable to open LAND Spatial Metadata file. Defaulting to missing values.')
1730 iret = nf90_inq_dimid(ftnGeo,'x',xDimId)
1731 call nwmCheck(diagFlag,iret,'ERROR: Unable to find x dimension in LAND spatial Metadata file')
1732 iret = nf90_inquire_dimension(ftnGeo,xDimId,len=numColLand)
1733 call nwmCheck(diagFlag,iret,'ERROR: Unable to retrieve number of columns in the LAND spatial metadata file')
1734 iret = nf90_inq_varid(ftnGeo,'x',xVarId)
1735 call nwmCheck(diagFlag,iret,'ERROR: Unable to find x variable in LAND spatial metadata file')
1736 iret = nf90_get_att(ftnGeo,xVarId,'resolution',resLand)
1737 call nwmCheck(diagFlag,iret,'ERROR: Unable to get x resolution in LAND spatial metadata file')
1738 iret = nf90_close(ftnGeo)
1739 call nwmCheck(diagFlag,iret,'ERROR: Unable to close the LAND spatial metadata file')
1742 ! Next get the number of columns on the high-resolution routing grid.
1743 ! This will be used to calculate the resolution of the routing grid in
1745 iret = nf90_open(trim(nlst(1)%geo_finegrid_flnm),NF90_NOWRITE,ncid=ftnGeo)
1746 if(iret .ne. 0) then
1747 call nwmCheck(diagFlag,iret,'ERROR: Unable to open Fulldom file')
1749 iret = nf90_inq_dimid(ftnGeo,'x',xDimId)
1750 call nwmCheck(diagFlag,iret,'ERROR: Unable to find x dimension in Fulldom file')
1751 iret = nf90_inquire_dimension(ftnGeo,xDimId,len=numColHydro)
1752 call nwmCheck(diagFlag,iret,'ERROR: Unable to retrieve number of columns in the Fulldom file')
1753 iret = nf90_close(ftnGeo)
1754 call nwmCheck(diagFlag,iret,'ERROR: Unable to close the Fulldom file')
1757 ! Calculate the aggregation factor and resolution of the hydro routing
1759 if(numColLand .ne. -9999) then
1760 aggFactor = float(numColHydro/numColLand)
1761 resHydro = resLand/aggFactor
1766 rtDomainDict%xRes = resHydro
1767 rtDomainDict%yRes = resHydro
1771 rtDomainDict%varNames(:) = [character(len=64) :: "zwattablrt","sfcheadsubrt","QSTRMVOLRT",&
1773 rtDomainDict%longName(:) = [character(len=64) :: "depth to saturation, rounded to highest saturated layer",&
1774 "surface head","channel inflow",&
1775 "accumulated value of the boundary flux, + into domain - out of domain",&
1776 "volumetric soil moisture"]
1777 rtDomainDict%units(:) = [character(len=64) :: "m","mm","mm","mm","m3 m-3"]
1778 rtDomainDict%scaleFactor(:) = [0.1,1.0,1.0,1.0,0.01]
1779 rtDomainDict%addOffset(:) = [0.0,0.0,0.0,0.0,0.0]
1780 rtDomainDict%outFlag(:) = [0,0,0,0,0]
1781 rtDomainDict%timeZeroFlag(:) = [1,1,1,1,1]
1782 rtDomainDict%missingReal(:) = [-9999.0,-9999.0,-9999.0,-9999.0,-9999.0]
1783 rtDomainDict%fillReal(:) = [-9999.0,-9999.0,-9999.0,-9999.0,-9999.0]
1784 rtDomainDict%validMinDbl(:) = [0.0d0, 0.0d0, 0.0d0, -1000000.0d0, 0.0d0]
1785 rtDomainDict%validMaxDbl(:) = [100.0d0, 1000000.0d0, 1000000.0d0, 1000000.0d0, 100.0d0]
1787 ! Loop through and calculate missing/fill/min/max values that will be placed
1788 ! into the NetCDF attributes after scale_factor/add_offset are applied.
1789 do i=1,numRtDomainVars
1790 rtDomainDict%fillComp(i) = &
1791 nint((rtDomainDict%fillReal(i) + rtDomainDict%addOffset(i)), 8) * &
1792 nint(one_dbl / rtDomainDict%scaleFactor(i), 8)
1793 rtDomainDict%missingComp(i) = &
1794 nint((rtDomainDict%missingReal(i) + rtDomainDict%addOffset(i)), 8) * &
1795 nint(one_dbl / rtDomainDict%scaleFactor(i), 8)
1796 rtDomainDict%validMinComp(i) = &
1797 nint((rtDomainDict%validMinDbl(i) + rtDomainDict%addOffset(i)) * &
1798 nint(one_dbl / rtDomainDict%scaleFactor(i), 8), 8)
1799 rtDomainDict%validMaxComp(i) = &
1800 nint((rtDomainDict%validMaxDbl(i) + rtDomainDict%addOffset(i)) * &
1801 nint(one_dbl / rtDomainDict%scaleFactor(i), 8), 8)
1804 end subroutine initRtDomainDict
1806 subroutine initLakeDict(lakeOutDict,procId,diagFlag)
1807 use config_base, only: nlst
1811 type(lakeMeta), intent(inout) :: lakeOutDict
1812 integer, intent(inout) :: diagFlag
1813 integer, intent(inout) :: procId
1816 integer :: ftnMeta,iret
1817 integer :: projVarId
1819 lakeOutDict%modelNdv = -9.E15
1821 ! NOTE !!!!! If you see PLC, this means OWP has no desire to output these,
1822 ! which means meta-data standards have yet to be determined
1823 ! for these variables. Fill in if it's desired to output....
1824 ! First establish global attributes for the channel output files
1825 lakeOutDict%title = "OUTPUT FROM " // trim(get_model_version())
1826 lakeOutDict%fType = 'timeSeries'
1827 lakeOutDict%initTime = '1970-01-01_00:00:00' ! This will be calculated in I/O code
1828 lakeOutDict%validTime = '1970-01-01_00:00:00' ! This will be calculated in I/O code
1829 lakeOutDict%lakeDim = 'lake_id'
1830 lakeOutDict%cdm = 'PLACEHOLDER'
1831 lakeOutDict%conventions = 'CF-1.6'
1833 ! Pull spatial metadata information about the modeling domain from the land
1834 ! spatial metadata file.
1835 if(procId .eq. 0) then
1836 iret = nf90_open(trim(nlst(1)%land_spatial_meta_flnm),NF90_NOWRITE,ncid=ftnMeta)
1837 if(iret .ne. 0) then
1838 ! Spatial metadata file not found for land grid.
1839 call postDiagMsg(diagFlag,'WARNING: Unable to open LAND spatial metadata file.')
1840 lakeOutDict%proj4 = ''
1841 lakeOutDict%esri = ''
1843 ! First pull metadata on coordinate system.
1844 iret = nf90_inq_varid(ftnMeta,'crs',projVarId)
1845 if(iret .ne. 0) then
1846 call postDiagMsg(diagFlag,'WARNING: Unable to find crs in LAND spatial metadata file')
1847 lakeOutDict%proj4 = ''
1848 lakeOutDict%esri = ''
1850 iret = nf90_get_att(ftnMeta,projVarId,'esri_pe_string',lakeOutDict%esri)
1851 if(iret .ne. 0) then
1852 call postDiagMsg(diagFlag,'WARNING: Unable to find esri_pe_string in LAND spatial metadata file.')
1853 lakeOutDict%esri = ''
1855 iret = nf90_get_att(ftnMeta,NF90_GLOBAL,'proj4',lakeOutDict%proj4)
1856 ! We are going to put a relaxed constraint on the proj4 string.
1857 if(iret .ne. 0) then
1858 call postDiagMsg(diagFlag,'WARNING: proj4 string not found. Defaulting to blank string.')
1859 lakeOutDict%proj4 = ''
1863 iret = nf90_close(ftnMeta)
1864 call nwmCheck(diagFlag,iret,'ERROR: Unable to close LAND spatial metadata file.')
1868 ! Next establish time attribues
1869 lakeOutDict%timeLName = 'valid output time'
1870 lakeOutDict%timeUnits = 'minutes since 1970-01-01 00:00:00 UTC'
1871 lakeOutDict%timeStName = 'time'
1872 lakeOutDict%rTimeLName = 'model initialization time'
1873 lakeOutDict%rTimeStName = 'forecast_reference_time'
1874 lakeOutDict%rTimeUnits = 'minutes since 1970-01-01 00:00:00 UTC'
1876 ! Establish elevation variable attributes
1877 lakeOutDict%elevLName = "Water Surface Elevation"
1878 lakeOutDict%elevUnits = "m" ! meters
1879 lakeOutDict%elevComment = "If reservoir_type = 4, " &
1880 // "water_sfc_elev is invalid because this value corresponds only to level pool"
1882 ! Establish feature_id attributes
1883 lakeOutDict%featureIdLName = "Lake ComID"
1884 lakeOutDict%featureIdComment = "ComID from NHDPlusV2 waterbody layer"
1885 lakeOutDict%cfRole = 'timeseries_id'
1887 ! Establish reservoir_type attributes
1888 lakeOutDict%reservoirTypeLName = "reservoir_type"
1889 lakeOutDict%reservoirTypeFlagValues = [1,2,3,4]
1890 lakeOutDict%reservoirTypeFlagMeanings = "Level_pool USGS-persistence USACE-persistence RFC-forecasts"
1892 ! Establish reservoir_assimilated_value attributes
1893 lakeOutDict%reservoirAssimilatedValueLName = "reservoir_assimilated_value"
1894 lakeOutDict%reservoirAssimilatedValueUnits = "m3 s-1"
1896 ! Establish reservoir_assimilated_source_file attributes
1897 lakeOutDict%reservoirAssimilatedSourceFileLName= "reservoir_assimilated_source_file"
1899 ! Esatablish lat/lon attributes
1900 lakeOutDict%latLName = "Lake latitude"
1901 lakeOutDict%latUnits = "degrees_north"
1902 lakeOutDict%latStName = "latitude"
1903 lakeOutDict%lonLName = "Lake longitude"
1904 lakeOutDict%lonUnits = "degrees_east"
1905 lakeOutDict%lonStName = "longitude"
1907 lakeOutDict%varNames(:) = [character(len=64) :: 'inflow','outflow']
1908 lakeOutDict%longName(:) = [character(len=64) :: 'Lake Inflow','Lake Outflow']
1909 lakeOutDict%units(:) = [character(len=64) :: 'm3 s-1','m3 s-1']
1910 lakeOutDict%coordNames(:) = [character(len=64) :: 'latitude longitude','latitude longitude']
1911 lakeOutDict%scaleFactor(:) = [0.01,0.01]
1912 lakeOutDict%addOffset(:) = [0.0,0.0]
1913 lakeOutDict%outFlag(:) = [0,0]
1914 lakeOutDict%timeZeroFlag(:) = [1,1]
1915 lakeOutDict%fillReal(:) = [-9999.0,-9999.0]
1916 lakeOutDict%missingReal(:) = [-9999.0,-9999.0]
1917 lakeOutDict%validMinDbl(:) = [-10000.0d0, -10000.0d0]
1918 lakeOutDict%validMaxDbl(:) = [10000.0d0, 10000.0d0]
1920 ! Loop through and calculate missing/fill/min/max values that will be placed
1921 ! into the NetCDF attributes after scale_factor/add_offset are applied.
1923 lakeOutDict%fillComp(i) = &
1924 nint((lakeOutDict%fillReal(i) + lakeOutDict%addOffset(i)), 8) * &
1925 nint(one_dbl / lakeOutDict%scaleFactor(i), 8)
1926 lakeOutDict%missingComp(i) = &
1927 nint((lakeOutDict%missingReal(i) + lakeOutDict%addOffset(i)), 8) * &
1928 nint(one_dbl / lakeOutDict%scaleFactor(i), 8)
1929 lakeOutDict%validMinComp(i) = &
1930 nint((lakeOutDict%validMinDbl(i) + lakeOutDict%addOffset(i)) * &
1931 nint(one_dbl / lakeOutDict%scaleFactor(i), 8), 8)
1932 lakeOutDict%validMaxComp(i) = &
1933 nint((lakeOutDict%validMaxDbl(i) + lakeOutDict%addOffset(i)) * &
1934 nint(one_dbl / lakeOutDict%scaleFactor(i), 8), 8)
1937 end subroutine initLakeDict
1939 subroutine initChrtGrdDict(chrtGrdDict,procId,diagFlag)
1940 use config_base, only: nlst
1944 type(chrtGrdMeta), intent(inout) :: chrtGrdDict
1945 integer, intent(inout) :: procId,diagFlag
1946 integer :: ftnMeta,projVarId,xVarId,yVarId,ftnGeo
1947 integer :: xDimId,numColLand,numColHydro
1948 real :: resLand,resHydro,aggFactor
1950 integer :: crsRealAttCnt,xRealAttCnt,yRealAttCnt
1951 integer :: crsCharAttCnt,xCharAttCnt,yCharAttCnt
1952 integer :: i, nCrsAtts,nxAtts,nyAtts
1954 integer :: floatFlag
1955 character(len=512) :: tmpAttName
1960 chrtGrdDict%modelNdv = -9.E15
1962 ! First establish global attributes.
1963 chrtGrdDict%title = "OUTPUT FROM " // trim(get_model_version())
1964 chrtGrdDict%initTime = "1970-01-01_00:00:00" ! This will be calculated in I/O code.
1965 chrtGrdDict%validTime = "1970-01-01_00:00:00" ! This will be calculated in I/O code.
1966 chrtGrdDict%decimation = 1
1967 chrtGrdDict%conventions = "CF-1.6"
1969 ! Next establish time attributes
1970 chrtGrdDict%timeLName = "valid output time"
1971 chrtGrdDict%timeUnits = "minutes since 1970-01-01 00:00:00 UTC"
1972 chrtGrdDict%timeStName = "time"
1973 chrtGrdDict%rTimeLName = "model initialization time"
1974 chrtGrdDict%rTimeUnits = "minutes since 1970-01-01 00:00:00 UTC"
1975 chrtGrdDict%rTimeStName = "forecast_reference_time"
1984 ! Pull spatial metadata information about the modeling domain from the
1986 if(procId .eq. 0) then
1987 iret = nf90_open(trim(nlst(1)%geo_finegrid_flnm),NF90_NOWRITE,ncid=ftnMeta)
1988 if(iret .ne. 0) then
1989 ! Spatial metadata file not found for hydro grid.
1990 call postDiagMsg(diagFlag,'WARNING: Unable to open Fulldom file. No crs variable or attributes will be created.')
1991 chrtGrdDict%nCrsCharAtts = 0
1992 chrtGrdDict%nCrsRealAtts = 0
1993 chrtGrdDict%nxCharAtts = 0
1994 chrtGrdDict%nxRealAtts = 0
1995 chrtGrdDict%nyCharAtts = 0
1996 chrtGrdDict%nyRealAtts = 0
1997 chrtGrdDict%proj4 = ''
1999 iret = nf90_get_att(ftnMeta,NF90_GLOBAL,'proj4',chrtGrdDict%proj4)
2000 if(iret .ne. 0) then
2001 call postDiagMsg(diagFlag,'WARNING: Unable to find proj4 global attribute. Defaulting to blank string.')
2002 chrtGrdDict%proj4 = ''
2006 ! Find the crs variable and pull out the attributes, their names, and
2007 ! their values. This will be translated to output files.
2008 iret = nf90_inq_varid(ftnMeta,'crs',projVarId)
2009 ! For now we are going to allow the code to move forward without
2010 ! finding this variable. In the future, we will probably restrict the
2011 ! code to ensure things are more seamless.
2013 if(iret .ne. 0) then
2014 call postDiagMsg(diagFlag,'WARNING: Unable to locate the crs variable. No crs variable or attributes will be created.')
2015 chrtGrdDict%nCrsCharAtts = 0
2016 chrtGrdDict%nCrsRealAtts = 0
2018 iret = nf90_inquire_variable(ftnMeta,projVarId,nAtts=nCrsAtts)
2019 call nwmCheck(diagFlag,iret,'ERROR: Unable to find crs number of attributes')
2021 iret = nf90_inq_attname(ftnMeta,projVarId,i,name=tmpAttName)
2022 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from crs variable.')
2023 iret = nf90_inquire_attribute(ftnMeta,projVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
2024 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from crs variable.')
2025 select case (xtypeTmp)
2047 if(charFlag .eq. 1) then
2048 chrtGrdDict%crsCharAttNames(crsCharAttCnt+1) = trim(tmpAttName)
2049 iret = nf90_get_att(ftnMeta,projVarId,trim(tmpAttName),chrtGrdDict%crsCharAttVals(crsCharAttCnt+1))
2050 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull crs attributes')
2051 crsCharAttCnt = crsCharAttCnt + 1
2053 chrtGrdDict%crsFloatAttNames(crsRealAttCnt+1) = trim(tmpAttName)
2054 chrtGrdDict%crsRealAttLen(crsRealAttCnt+1) = tmpLen
2055 iret = nf90_get_att(ftnMeta,projVarId,trim(tmpAttName),chrtGrdDict%crsRealAttVals(crsRealAttCnt+1,1:tmpLen))
2056 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull crs attributes')
2057 crsRealAttCnt = crsRealAttCnt + 1
2062 chrtGrdDict%nCrsRealAtts = crsRealAttCnt
2063 chrtGrdDict%nCrsCharAtts = crsCharAttCnt
2066 ! Next pull the attributes from the x/y dimensions
2069 iret = nf90_inq_varid(ftnMeta,'x',xVarId)
2070 if(iret .ne. 0) then
2071 call postDiagMsg(diagFlag,'WARNING: Unable to locate the x variable. No x variable or attributes will be created.')
2072 chrtGrdDict%nxCharAtts = 0
2073 chrtGrdDict%nxRealAtts = 0
2075 iret = nf90_inquire_variable(ftnMeta,xVarId,nAtts=nxAtts)
2076 call nwmCheck(diagFlag,iret,'ERROR: Unable to find x number of attributes')
2078 iret = nf90_inq_attname(ftnMeta,xVarId,i,name=tmpAttName)
2079 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from x variable.')
2080 iret = nf90_inquire_attribute(ftnMeta,xVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
2081 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from x variable.')
2082 select case (xtypeTmp)
2104 if(charFlag .eq. 1) then
2105 chrtGrdDict%xCharAttNames(xCharAttCnt+1) = trim(tmpAttName)
2106 iret = nf90_get_att(ftnMeta,xVarId,trim(tmpAttName),chrtGrdDict%xCharAttVals(xCharAttCnt+1))
2107 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull x attributes')
2108 xCharAttCnt = xCharAttCnt + 1
2110 chrtGrdDict%xFloatAttNames(xRealAttCnt+1) = trim(tmpAttName)
2111 chrtGrdDict%xRealAttLen(xRealAttCnt+1) = tmpLen
2112 iret = nf90_get_att(ftnMeta,xVarId,trim(tmpAttName),chrtGrdDict%xRealAttVals(xRealAttCnt+1,1:tmpLen))
2113 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull x attributes')
2114 xRealAttCnt = xRealAttCnt + 1
2119 chrtGrdDict%nxRealAtts = xRealAttCnt
2120 chrtGrdDict%nxCharAtts = xCharAttCnt
2125 iret = nf90_inq_varid(ftnMeta,'y',yVarId)
2126 if(iret .ne. 0) then
2127 call postDiagMsg(diagFlag,'WARNING: Unable to locate the y variable. No y variable or attributes will be created.')
2128 chrtGrdDict%nyCharAtts = 0
2129 chrtGrdDict%nyRealAtts = 0
2131 iret = nf90_inquire_variable(ftnMeta,yVarId,nAtts=nyAtts)
2132 call nwmCheck(diagFlag,iret,'ERROR: Unable to find y number of attributes')
2134 iret = nf90_inq_attname(ftnMeta,yVarId,i,name=tmpAttName)
2135 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from y variable.')
2136 iret = nf90_inquire_attribute(ftnMeta,yVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
2137 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from y variable.')
2138 select case (xtypeTmp)
2160 if(charFlag .eq. 1) then
2161 chrtGrdDict%yCharAttNames(yCharAttCnt+1) = trim(tmpAttName)
2162 iret = nf90_get_att(ftnMeta,yVarId,trim(tmpAttName),chrtGrdDict%yCharAttVals(yCharAttCnt+1))
2163 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull y attributes')
2164 yCharAttCnt = yCharAttCnt + 1
2166 chrtGrdDict%yFloatAttNames(yRealAttCnt+1) = trim(tmpAttName)
2167 chrtGrdDict%yRealAttLen(yRealAttCnt+1) = tmpLen
2168 iret = nf90_get_att(ftnMeta,yVarId,trim(tmpAttName),chrtGrdDict%yRealAttVals(yRealAttCnt+1,1:tmpLen))
2169 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull y attributes')
2170 yRealAttCnt = yRealAttCnt + 1
2175 chrtGrdDict%nyRealAtts = yRealAttCnt
2176 chrtGrdDict%nyCharAtts = yCharAttCnt
2179 iret = nf90_close(ftnMeta)
2180 call nwmCheck(diagFlag,iret,'ERROR: Unable to close LAND spatial metadata file.')
2184 ! Next get the number of columns on the land grid. This will be used to
2185 ! calculate the resolution of the routing grid in meters and aggfactor.
2186 iret = nf90_open(trim(nlst(1)%land_spatial_meta_flnm),NF90_NOWRITE,ncid=ftnGeo)
2187 if(iret .ne. 0) then
2188 call postDiagMsg(diagFlag,'WARNING: Unable to open LAND Spatial Metadata file. Defaulting to missing values.')
2192 iret = nf90_inq_dimid(ftnGeo,'x',xDimId)
2193 call nwmCheck(diagFlag,iret,'ERROR: Unable to find x dimension in LAND spatial Metadata file')
2194 iret = nf90_inquire_dimension(ftnGeo,xDimId,len=numColLand)
2195 call nwmCheck(diagFlag,iret,'ERROR: Unable to retrieve number of columns in the LAND spatial metadata file')
2196 iret = nf90_inq_varid(ftnGeo,'x',xVarId)
2197 call nwmCheck(diagFlag,iret,'ERROR: Unable to find x variable in LAND spatial metadata file')
2198 iret = nf90_get_att(ftnGeo,xVarId,'resolution',resLand)
2199 call nwmCheck(diagFlag,iret,'ERROR: Unable to get x resolution in LAND spatial metadata file')
2200 iret = nf90_close(ftnGeo)
2201 call nwmCheck(diagFlag,iret,'ERROR: Unable to close the LAND spatial metadata file')
2204 ! Next get the number of columns on the high-resolution routing grid.
2205 ! This will be used to calculate the resolution of the routing grid in
2207 iret = nf90_open(trim(nlst(1)%geo_finegrid_flnm),NF90_NOWRITE,ncid=ftnGeo)
2208 if(iret .ne. 0) then
2209 call nwmCheck(diagFlag,iret,'ERROR: Unable to open Fulldom file')
2211 iret = nf90_inq_dimid(ftnGeo,'x',xDimId)
2212 call nwmCheck(diagFlag,iret,'ERROR: Unable to find x dimension in Fulldom file')
2213 iret = nf90_inquire_dimension(ftnGeo,xDimId,len=numColHydro)
2214 call nwmCheck(diagFlag,iret,'ERROR: Unable to retrieve number of columns in the Fulldom file')
2215 iret = nf90_close(ftnGeo)
2216 call nwmCheck(diagFlag,iret,'ERROR: Unable to close the Fulldom file')
2219 ! Calculate the aggregation factor and resolution of the hydro routing
2221 if(numColLand .ne. -9999) then
2222 aggFactor = float(numColHydro/numColLand)
2223 resHydro = resLand/aggFactor
2228 chrtGrdDict%xRes = resHydro
2229 chrtGrdDict%yRes = resHydro
2233 chrtGrdDict%varNames(:) = [character(len=64) :: "streamflow"]
2234 chrtGrdDict%longName(:) = [character(len=64) :: "River Flow"]
2235 chrtGrdDict%units(:) = [character(len=64) :: "m3 s-1"]
2236 chrtGrdDict%scaleFactor(:) = [0.1]
2237 chrtGrdDict%addOffset(:) = [0.0]
2238 chrtGrdDict%outFlag(:) = [0]
2239 chrtGrdDict%timeZeroFlag(:) = [1]
2240 chrtGrdDict%missingReal(:) = [-9999.0]
2241 chrtGrdDict%fillReal(:) = [-9999.0]
2242 chrtGrdDict%validMinDbl(:) = [0.0d0]
2243 chrtGrdDict%validMaxDbl(:) = [50000.0d0]
2245 ! Loop through and calculate missing/fill/min/max values that will be placed
2246 ! into the NetCDF attributes after scale_factor/add_offset are applied.
2248 chrtGrdDict%fillComp(i) = &
2249 nint((chrtGrdDict%fillReal(i) + chrtGrdDict%addOffset(i)), 8) * &
2250 nint(one_dbl / chrtGrdDict%scaleFactor(i), 8)
2251 chrtGrdDict%missingComp(i) = &
2252 nint((chrtGrdDict%missingReal(i) + chrtGrdDict%addOffset(i)), 8) * &
2253 nint(one_dbl / chrtGrdDict%scaleFactor(i), 8)
2254 chrtGrdDict%validMinComp(i) = &
2255 nint((chrtGrdDict%validMinDbl(i) + chrtGrdDict%addOffset(i)) * &
2256 nint(one_dbl / chrtGrdDict%scaleFactor(i), 8), 8)
2257 chrtGrdDict%validMaxComp(i) = &
2258 nint((chrtGrdDict%validMaxDbl(i) + chrtGrdDict%addOffset(i)) * &
2259 nint(one_dbl / chrtGrdDict%scaleFactor(i), 8), 8)
2262 end subroutine initChrtGrdDict
2264 subroutine initLsmOutDict(lsmOutDict,procId,diagFlag)
2265 use config_base, only: nlst
2269 type(lsmMeta), intent(inout) :: lsmOutDict
2270 integer, intent(inout) :: procId
2271 integer, intent(inout) :: diagFlag
2272 integer :: ftnMeta,projVarId,xVarId,yVarId
2274 integer :: crsRealAttCnt,xRealAttCnt,yRealAttCnt
2275 integer :: crsCharAttCnt,xCharAttCnt,yCharAttCnt
2276 integer :: i, nCrsAtts,nxAtts,nyAtts
2278 integer :: floatFlag
2279 character(len=512) :: tmpAttName
2284 lsmOutDict%numSoilLayers = nlst(1)%nsoil
2285 lsmOutDict%act_lev = nlst(1)%act_lev
2287 lsmOutDict%modelNdv = 9.9692099683868690E36
2288 lsmOutDict%modelNdvInt = -2147483647
2290 ! First establish global attributes.
2291 lsmOutDict%title = "OUTPUT FROM " // trim(get_model_version())
2292 lsmOutDict%initTime = "1970-01-01_00:00:00" ! This will be calculated in I/O code.
2293 lsmOutDict%validTime = "1970-01-01_00:00:00" ! This will be calculated in I/O code.
2294 lsmOutDict%conventions = "CF-1.6"
2296 ! Next establish time attributes
2297 lsmOutDict%timeLName = "valid output time"
2298 lsmOutDict%timeUnits = "minutes since 1970-01-01 00:00:00 UTC"
2299 lsmOutDict%timeStName = "time"
2300 lsmOutDict%rTimeLName = "model initialization time"
2301 lsmOutDict%rTimeUnits = "minutes since 1970-01-01 00:00:00 UTC"
2302 lsmOutDict%rTimeStName = "forecast_reference_time"
2311 ! Pull spatial metadata information about the modeling domain from the land
2312 ! spatial metadata file.
2313 if(procId .eq. 0) then
2314 iret = nf90_open(trim(nlst(1)%land_spatial_meta_flnm),NF90_NOWRITE,ncid=ftnMeta)
2315 if(iret .ne. 0) then
2316 ! Spatial metadata file not found for land grid.
2317 call postDiagMsg(diagFlag,'WARNING: Unable to open LAND spatial metadata file. No crs variable or attributes will be created.')
2318 lsmOutDict%nCrsCharAtts = 0
2319 lsmOutDict%nCrsRealAtts = 0
2320 lsmOutDict%nxCharAtts = 0
2321 lsmOutDict%nxRealAtts = 0
2322 lsmOutDict%nyCharAtts = 0
2323 lsmOutDict%nyRealAtts = 0
2324 lsmOutDict%proj4 = ''
2326 iret = nf90_get_att(ftnMeta,NF90_GLOBAL,'proj4',lsmOutDict%proj4)
2327 if(iret .ne. 0) then
2328 call postDiagMsg(diagFlag,'WARNING: Unable to find proj4 global attribute. Defaulting to blank string.')
2329 lsmOutDict%proj4 = ''
2333 ! Find the crs variable and pull out the attributes, their names, and
2334 ! their values. This will be translated to output files.
2335 iret = nf90_inq_varid(ftnMeta,'crs',projVarId)
2336 ! For now we are going to allow the code to move forward without
2337 ! finding this variable. In the future, we will probably restrict the
2338 ! code to ensure things are more seamless.
2340 ! code to ensure things are more seamless.
2341 if(iret .ne. 0) then
2342 call postDiagMsg(diagFlag,'WARNING: Unable to locate the crs variable. No crs variable or attributes will be created.')
2343 lsmOutDict%nCrsCharAtts = 0
2344 lsmOutDict%nCrsRealAtts = 0
2346 iret = nf90_inquire_variable(ftnMeta,projVarId,nAtts=nCrsAtts)
2347 call nwmCheck(diagFlag,iret,'ERROR: Unable to find crs number of attributes')
2349 iret = nf90_inq_attname(ftnMeta,projVarId,i,name=tmpAttName)
2350 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from crs variable.')
2351 iret = nf90_inquire_attribute(ftnMeta,projVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
2352 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from crs variable.')
2353 select case (xtypeTmp)
2375 if(charFlag .eq. 1) then
2376 lsmOutDict%crsCharAttNames(crsCharAttCnt+1) = trim(tmpAttName)
2377 iret = nf90_get_att(ftnMeta,projVarId,trim(tmpAttName),lsmOutDict%crsCharAttVals(crsCharAttCnt+1))
2378 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull crs attributes')
2379 crsCharAttCnt = crsCharAttCnt + 1
2381 lsmOutDict%crsFloatAttNames(crsRealAttCnt+1) = trim(tmpAttName)
2382 lsmOutDict%crsRealAttLen(crsRealAttCnt+1) = tmpLen
2383 iret = nf90_get_att(ftnMeta,projVarId,trim(tmpAttName),lsmOutDict%crsRealAttVals(crsRealAttCnt+1,1:tmpLen))
2384 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull crs attributes')
2385 crsRealAttCnt = crsRealAttCnt + 1
2390 lsmOutDict%nCrsRealAtts = crsRealAttCnt
2391 lsmOutDict%nCrsCharAtts = crsCharAttCnt
2395 ! Next pull the attributes from the x/y dimensions
2398 iret = nf90_inq_varid(ftnMeta,'x',xVarId)
2399 if(iret .ne. 0) then
2400 call postDiagMsg(diagFlag,'WARNING: Unable to locate the x variable. No x variable or attributes will be created.')
2401 lsmOutDict%nxCharAtts = 0
2402 lsmOutDict%nxRealAtts = 0
2404 iret = nf90_inquire_variable(ftnMeta,xVarId,nAtts=nxAtts)
2405 call nwmCheck(diagFlag,iret,'ERROR: Unable to find x number of attributes')
2407 iret = nf90_inq_attname(ftnMeta,xVarId,i,name=tmpAttName)
2408 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from x variable.')
2409 iret = nf90_inquire_attribute(ftnMeta,xVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
2410 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from x variable.')
2411 select case (xtypeTmp)
2433 if(charFlag .eq. 1) then
2434 lsmOutDict%xCharAttNames(xCharAttCnt+1) = trim(tmpAttName)
2435 iret = nf90_get_att(ftnMeta,xVarId,trim(tmpAttName),lsmOutDict%xCharAttVals(xCharAttCnt+1))
2436 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull x attributes')
2437 xCharAttCnt = xCharAttCnt + 1
2439 lsmOutDict%xFloatAttNames(xRealAttCnt+1) = trim(tmpAttName)
2440 lsmOutDict%xRealAttLen(xRealAttCnt+1) = tmpLen
2441 iret = nf90_get_att(ftnMeta,xVarId,trim(tmpAttName),lsmOutDict%xRealAttVals(xRealAttCnt+1,1:tmpLen))
2442 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull x attributes')
2443 xRealAttCnt = xRealAttCnt + 1
2448 lsmOutDict%nxRealAtts = xRealAttCnt
2449 lsmOutDict%nxCharAtts = xCharAttCnt
2455 iret = nf90_inq_varid(ftnMeta,'y',yVarId)
2456 if(iret .ne. 0) then
2457 call postDiagMsg(diagFlag,'WARNING: Unable to locate the y variable. No y variable or attributes will be created.')
2458 lsmOutDict%nyCharAtts = 0
2459 lsmOutDict%nyRealAtts = 0
2461 iret = nf90_inquire_variable(ftnMeta,yVarId,nAtts=nyAtts)
2462 call nwmCheck(diagFlag,iret,'ERROR: Unable to find y number of attributes')
2464 iret = nf90_inq_attname(ftnMeta,yVarId,i,name=tmpAttName)
2465 call nwmCheck(diagFlag,iret,'ERROR: Unable to extract attribute names from y variable.')
2466 iret = nf90_inquire_attribute(ftnMeta,yVarId,trim(tmpAttName),xtype=xtypeTmp,len=tmpLen)
2467 call nwmCheck(diagFlag,iret,'ERROR: Unable to find attribute types from y variable.')
2468 select case (xtypeTmp)
2490 if(charFlag .eq. 1) then
2491 lsmOutDict%yCharAttNames(yCharAttCnt+1) = trim(tmpAttName)
2492 iret = nf90_get_att(ftnMeta,yVarId,trim(tmpAttName),lsmOutDict%yCharAttVals(yCharAttCnt+1))
2493 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull y attributes')
2494 yCharAttCnt = yCharAttCnt + 1
2496 lsmOutDict%yFloatAttNames(yRealAttCnt+1) = trim(tmpAttName)
2497 lsmOutDict%yRealAttLen(yRealAttCnt+1) = tmpLen
2498 iret = nf90_get_att(ftnMeta,yVarId,trim(tmpAttName),lsmOutDict%yRealAttVals(yRealAttCnt+1,1:tmpLen))
2499 call nwmCheck(diagFlag,iret,'ERROR: Unable to pull y attributes')
2500 yRealAttCnt = yRealAttCnt + 1
2505 lsmOutDict%nyRealAtts = yRealAttCnt
2506 lsmOutDict%nyCharAtts = yCharAttCnt
2510 iret = nf90_close(ftnMeta)
2511 call nwmCheck(diagFlag,iret,'ERROR: Unable to close LAND spatial metadata file.')
2516 lsmOutDict%varNames(:) = [character(len=64) :: "stc1","smc1","sh2ox1","stc2",&
2517 "smc2","sh2ox2","stc3","smc3","sh2ox3","stc4",&
2518 "smc4","sh2ox4","infxsrt","sfcheadrt"]
2519 lsmOutDict%longName(:) = [character(len=64) :: "Soil temperature in the top layer",&
2520 "Soil moisture in the top layer",&
2521 "Volumetric soil moisture in the top layer",&
2522 "Soil temperature in the second layer",&
2523 "Soil moisture in the second layer",&
2524 "Volumetric soil moisture in the second layer",&
2525 "Soil temperature in the third layer",&
2526 "Soil moisture in the third layer",&
2527 "Volumetric soil moisture in the third layer",&
2528 "Soil temperature in the fourth layer",&
2529 "Soil moisture in the fourth layer",&
2530 "Volumetric soil moisture in the fourth layer",&
2531 "Infiltration excess","Surface head"]
2532 lsmOutDict%units(:) = [character(len=64) :: "K","fraction","fraction",&
2533 "K","fraction","fraction","K","fraction",&
2534 "fraction","K","fraction","fraction",&
2536 lsmOutDict%scaleFactor(:) = [0.1,0.01,0.01,0.1,0.01,0.01,0.1,0.01,0.01,&
2537 0.1,0.01,0.01,1.0,1.0]
2538 lsmOutDict%addOffset(:) = [0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,&
2540 lsmOutDict%timeZeroFlag(:) = [1,1,1,1,1,1,1,1,1,1,1,1,1,1]
2541 lsmOutDict%numLev(:) = [1,1,1,1,1,1,1,1,1,1,1,1,1,1]
2542 lsmOutDict%missingReal(:) = [-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,&
2543 -9999.0,-9999.0,-9999.0,-9999.0]
2544 lsmOutDict%fillReal(:) = [-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,-9999.0,&
2545 -9999.0,-9999.0,-9999.0,-9999.0]
2546 lsmOutDict%validMinDbl(:) = [150.0d0, 0.0d0, 0.0d0, 150.0d0, 0.0d0, 0.0d0, 150.0d0, 0.0d0, 0.0d0, &
2547 150.0d0, 0.0d0, 0.0d0, 0.0d0, 0.0d0]
2548 lsmOutDict%validMaxDbl(:) = [400.0d0, 1.0d0, 1.0d0, 400.0d0, 1.0d0, 1.0d0, 400.0d0, 1.0d0, 1.0d0, &
2549 400.0d0, 1.0d0, 1.0d0, 100000.0d0, 100000.0d0]
2550 ! Loop through and calculate missing/fill/min/max values that will be placed
2551 ! into the NetCDF attributes after scale_factor/add_offset are applied.
2553 lsmOutDict%fillComp(i) = &
2554 nint((lsmOutDict%fillReal(i) + lsmOutDict%addOffset(i)), 8) * &
2555 nint(one_dbl / lsmOutDict%scaleFactor(i), 8)
2556 lsmOutDict%missingComp(i) = &
2557 nint((lsmOutDict%missingReal(i) + lsmOutDict%addOffset(i)), 8) * &
2558 nint(one_dbl / lsmOutDict%scaleFactor(i), 8)
2559 lsmOutDict%validMinComp(i) = &
2560 nint((lsmOutDict%validMinDbl(i) + lsmOutDict%addOffset(i)) * &
2561 nint(one_dbl / lsmOutDict%scaleFactor(i), 8), 8)
2562 lsmOutDict%validMaxComp(i) = &
2563 nint((lsmOutDict%validMaxDbl(i) + lsmOutDict%addOffset(i)) * &
2564 nint(one_dbl / lsmOutDict%scaleFactor(i), 8), 8)
2567 end subroutine initLsmOutDict
2569 subroutine initChanObsDict(chObsDict,diagFlag,procId)
2570 use config_base, only: nlst
2574 type(chObsMeta), intent(inout) :: chObsDict
2575 integer, intent(inout) :: diagFlag
2576 integer, intent(inout) :: procId
2579 integer :: ftnMeta, iret
2580 integer :: projVarId
2582 chObsDict%modelNdv = -9.E15
2585 ! Pull spatial metadata information about the modeling domain from the land
2586 ! spatial metadata file.
2587 if(procId .eq. 0) then
2588 iret = nf90_open(trim(nlst(1)%land_spatial_meta_flnm),NF90_NOWRITE,ncid=ftnMeta)
2589 if(iret .ne. 0) then
2590 ! Spatial metadata file not found for land grid.
2591 call postDiagMsg(diagFlag,'WARNING: Unable to open LAND spatial metadata file.')
2592 chObsDict%proj4 = ''
2595 ! First pull metadata on coordinate system.
2596 iret = nf90_inq_varid(ftnMeta,'crs',projVarId)
2597 if(iret .ne. 0) then
2598 call postDiagMsg(diagFlag,'WARNING: Unable to find crs in LAND spatial metadata file')
2599 chObsDict%proj4 = ''
2602 iret = nf90_get_att(ftnMeta,projVarId,'esri_pe_string',chObsDict%esri)
2603 if(iret .ne. 0) then
2604 call postDiagMsg(diagFlag,'WARNING: Unable to find esri_pe_string in LAND spatial metadata file.')
2607 iret = nf90_get_att(ftnMeta,NF90_GLOBAL,'proj4',chObsDict%proj4)
2608 ! We are going to put a relaxed constraint on the proj4 string.
2609 if(iret .ne. 0) then
2610 call postDiagMsg(diagFlag,'WARNING: proj4 string not found. Defaulting to blank string.')
2611 chObsDict%proj4 = ''
2615 iret = nf90_close(ftnMeta)
2616 call nwmCheck(diagFlag,iret,'ERROR: Unable to close LAND spatial metadata file.')
2620 chObsDict%title = "OUTPUT FROM " // trim(get_model_version())
2621 chObsDict%fType = 'timeSeries'
2622 chObsDict%initTime = '1970-01-01_00:00:00' ! This will be calculated in I/O code
2623 chObsDict%validTime = '1970-01-01_00:00:00' ! This will be calculated in I/O code
2624 chObsDict%stDim = 'feature_id'
2625 chObsDict%stOrder = 1
2626 chObsDict%cdm = 'Station'
2627 chObsDict%conventions = 'CF-1.6'
2629 ! Next establish time attribues
2630 chObsDict%timeLName = 'valid output time'
2631 chObsDict%timeUnits = 'minutes since 1970-01-01 00:00:00 UTC'
2632 chObsDict%timeStName = 'time'
2633 chObsDict%rTimeLName = 'model initialization time'
2634 chObsDict%rTimeStName = 'forecast_reference_time'
2635 chObsDict%rTimeUnits = 'minutes since 1970-01-01 00:00:00 UTC'
2637 ! Esatablish lat/lon attributes
2638 chObsDict%latLName = "Feature latitude"
2639 chObsDict%latUnits = "degrees_north"
2640 chObsDict%latStName = "latitude"
2641 chObsDict%lonLName = "Feature longitude"
2642 chObsDict%lonUnits = "degrees_east"
2643 chObsDict%lonStName = "longitude"
2645 ! Establish streamflw order attributes
2646 chObsDict%orderLName = "Streamflow Order"
2647 chObsDict%orderStName = "order"
2649 ! Establish point elevation attributes
2650 chObsDict%elevLName = "Feature Elevation"
2651 chObsDict%elevUnits = "meters"
2652 chObsDict%elevStName = "Elevation"
2654 ! Next establish feature_id attributes
2655 chObsDict%featureIdLName = 'Reach ID'
2656 chObsDict%featureIdComment = 'NHDPlusv2 ComIDs within CONUS, arbitrary Reach IDs outside of CONUS'
2657 chObsDict%cfRole = 'timeseries_id'
2659 chObsDict%varNames(:) = [character(len=64) :: "streamflow"]
2660 chObsDict%longName(:) = [character(len=64) :: "River Flow"]
2661 chObsDict%units(:) = [character(len=64) :: "m3 s-1"]
2662 chObsDict%coordNames(:) = [character(len=64) :: "latitude longitude"]
2663 chObsDict%scaleFactor(:) = [0.01]
2664 chObsDict%addOffset(:) = [0.0]
2665 ! Initialize all output flags to 0. Modify (if absolutely necessary) in the
2666 ! output subroutine.
2667 chObsDict%outFlag(:) = [0]
2668 chObsDict%timeZeroFlag(:) = [1]
2669 chObsDict%fillReal(:) = [-9999.0]
2670 chObsDict%missingReal(:) = [-9999.0]
2671 chObsDict%validMinDbl(:) = [0.0d0]
2672 chObsDict%validMaxDbl(:) = [50000.0d0]
2673 ! Loop through and calculate missing/fill/min/max values that will be placed
2674 ! into the NetCDF attributes after scale_factor/add_offset are applied.
2676 chObsDict%fillComp(i) = &
2677 nint((chObsDict%fillReal(i) + chObsDict%addOffset(i)), 8) * &
2678 nint(one_dbl / chObsDict%scaleFactor(i), 8)
2679 chObsDict%missingComp(i) = &
2680 nint((chObsDict%missingReal(i) + chObsDict%addOffset(i)), 8) * &
2681 nint(one_dbl / chObsDict%scaleFactor(i), 8)
2682 chObsDict%validMinComp(i) = &
2683 nint((chObsDict%validMinDbl(i) + chObsDict%addOffset(i)) * &
2684 nint(one_dbl / chObsDict%scaleFactor(i), 8), 8)
2685 chObsDict%validMaxComp(i) = &
2686 nint((chObsDict%validMaxDbl(i) + chObsDict%addOffset(i)) * &
2687 nint(one_dbl / chObsDict%scaleFactor(i), 8), 8)
2690 end subroutine initChanObsDict
2692 subroutine initGwDict(gwOutDict)
2695 type(gwMeta), intent(inout) :: gwOutDict
2698 gwOutDict%modelNdv = -9.E15
2700 gwOutDict%title = "OUTPUT FROM " // trim(get_model_version())
2701 gwOutDict%fType = 'timeSeries'
2702 !gwOutDict%proj4 = '+proj=longlat +datum=NAD83 +no_defs'
2703 gwOutDict%initTime = '1970-01-01_00:00:00' ! This will be calculated in I/O code
2704 gwOutDict%validTime = '1970-01-01_00:00:00' ! This will be calculated in I/O code
2705 gwOutDict%gwDim = 'gw_id'
2706 !gwOutDict%cdm = 'PLACEHOLDER'
2707 !gwOutDict%esri = 'GEOGCS[GCS_North_American_1983,DATUM[D_North_American_1983,&
2708 ! &SPHEROID[GRS_1980,6378137.0,298.257222101]],&
2709 ! &PRIMEM[Greenwich,0.0],UNIT[Degree,0.017453292519943295]]'
2710 gwOutDict%conventions = 'CF-1.6'
2712 ! Next establish time attribues
2713 gwOutDict%timeLName = 'valid output time'
2714 gwOutDict%timeUnits = 'minutes since 1970-01-01 00:00:00 UTC'
2715 gwOutDict%timeStName = 'time'
2716 gwOutDict%rTimeLName = 'model initialization time'
2717 gwOutDict%rTimeStName = 'forecast_reference_time'
2718 gwOutDict%rTimeUnits = 'minutes since 1970-01-01 00:00:00 UTC'
2720 ! Establish elevation variable attributes
2721 !gwOutDict%elevLName = "Water Surface Elevation"
2722 !gwOutDict%elevUnits = "meters"
2724 ! Establish feature_id attributes
2725 gwOutDict%featureIdLName = "Groundwater Bucket ID"
2726 gwOutDict%featureIdComment = "Groundwater Bucket ID"
2727 gwOutDict%cfRole = 'timeseries_id'
2729 ! Esatablish lat/lon attributes
2730 !gwOutDict%latLName = "Groundwater Bucket latitude"
2731 !gwOutDict%latUnits = "degrees_north"
2732 !gwOutDict%latStName = "latitude"
2733 !gwOutDict%lonLName = "Groundwater Bucket longitude"
2734 !gwOutDict%lonUnits = "degrees_east"
2735 !gwOutDict%lonStName = "longitude"
2737 gwOutDict%varNames(:) = [character(len=64) :: 'inflow','outflow','loss','depth']
2738 gwOutDict%longName(:) = [character(len=64) :: 'Bucket Inflow','Bucket Outflow','Bucket Loss','Bucket Depth']
2739 gwOutDict%units(:) = [character(len=64) :: 'm3 s-1','m3 s-1','m3 s-1','mm']
2740 !gwOutDict%coordNames(:) = [character(len=64) :: 'latitude longitude','latitude longitude','latitude longitude','latitude longitude']
2741 gwOutDict%scaleFactor(:) = [0.001,0.001,0.001,0.1]
2742 gwOutDict%addOffset(:) = [0.0,0.0,0.0,0.0]
2743 gwOutDict%outFlag(:) = [0,0,0,0]
2744 gwOutDict%timeZeroFlag(:) = [0,0,0,1]
2745 gwOutDict%fillReal(:) = [-9999.0,-9999.0,-9999.0,-9999.0]
2746 gwOutDict%missingReal(:) = [-9999.0,-9999.0,-9999.0,-9999.0]
2747 gwOutDict%validMinDbl(:) = [0.0d0, 0.0d0, 0.0d0, 0.0d0]
2748 gwOutDict%validMaxDbl(:) = [50000.0d0, 50000.0d0, 50000.0d0, 10000.0d0]
2750 ! Loop through and calculate missing/fill/min/max values that will be placed
2751 ! into the NetCDF attributes after scale_factor/add_offset are applied.
2753 gwOutDict%fillComp(i) = &
2754 nint((gwOutDict%fillReal(i) + gwOutDict%addOffset(i)), 8) * &
2755 nint(one_dbl / gwOutDict%scaleFactor(i), 8)
2756 gwOutDict%missingComp(i) = &
2757 nint((gwOutDict%missingReal(i) + gwOutDict%addOffset(i)), 8) * &
2758 nint(one_dbl / gwOutDict%scaleFactor(i), 8)
2759 gwOutDict%validMinComp(i) = &
2760 nint((gwOutDict%validMinDbl(i) + gwOutDict%addOffset(i)) * &
2761 nint(one_dbl / gwOutDict%scaleFactor(i), 8), 8)
2762 gwOutDict%validMaxComp(i) = &
2763 nint((gwOutDict%validMaxDbl(i) + gwOutDict%addOffset(i)) * &
2764 nint(one_dbl / gwOutDict%scaleFactor(i), 8), 8)
2767 end subroutine initGwDict
2769 subroutine postDiagMsg(diagFlag,diagMsg)
2772 ! Subroutine arguments.
2773 integer, intent(in) :: diagFlag
2774 character(len=*), intent(in) :: diagMsg
2776 ! Only write out message if the diagnostic WRF_HYDRO_D flag was
2778 if (diagFlag .eq. 1) then
2779 print*, trim(diagMsg)
2782 end subroutine postDiagMsg
2784 subroutine nwmCheck(diagFlag,iret,msg)
2787 ! Subroutine arguments.
2788 integer, intent(in) :: diagFlag,iret
2789 character(len=*), intent(in) :: msg
2791 ! Check status. If status of command is not 0, then post the error message
2792 ! if WRF_HYDRO_D was set to be 1.
2793 if (iret .ne. 0) then
2794 call hydro_stop(trim(msg))
2797 end subroutine nwmCheck
2799 end module module_NWM_io_dict