FV3 Bundle
ncdw_chaninfo.F90
Go to the documentation of this file.
1 ! nc_diag_write - NetCDF Layer Diag Writing Module
2 ! Copyright 2015 Albert Huang - SSAI/NASA for NASA GSFC GMAO (610.1).
3 !
4 ! Licensed under the Apache License, Version 2.0 (the "License");
5 ! you may not use this file except in compliance with the License.
6 ! You may obtain a copy of the License at
7 !
8 ! http://www.apache.org/licenses/LICENSE-2.0
9 !
10 ! Unless required by applicable law or agreed to in writing, software
11 ! distributed under the License is distributed on an "AS IS" BASIS,
12 ! WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or
13 ! implied. See the License for the specific language governing
14 ! permissions and limitations under the License.
15 !
16 ! chaninfo module - ncdw_chaninfo
17 !
19  ! Module that provides chaninfo variable storage support.
20  !
21  ! This module has all of the subroutines needed to store chaninfo
22  ! data. It includes the chaninfo storing subroutine
23  ! (nc_diag_chaninfo), subroutines for controlling chaninfo data
24  ! (dimension setting, loading definitions, saving definitions,
25  ! saving data, etc.), and preallocation subroutines.
26  !
27  ! Background:
28  ! chaninfo is a fixed storage variable, with dimensions of
29  ! 1 x nchans, where nchans is a known integer.
30  !
31  ! Because we can know nchans, we can constrain the dimensions and
32  ! make a few assumptions:
33  !
34  ! -> nchans won't change for the duration of the file being open;
35  ! -> nchans will be the same for all chaninfo variables, for any
36  ! type involved;
37  ! -> because everything is fixed, we can store variables block
38  ! by block!
39  !
40  ! Because Fortran is a strongly typed language, we can't do silly
41  ! tricks in C, like allocating some memory to a void pointer and
42  ! just storing our byte, short, int, long, float, or double
43  ! numeric data there, and later casting it back...
44  !
45  ! (e.g. void **data_ref; data_ref = malloc(sizeof(void *) * 1000);
46  ! float *f = malloc(sizeof(float)); *f = 1.2345;
47  ! data_ref[0] = f; ...)
48  !
49  ! No frets - we can work around this issue with some derived types
50  ! and arrays! We create an array for each type we want to support.
51  ! Since we're using kinds.F90, we support the following types:
52  ! i_byte, i_short, i_long, r_single, r_double, character(len=*)
53  !
54  ! The derived type used, diag_chaninfo, has these variables to
55  ! help us keep track of everything:
56  !
57  ! -> ci_* - these arrays have the types listed above, plus string
58  ! support! These arrays are simply arrays that we throw our
59  ! data in. However, don't mistaken "throw in" with
60  ! "disorganized" - chaninfo uses a very structured format for
61  ! these variables! Keep reading to find out how we structure
62  ! it...
63  !
64  ! -> nchans - the number of channels to use. Remember that
65  ! chaninfo variables have dimensions 1 x nchans - basically, we
66  ! need to store nchans values. We'll need this a LOT to do
67  ! consistency checks, and to keep track of everything!
68  !
69  ! -> names - all of the chaninfo variable names! We'll be using
70  ! this array to store and lookup chaninfo variables, as well as
71  ! storing them!
72  !
73  ! -> types - all of the chaninfo variable types! These are byte
74  ! integers that get compared to our NLAYER_* type constants
75  ! (see: ncdw_types.F90).
76  !
77  ! -> var_usage - the amount of entries we've stored in our
78  ! chaninfo variable! For instance, if we called
79  ! nc_diag_chaninfo("myvar", 1) three times, for that particular
80  ! var_usage(i), we would have an entry of 3.
81  !
82  ! -> var_rel_pos - the star of the show! This is an abbreviation
83  ! of "variable relative positioning". Recall that we store
84  ! our variable data in ci_* specific type arrays. We know
85  ! the nchans amount, and we know the type. This variable stores
86  ! the "block" that our data is in within the type array.
87  !
88  ! Therefore, we can use the equation to find our starting
89  ! position: 1 + [(REL_VAL - 1) * nchans]
90  !
91  ! For instance, if var_rel_pos(1) for variable names(1) = "bla"
92  ! is set to 2, nchans is 10, and the type is NLAYER_FLOAT, we
93  ! can deduce that in ci_rsingle, our data can be found starting
94  ! at 1 + (1 * 10) = 11. This makes sense, as seen with our mini
95  ! diagram below:
96  !
97  ! ci_rsingle:
98  ! / ci_rsingle index \
99  ! 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20
100  ! [ x, x, x, x, x, x, x, x, x, x, y, y, y, y, y, y, y, y, y, y]
101  ! \ ci_rsingle array /
102  !
103  ! Indeed, our second block does start at index 11!
104  ! As a bonus, since our data is in blocks, things can be super
105  ! fast since we're just cutting our big array into small ones!
106  !
107  ! -> acount_v: Finally, we have dynamic allocation. We have in our
108  ! type a variable called acount_v. This tells us how many
109  ! variables are stored in each type. Using the same equation
110  ! above, and combining with var_usage, we can figure out where
111  ! to put our data!
112  !
113  ! Assume var_usage(i) = 2, block starts at index 11 with the
114  ! equation above.
115  !
116  ! Again, with our fun little diagram:
117  !
118  ! ci_rsingle:
119  ! / ci_rsingle index \
120  ! 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20
121  ! [ x, x, x, x, x, x, x, x, x, x, y, y, Y, y, y, y, y, y, y, y]
122  ! [ BLOCK 1 SEEK = 1->10->11 ][var_u=2|--block 2 area 11->20]
123  ! \ ci_rsingle array /
124  !
125  ! The capital Y marks the place we store our data!
126  !
127  ! For the non-data variables (e.g. variable names, types, etc.),
128  ! they are indexed by variable insertion order. This allows for
129  ! easy lookup by looking up the variable name, and using the
130  ! resulting index for fetching other information.
131  !
132  ! Example:
133  ! names: [ 'asdf', 'ghjk', 'zxcv' ]
134  ! types: [ BYTE, FLOAT, BYTE ]
135  ! var_rel_pos: [ 1, 1, 2 ]
136  !
137  ! Lookup: "ghjk", result index = 2
138  !
139  ! Therefore, the "ghjk" variable type is types(2) = FLOAT, and
140  ! the var_rel_pos for "ghjk" variable is var_rel_pos(2) = 1.
141  !
142  ! These variables are allocated and reallocated, as needed.
143  !
144  ! For the variable metadata fields (variable names, types,
145  ! relative indicies, etc.), these are reallocated incrementally
146  ! when a new variable is added.
147  !
148  ! For the data storage fields, these are reallocated incrementally
149  ! when new data is added.
150  !
151  ! Initial allocation and subsequent reallocation is done by
152  ! chunks. Allocating one element and/or reallocating and adding
153  ! just one element is inefficient, since it's likely that much
154  ! more data (and variables) will be added. Thus, allocation and
155  ! reallocation is done by (re-)allocating exponentially increasing
156  ! chunk sizes. See nc_diag_chaninfo_allocmulti help for more
157  ! details.
158  !
159  ! Because (re-)allocation is done in chunks, we keep a count of
160  ! how much of the memory we're using so that we know when it's
161  ! time to (re-)allocate. Once we need to (re-)allocate, we
162  ! perform it, and then update our total memory counter to keep
163  ! track of the memory already allocated.
164  !
165  ! With all of these variables (and a few more state variables),
166  ! we can reliably store our chaninfo data quickly and
167  ! efficiently!
168  !
169 
170  ! Load our numerical types from kinds
171  ! Note that i_llong is not a type we store - it's just for keeping
172  ! track of numeric indexes. (Maybe this is too excessive...)
173  use ncd_kinds, only: i_byte, i_short, i_long, i_llong, r_single, &
174  r_double
175 
176  ! Load state variables! We need to know:
177  ! init_done - ...whether a file is currently loaded or
178  ! not.
179  ! ncid - ...the current NCID of our file.
180  ! append_only - ...whether we are in append mode or not.
181  ! enable_trim - ...whether we need to automatically trim
182  ! our strings for chaninfo string storage or
183  ! not.
184  ! diag_chaninfo_store - ...chaninfo variable information.
185  ! We pretty much do everything related to
186  ! chaninfo here, so we're using everything
187  ! inside this derived type!
188  use ncdw_state, only: init_done, ncid, append_only, &
189  enable_trim, &
191 
192  ! Load types! We need:
193  ! NLAYER_* - nc_diag types.
194  ! NLAYER_FILL_* - nc_diag type fill. This is pretty much
195  ! equivalent to NF90_FILL_*.
196  ! NLAYER_COMPRESSION - zlib (a la gzip) compression level setting.
197  ! NLAYER_DEFAULT_ENT - default starting number of element entries.
198  ! This is for the initial allocation of
199  ! space for data storage arrays, e.g.
200  ! the ci_* data arrays within diag_chaninfo.
201  ! NLAYER_MULTI_BASE - the base number to use when exponentiating
202  ! to allocate or reallocate data storage
203  ! arrays.
209 
210  ! Load our varattr adder! We need this to store our new shiny
211  ! variable in the varattr database so we can add variable attributes
212  ! to our variables.
214 
215  ! Load our function - given an array of strings, find
216  ! max(len_trim(str_array)) - aka the maximum for len_trim()s on each
217  ! variable.
219 
220  use ncdw_climsg, only: &
221 #ifdef enable_action_msgs
222  nclayer_enable_action, nclayer_actionm, &
223 #endif
224  nclayer_error, nclayer_warning, nclayer_info, nclayer_check
225 
226  ! Load our fun reallocation subroutine - we need this to reallocate
227  ! a few things in our preallocation subroutines:
228  use ncdw_realloc, only: nc_diag_realloc
229 
230  ! Load our chaninfo resizing subroutines - these resize our data
231  ! storage arrays automatically when needed!
236 
237  use netcdf, only: nf90_inquire, nf90_inq_dimid, &
238  nf90_inquire_dimension, nf90_inquire_variable, nf90_def_dim, &
239  nf90_def_var, nf90_get_var, nf90_put_var, &
240  nf90_def_var_deflate, nf90_def_var_chunking, &
241  nf90_byte, nf90_short, nf90_int, nf90_float, nf90_double, &
242  nf90_char, &
243  nf90_ebaddim, nf90_noerr, nf90_max_name, nf90_chunked
244 
245  implicit none
246 
247  ! Add a single chaninfo value to a new or existing chaninfo
248  ! variable.
249  !
250  ! Given the chaninfo variable name and value, add or update the
251  ! variable with the corresponding value.
252  !
253  ! If the variable doesn't already exist, this will automatically
254  ! create it and store the value into it.
255  !
256  ! If the variable does exist, it will simply append to the
257  ! variable's existing values.
258  !
259  ! Values are inserted in the order of the calls made. As such,
260  ! this subroutine is best designed to be used in a loop, where
261  ! for every channel iteration, a value is added using this
262  ! subroutine.
263  !
264  ! chaninfo is stored element by element - no arrays are accepted,
265  ! only scalar values. The best way to call chaninfo is in a loop,
266  ! where each channel is being accessed and stored.
267  !
268  ! Once a value has been added, it may not be removed. Make sure you
269  ! are certain that the value should be added!
270  !
271  ! The number of values may not exceed the number of channels
272  ! (nchans). If more values are added and nchans is exceeded, an
273  ! error will occur.
274  !
275  ! Data locking and definition locking will also affect adding
276  ! chaninfo variables and value. If data locking is in effect, any
277  ! variable or value adding will not work. If definition locking is
278  ! in effect, adding variable values to existing variables will still
279  ! work, but adding new variables will not.
280  !
281  ! For strings, if the length of the string changes when trimming is
282  ! disabled, or when the definitions have been locked, an error will
283  ! occur as well.
284  !
285  ! To see more details about what checks are made, see the
286  ! corresponding called subroutine documentation for details.
287  !
288  ! Valid data types (represented below as data_types):
289  ! integer(i_byte), integer(i_short), integer(i_long),
290  ! real(r_single), real(r_double), character(len=*)
291  !
292  ! Args:
293  ! name (character(len=*)): the name of the chaninfo variable to
294  ! add or update.
295  ! value (data_types): the value to add to chaninfo.
296  !
297  ! Raises:
298  ! If data writing is locked, this will result in an error.
299  !
300  ! If the variable doesn't exist yet, and definitions are locked,
301  ! this will result in an error.
302  !
303  ! If the amount of data in the chaninfo variable is already at
304  ! or exceeding nchans, this will result in an error.
305  !
306  ! For string data, if the string length changes and the
307  ! definitions have already been locked, this will result in an
308  ! error.
309  !
310  ! Also, for string data, if the string length changes and
311  ! trimming is turned off, this will also result in an error.
312  !
313  ! The following errors will trigger indirectly from other
314  ! subroutines called here:
315  !
316  ! If nchans has not been set yet, this will result in an error.
317  !
318  ! If there is no file open (or the file is already closed),
319  ! this will result in an error.
320  !
321  ! Other errors may result from invalid data storage, NetCDF
322  ! errors, or even a bug. See the called subroutines'
323  ! documentation for details.
324  !
326  module procedure nc_diag_chaninfo_byte, &
330  end interface nc_diag_chaninfo
331 
332  contains
333  ! Set the number of channels (nchans) for chaninfo to use for
334  ! variable storage and configuration.
335  !
336  ! This set the number of channels (nchans) for all of the future
337  ! chaninfo variables that will be added. nchans will be used
338  ! as the number of elements to use for every chaninfo variable
339  ! added. It will also be used as a bounds check for variable
340  ! data amounts.
341  !
342  ! Args:
343  ! nchans (integer(i_long)): the number of channels to use
344  ! for chaninfo.
345  !
346  ! Raises:
347  ! If nchans was already set, this will result in an error.
348  ! (You can't change nchans arbitarily - otherwise, variable
349  ! data amounts could become invalid!)
350  !
351  ! If the nchans specified is invalid (<1), this will result
352  ! in an error. If you have no chaninfo variables to write,
353  ! don't call this subroutine at all. No chaninfo variables
354  ! will be processed or written if you don't set anything!
355  !
356  ! If there is no file open (or the file is already closed),
357  ! this will result in an error.
358  !
359  ! Although unlikely, other errors may indirectly occur.
360  ! They may be general storage errors, or even a bug.
361  ! See the called subroutines' documentation for details.
362  !
363  subroutine nc_diag_chaninfo_dim_set(nchans)
364  integer(i_long), intent(in) :: nchans
365 #ifdef ENABLE_ACTION_MSGS
366  character(len=1000) :: action_str
367 
368  if (nclayer_enable_action) then
369  write(action_str, "(A, I0, A)") "nc_diag_chaninfo_dim_set(nchans = ", nchans, ")"
370  call nclayer_actionm(trim(action_str))
371  end if
372 #endif
373  ! Check if everything is initialized / file is open
374  if (init_done .AND. allocated(diag_chaninfo_store)) then
375  ! nchans can't be less than 1!
376  if (nchans < 1) then
377  call nclayer_error("Critical error - specified a nchan < 1!")
378  end if
379 
380  ! Is nchans already set?
381  if (diag_chaninfo_store%nchans /= -1) &
382  call nclayer_error("nchans already set!")
383 
384  ! Set nchans
385  diag_chaninfo_store%nchans = nchans
386  else
387  call nclayer_error("NetCDF4 layer not initialized yet!")
388  end if
389  end subroutine nc_diag_chaninfo_dim_set
390 
391  ! Set the allocation multiplier for chaninfo variable storage
392  ! allocation and reallocation.
393  !
394  ! This sets the allocation multiplier (exponentiator?) for
395  ! chaninfo variable storage allocation and reallocation.
396  !
397  ! Reallocation looks like this:
398  ! new_size = old_size + addl_num_entries +
399  ! (NLAYER_DEFAULT_ENT * (NLAYER_MULTI_BASE **
400  ! diag_chaninfo_store%alloc_multi))
401  !
402  ! NLAYER_DEFAULT_ENT and NLAYER_MULTI_BASE are constants defined
403  ! in ncdw_types. The alloc_multi part is set with this
404  ! subroutine.
405  !
406  ! As reallocation occurs, the alloc_multi continues to increase
407  ! by one, causing subsequent reallocations to allocate
408  ! exponentially more memory.
409  !
410  ! You can use this subroutine to increase the initial amount of
411  ! memory allocated/reallocated, or you can use it to prevent
412  ! the reallocating counter from increasing by calling this
413  ! every so often.
414  !
415  ! If this is not set, it will be initially set to 0 and will
416  ! increase from there.
417  !
418  ! Args:
419  ! multiplier (integer(i_long)): the multiplier to use when
420  ! allocating or reallocating.
421  !
422  ! Raises:
423  ! If there is no file open (or the file is already closed),
424  ! this will result in an error.
425  !
426  ! Although unlikely, other errors may indirectly occur.
427  ! They may be general storage errors, or even a bug.
428  ! See the called subroutines' documentation for details.
429  !
430  subroutine nc_diag_chaninfo_allocmulti(multiplier)
431  integer(i_long), intent(in) :: multiplier
432 #ifdef ENABLE_ACTION_MSGS
433  character(len=1000) :: action_str
434 
435  if (nclayer_enable_action) then
436  write(action_str, "(A, I0, A)") "nc_diag_chaninfo_allocmulti(multiplier = ", multiplier, ")"
437  call nclayer_actionm(trim(action_str))
438  end if
439 #endif
440  if (init_done) then
441  ! # of times we needed to realloc simple metadata
442  ! also the multiplier factor for allocation (2^x)
443  diag_chaninfo_store%alloc_multi = multiplier
444  end if
445  end subroutine nc_diag_chaninfo_allocmulti
446 
447  ! Load chaninfo variable definitions from an existing, already
448  ! open NetCDF file.
449  !
450  ! This will scan the currently open NetCDF file for chaninfo
451  ! variables. If any exist, the metadata and last position will
452  ! get loaded into the chaninfo variable data buffer.
453  !
454  ! Basically, this scans for the "nchans" dimension. If it
455  ! exists, we set our internal nchans to that dimension's value.
456  ! Then we fetch the dimension names for all variables, and try
457  ! to match them to "nchans". (This is slow... see TODO.txt for
458  ! a better idea!)
459  !
460  ! Once we find our chaninfo variable(s), we scan them for NetCDF
461  ! fill bytes, starting at the end of the variable. When we find
462  ! a spot that does NOT have a fill byte, we set our relative
463  ! index at that spot, and set everything up to resume at that
464  ! position.
465  !
466  ! For string data, we also our maximum string length constraint
467  ! so that we still keep additional variable data within bounds.
468  !
469  ! This is an internal subroutine, and is NOT meant to be called
470  ! outside of nc_diag_write. Calling this subroutine in your
471  ! program may result in unexpected behavior and/or data
472  ! corruption!
473  !
474  ! Args:
475  ! None
476  !
477  ! Raises:
478  ! If the chaninfo variable uses an unsupported type (e.g.
479  ! not one of the types listed above), this will result in
480  ! an error.
481  !
482  ! If there is no file open (or the file is already closed),
483  ! this will result in an error. (NetCDF error here, since
484  ! init_done is not being checked... see TODO.txt)
485  !
486  ! Other errors may result from invalid data storage, NetCDF
487  ! errors, or even a bug. See the called subroutines'
488  ! documentation for details.
489  !
490  subroutine nc_diag_chaninfo_load_def
491  integer(i_long) :: ndims, nvars, var_index, type_index
492  integer(i_long) :: rel_index, i, j
493 
494  ! Temporary variables used when scanning variables and dimensions
495  ! from our NetCDF file
496  character(len=NF90_MAX_NAME) :: tmp_var_name
497  integer(i_long) :: tmp_var_type, tmp_var_ndims
498 
499  integer(i_long), dimension(:), allocatable :: tmp_var_dimids, tmp_var_dim_sizes
500  character(len=NF90_MAX_NAME) , allocatable :: tmp_var_dim_names(:)
501 
502  ! Is this a nchans var?
503  logical :: is_nchans_var
504 
505  ! Data buffers - we need these to fetch our data and see where
506  ! we left off...
507  integer(i_byte), dimension(:), allocatable :: byte_buffer
508  integer(i_short), dimension(:), allocatable :: short_buffer
509  integer(i_long), dimension(:), allocatable :: long_buffer
510 
511  real(r_single), dimension(:), allocatable :: rsingle_buffer
512  real(r_double), dimension(:), allocatable :: rdouble_buffer
513 
514  character(1), dimension(:,:), allocatable :: string_buffer
515 
516  ! Dimension checking NetCDF error storage
517  integer(i_long) :: dim_nc_err
518 
519  ! Get top level info about the file!
520  call nclayer_check(nf90_inquire(ncid, ndimensions = ndims, &
521  nvariables = nvars))
522 
523  ! Fetch nchans first!
524  dim_nc_err = nf90_inq_dimid(ncid, "nchans", diag_chaninfo_store%nchans_dimid)
525 
526  ! Check if we found anything!
527  ! If we got NF90_EBADDIM, then exit.
528  if (dim_nc_err == nf90_ebaddim) then
529  return
530  else if (dim_nc_err /= nf90_noerr) then
531  ! If an error besides not finding the dimension occurs,
532  ! raise an exception.
533  call nclayer_check(dim_nc_err)
534  end if
535 
536  ! Then grab nchans value...
537  call nclayer_check(nf90_inquire_dimension(ncid, diag_chaninfo_store%nchans_dimid, &
538  len = diag_chaninfo_store%nchans))
539 
540  ! Now search for variables that use nchans!
541  ! Loop through each variable!
542  do var_index = 1, nvars
543  ! Grab number of dimensions and attributes first
544  call nclayer_check(nf90_inquire_variable(ncid, var_index, name = tmp_var_name, ndims = tmp_var_ndims))
545 
546  ! Allocate temporary variable dimids storage!
547  allocate(tmp_var_dimids(tmp_var_ndims))
548  allocate(tmp_var_dim_names(tmp_var_ndims))
549  allocate(tmp_var_dim_sizes(tmp_var_ndims))
550 
551  ! Grab the actual dimension IDs and attributes
552  call nclayer_check(nf90_inquire_variable(ncid, var_index, dimids = tmp_var_dimids, &
553  xtype = tmp_var_type))
554 
555  if ((tmp_var_ndims == 1) .OR. &
556  ((tmp_var_ndims == 2) .AND. (tmp_var_type == nf90_char))) then
557  ! Reset our is_nchans_var switch to FALSE!
558  is_nchans_var = .false.
559 
560  ! Fetch all dimension names for the dimensions in the
561  ! variable, and check if the variable is a nchans
562  ! variable. We do so by (slowly) checking all
563  ! dimension names and seeing if they match "nchans".
564  ! If they do, is_nchans_var is set to TRUE.
565  do i = 1, tmp_var_ndims
566  call nclayer_check(nf90_inquire_dimension(ncid, tmp_var_dimids(i), tmp_var_dim_names(i), &
567  tmp_var_dim_sizes(i)))
568 
569  if (tmp_var_dim_names(i) == "nchans") is_nchans_var = .true.
570  end do
571 
572  if (is_nchans_var) then
573  ! Expand variable metadata first!
574  ! Make sure we have enough variable metadata storage
575  ! (and if not, reallocate!)
577 
578  ! Add to the total!
579  diag_chaninfo_store%total = diag_chaninfo_store%total + 1
580 
581  ! Store name and type!
582  diag_chaninfo_store%names(diag_chaninfo_store%total) = trim(tmp_var_name)
583 
584  ! Reset relative index to zero...
585  rel_index = 0
586 
587  ! For the rest of the code, we basically do the following:
588  ! -> We allocate a temporary data storage variable.
589  ! -> We set the NLAYER variable type for the variable.
590  ! -> We fetch all of the data for the variable.
591  ! -> We search, starting at the end of the variable, for
592  ! fill bytes. We keep going if we see filler bytes, and
593  ! stop when we encounter a non-fill byte.
594  ! -> Since the place we stop is where we last stored a value,
595  ! we set our relative index to the stopped index variable.
596  ! -> We deallocate our temporary data storage variable.
597  ! -> We set our type_index to update our data storage array count.
598 
599  if (tmp_var_type == nf90_byte) then
601  call nc_diag_chaninfo_resize_byte(int8(diag_chaninfo_store%nchans), .false.)
602  allocate(byte_buffer(diag_chaninfo_store%nchans))
603  call nclayer_check(nf90_get_var(ncid, var_index, byte_buffer))
604 
605  do j = diag_chaninfo_store%nchans, 1, -1
606  if (byte_buffer(j) /= nlayer_fill_byte) then
607  exit
608  end if
609  end do
610 
611  rel_index = j
612 
613  deallocate(byte_buffer)
614 
615  type_index = 1
616  else if (tmp_var_type == nf90_short) then
618  call nc_diag_chaninfo_resize_short(int8(diag_chaninfo_store%nchans), .false.)
619  allocate(short_buffer(diag_chaninfo_store%nchans))
620  call nclayer_check(nf90_get_var(ncid, var_index, short_buffer))
621 
622  do j = diag_chaninfo_store%nchans, 1, -1
623  if (short_buffer(j) /= nlayer_fill_short) then
624  exit
625  end if
626  end do
627 
628  rel_index = j
629 
630  deallocate(short_buffer)
631 
632  type_index = 2
633  else if (tmp_var_type == nf90_int) then
635  call nc_diag_chaninfo_resize_long(int8(diag_chaninfo_store%nchans), .false.)
636  allocate(long_buffer(diag_chaninfo_store%nchans))
637  call nclayer_check(nf90_get_var(ncid, var_index, long_buffer))
638 
639  do j = diag_chaninfo_store%nchans, 1, -1
640  if (long_buffer(j) /= nlayer_fill_long) then
641  exit
642  end if
643  end do
644 
645  rel_index = j
646 
647  deallocate(long_buffer)
648 
649  type_index = 3
650  else if (tmp_var_type == nf90_float) then
652  call nc_diag_chaninfo_resize_rsingle(int8(diag_chaninfo_store%nchans), .false.)
653  allocate(rsingle_buffer(diag_chaninfo_store%nchans))
654  call nclayer_check(nf90_get_var(ncid, var_index, rsingle_buffer))
655 
656  do j = diag_chaninfo_store%nchans, 1, -1
657  if (rsingle_buffer(j) /= nlayer_fill_float) then
658  exit
659  end if
660  end do
661 
662  rel_index = j
663 
664  deallocate(rsingle_buffer)
665 
666  type_index = 4
667  else if (tmp_var_type == nf90_double) then
669  call nc_diag_chaninfo_resize_rdouble(int8(diag_chaninfo_store%nchans), .false.)
670  allocate(rdouble_buffer(diag_chaninfo_store%nchans))
671  call nclayer_check(nf90_get_var(ncid, var_index, rdouble_buffer))
672 
673  do j = diag_chaninfo_store%nchans, 1, -1
674  if (rdouble_buffer(j) /= nlayer_fill_double) then
675  exit
676  end if
677  end do
678 
679  rel_index = j
680 
681  deallocate(rdouble_buffer)
682 
683  type_index = 5
684  else if (tmp_var_type == nf90_char) then
686  call nc_diag_chaninfo_resize_string(int8(diag_chaninfo_store%nchans), .false.)
687  allocate(string_buffer(diag_chaninfo_store%nchans, tmp_var_dim_sizes(1)))
688  call nclayer_check(nf90_get_var(ncid, var_index, string_buffer))
689 
690  do j = diag_chaninfo_store%nchans, 1, -1
691  if (string_buffer(j, 1) /= nlayer_fill_char) then
692  exit
693  end if
694  end do
695 
696  rel_index = j
697 
698  deallocate(string_buffer)
699 
700  ! Set max string length constraint
701  diag_chaninfo_store%max_str_lens(diag_chaninfo_store%total) = tmp_var_dim_sizes(1)
702 
703  type_index = 6
704  else
705  ! The type is not supported by chaninfo - error!
706  call nclayer_error("NetCDF4 type invalid!")
707  end if
708 
709  print *, trim(tmp_var_name), "rel index", rel_index
710 
711  ! Now add a relative position... based on the next position!
712 
713  ! First, increment the number of variables stored for this type:
714  diag_chaninfo_store%acount_v(type_index) = diag_chaninfo_store%acount_v(type_index) + 1
715 
716  ! Then, set the next variable's relative positioning,
717  ! based on the number of variables stored for this type.
718  diag_chaninfo_store%var_rel_pos(diag_chaninfo_store%total) = diag_chaninfo_store%acount_v(type_index)
719 
720  ! Initialize the amount of memory used to 0.
721  diag_chaninfo_store%var_usage(diag_chaninfo_store%total) = 0
722 
723  ! Set relative index!
724  diag_chaninfo_store%rel_indexes(diag_chaninfo_store%total) = rel_index
725 
726  ! Set variable ID! Note that var_index here is the actual variable ID.
727  diag_chaninfo_store%var_ids(diag_chaninfo_store%total) = var_index
728  end if
729 
730  !call nc_diag_cat_metadata_add_var(trim(tmp_var_name), tmp_var_type, tmp_var_ndims, tmp_var_dim_names)
731  end if
732 
733  ! Deallocate
734  deallocate(tmp_var_dimids)
735  deallocate(tmp_var_dim_names)
736  deallocate(tmp_var_dim_sizes)
737  end do
738 
739  ! Set our definition lock!
740  diag_chaninfo_store%def_lock = .true.
741  end subroutine nc_diag_chaninfo_load_def
742 
743  ! Write out chaninfo variable dimensions and variable
744  ! definitions to NetCDF via the NetCDF API.
745  !
746  ! Commit the current variables and make them known to NetCDF to
747  ! allow chaninfo variable data writing. If successfully written,
748  ! this will always set the definition lock flag to prevent any
749  ! further changes.
750  !
751  ! Definitions are only written once for every file opened, and
752  ! can not be modified or written again within the opened file.
753  ! This is enforced with a definition lock (def_lock) that is
754  ! set here and checked everywhere.
755  !
756  ! If definitions are already locked, no additional definitions
757  ! will be created. Depending on conditions, the following may
758  ! occur:
759  !
760  ! -> If the internal argument is defined and set to TRUE, no
761  ! error will be triggered. This is used internally by
762  ! nc_diag_write to prevent errors from occuring when the
763  ! lock may have already been set elsewhere.
764  !
765  ! -> Otherwise, an error will be triggered, since the
766  ! definition write occurred when the definitions were
767  ! already written and locked.
768  !
769  ! The inner workings:
770  !
771  ! -> First and foremost, it performs sanity checks to ensure
772  ! that we have a file loaded. If the check fails, an error
773  ! occurs.
774  !
775  ! -> It then checks to make sure we have chaninfo variables to
776  ! write in the first place. If we don't have any, we simply
777  ! return.
778  !
779  ! -> We then do another sanity check to ensure that nchans is
780  ! defined. We probably shouldn't have any variables in the
781  ! first place if nchans isn't defined, but it doesn't hurt
782  ! to check! (If this check fails, we probably have a
783  ! serious bug...)
784  !
785  ! -> If necessary (aka not in append mode, where this might
786  ! already exist), define the nchans dimension in NetCDF.
787  !
788  ! -> For every variable, fetch the type and name of the
789  ! variable. If the variable is a string type, we also
790  ! figure out the maximum string length, and create an
791  ! extra dimension for that as well. Finally, we can go and
792  ! define the variable itself to NetCDF, with the variable's
793  ! respective dimensions (and NetCDF dimension IDs).
794  !
795  ! -> We then add the variable to the varattr list to allow
796  ! variable attributes for the chaninfo variable.
797  !
798  ! -> If we're not in append mode, we set the appropriate
799  ! chunking and compression settings for the variable to
800  ! make storing the data more efficient.
801  !
802  ! -> After we've gone through all of the chaninfo variables,
803  ! we lock the definitions. That's it!
804  !
805  ! This is an internal subroutine, and is NOT meant to be called
806  ! outside of nc_diag_write. Calling this subroutine in your
807  ! program may result in unexpected behavior and/or data
808  ! corruption!
809  !
810  ! Args:
811  ! internal (logical, optional): whether or not to disable
812  ! triggering an error when a definition lock is
813  ! detected. This flag is used internally for the final
814  ! nc_diag_write, where this flag is purposely set to
815  ! avoid any errors with definition locking, since the
816  ! lock could have already been set earlier by
817  ! nc_diag_lock_def or others.
818  !
819  ! Raises:
820  ! If definitions are already locked, and the internal
821  ! argument is not set or is not TRUE, this will result in an
822  ! error.
823  !
824  ! If the nchans dimension hasn't been defined yet, this will
825  ! result in an error.
826  !
827  ! If there is no file open (or the file is already closed),
828  ! this will result in an error.
829  !
830  ! Other errors may result from invalid data storage, NetCDF
831  ! errors, or even a bug. See the called subroutines'
832  ! documentation for details.
833  !
834  subroutine nc_diag_chaninfo_write_def(internal)
835  logical, intent(in), optional :: internal
836 
837  ! Just write the definitions out!
838  integer(i_llong) :: curdatindex
839  integer(i_byte) :: data_type
840  integer(i_long) :: data_type_index
841  character(len=100) :: data_name
842  integer(i_long) :: nc_data_type
843 
844  integer(i_long) :: tmp_dim_id
845  character(len=120) :: data_dim_name
846 
847  character(len=:), allocatable :: string_arr(:)
848 
849 #ifdef ENABLE_ACTION_MSGS
850  character(len=1000) :: action_str
851 
852  if (nclayer_enable_action) then
853  if (present(internal)) then
854  write(action_str, "(A, L, A)") "nc_diag_chaninfo_write_def(internal = ", internal, ")"
855  else
856  write(action_str, "(A)") "nc_diag_chaninfo_write_def(internal = (not specified))"
857  end if
858  call nclayer_actionm(trim(action_str))
859  end if
860 #endif
861  ! Ensure that we have a file open and that things are loaded!
862  if (init_done .AND. allocated(diag_chaninfo_store)) then
863  ! Ensure that we have at least one variable to store!
864  ! Otherwise, just return and do nothing.
865  if (diag_chaninfo_store%total > 0) then
866  ! Make sure nchans is defined before doing anything!
867  if (diag_chaninfo_store%nchans /= -1) then
868  ! Finally, make sure definitions are not locked!
869  if (.NOT. diag_chaninfo_store%def_lock) then
870  ! First, set the dimensions... if necessary!
871  if (.NOT. append_only) &
872  call nclayer_check(nf90_def_dim(ncid, "nchans", diag_chaninfo_store%nchans, diag_chaninfo_store%nchans_dimid))
873 
874  ! Once we have the dimension, we can start writing
875  ! variable definitions!
876  do curdatindex = 1, diag_chaninfo_store%total
877  ! Fetch variable name and type:
878  data_name = diag_chaninfo_store%names(curdatindex)
879  data_type = diag_chaninfo_store%types(curdatindex)
880 
881  ! Figure out where our data is stored, given var_rel_pos
882  ! and nchans... (see equation/discussion above for more
883  ! details!)
884  data_type_index = 1 + &
885  ((diag_chaninfo_store%var_rel_pos(curdatindex) - 1) * diag_chaninfo_store%nchans)
886 
887  call nclayer_info("chaninfo: defining " // trim(data_name))
888 
889  ! Map our NLAYER type to the NF90 NetCDF native type!
890  if (data_type == nlayer_byte) nc_data_type = nf90_byte
891  if (data_type == nlayer_short) nc_data_type = nf90_short
892  if (data_type == nlayer_long) nc_data_type = nf90_int
893  if (data_type == nlayer_float) nc_data_type = nf90_float
894  if (data_type == nlayer_double) nc_data_type = nf90_double
895  if (data_type == nlayer_string) nc_data_type = nf90_char
896 
897 #ifdef _DEBUG_MEM_
898  print *, "chaninfo part 1"
899 #endif
900 
901  ! If our variable type is a string, we need to compute the maximum
902  ! string length.
903  !
904  ! If we're trimming, we take the maximum of the length of strings
905  ! in the variable, and use that as our maximum string length.
906  !
907  ! Otherwise, we simply use the previously defined fixed length,
908  ! which is already stored as the maximum string length from the
909  ! initial string add.
910  !
911  ! Once we know our maximum string length, we add that as a
912  ! dimension, and use it (along with our nchans dimension) to
913  ! create our string chaninfo variable!
914 
915  if (data_type == nlayer_string) then
916  ! Figure out the dimension name for this chaninfo variable
917  write (data_dim_name, "(A, A)") trim(data_name), "_maxstrlen"
918 
919  ! Assume that the maximum string length is 10000
920  ! Allocate an array of 10000, with a size of the
921  ! variable's var_usage
922  allocate(character(10000) :: string_arr(diag_chaninfo_store%var_usage(curdatindex)))
923 
924  ! Fetch the strings from our variable storage
925  string_arr = diag_chaninfo_store%ci_string(data_type_index:(data_type_index + &
926  diag_chaninfo_store%var_usage(curdatindex) - 1))
927 
928  ! If trimming is enabled, we haven't found our max_str_len yet.
929  ! Go find it!
930  if (enable_trim) then
931  ! Save the max string len
932  diag_chaninfo_store%max_str_lens(curdatindex) = &
933  max_len_string_array(string_arr, diag_chaninfo_store%var_usage(curdatindex))
934  end if
935 
936  ! Add our custom string dimension to NetCDF, if necessary
937  if (.NOT. append_only) &
938  call nclayer_check(nf90_def_dim(ncid, data_dim_name, &
939  diag_chaninfo_store%max_str_lens(curdatindex), &
940  tmp_dim_id))
941 #ifdef _DEBUG_MEM_
942  print *, "Defining char var type..."
943 #endif
944  ! Add our string variable to NetCDF!
945  if (.NOT. append_only) &
946  call nclayer_check(nf90_def_var(ncid, diag_chaninfo_store%names(curdatindex), &
947  nc_data_type, (/ tmp_dim_id, diag_chaninfo_store%nchans_dimid /), &
948  diag_chaninfo_store%var_ids(curdatindex)))
949 #ifdef _DEBUG_MEM_
950  print *, "Done defining char var type..."
951 #endif
952  ! Deallocate temp string array
953  deallocate(string_arr)
954  else
955  ! Nothing fancy here!
956  ! Just add our non-string variable to NetCDF!
957  if (.NOT. append_only) &
958  call nclayer_check(nf90_def_var(ncid, diag_chaninfo_store%names(curdatindex), &
959  nc_data_type, diag_chaninfo_store%nchans_dimid, &
960  diag_chaninfo_store%var_ids(curdatindex)))
961  end if
962 
963 #ifdef _DEBUG_MEM_
964  print *, "chaninfo part 2"
965 #endif
966 
967  ! Make our variable known to varattr - add it to the varattr database!
968  call nc_diag_varattr_add_var(diag_chaninfo_store%names(curdatindex), &
969  diag_chaninfo_store%types(curdatindex), &
970  diag_chaninfo_store%var_ids(curdatindex))
971 
972  ! If we are not appending, make sure to also set chunking and
973  ! compression for efficiency + optimization!
974  if (.NOT. append_only) then
975  ! If we're storing a string, we need to specify both dimensions
976  ! for our chunking parameters. Otherwise, we just need to
977  ! specify nchans...
978  if (data_type == nlayer_string) then
979  call nclayer_check(nf90_def_var_chunking(ncid, diag_chaninfo_store%var_ids(curdatindex), &
980  nf90_chunked, (/ diag_chaninfo_store%max_str_lens(curdatindex), diag_chaninfo_store%nchans /)))
981  else
982  call nclayer_check(nf90_def_var_chunking(ncid, diag_chaninfo_store%var_ids(curdatindex), &
983  nf90_chunked, (/ diag_chaninfo_store%nchans /)))
984  end if
985 
986  ! Enable zlib (gzip-like) compression based on our level settings
987  call nclayer_check(nf90_def_var_deflate(ncid, diag_chaninfo_store%var_ids(curdatindex), &
988  1, 1, nlayer_compression))
989  end if
990  end do
991 
992  ! Lock the definitions!
993  diag_chaninfo_store%def_lock = .true.
994  else
995  ! Show an error message if we didn't suppress errors on purpose
996  if(.NOT. present(internal)) &
997  call nclayer_error("Can't write definitions - definitions have already been written and locked!")
998  end if
999  else
1000  call nclayer_error("Can't write definitions - number of chans not set yet!")
1001  end if
1002 
1003  ! End: if (diag_chaninfo_store%total > 0)
1004  end if
1005  else
1006  call nclayer_error("Can't write definitions - NetCDF4 layer not initialized yet!")
1007  end if
1008  end subroutine nc_diag_chaninfo_write_def
1009 
1010  ! Write all of the currently stored chaninfo data to NetCDF via
1011  ! the NetCDF APIs ("put").
1012  !
1013  ! This will go through all of the variables stored in chaninfo,
1014  ! and write their data to NetCDF.
1015  !
1016  ! Buffer flushing mode is enabled if flush_data_only is set and
1017  ! is TRUE. Otherwise, this will operate normally.
1018  !
1019  ! For buffer flushing mode, data locking will not be performed.
1020  ! Instead, it "flushes" the variable storage buffer. For all
1021  ! of the variables stored, it increments the relative index of
1022  ! the variable with the amount of data currently stored in the
1023  ! variable.
1024  !
1025  ! (Essentially, new_rel_index = old_rel_index + var_data_count)
1026  !
1027  ! Recall that the relative index stores the position of the last
1028  ! data entered for the variable. This is set by write_data, as
1029  ! well as load_def for the data append mode. In turn, write_data
1030  ! also uses it to store at the correct position.
1031  !
1032  ! We also reset the var_usage, or the variable memory usage
1033  ! counter, back to zero to allow data storage to start at the
1034  ! beginning again. We use var_usage in write_data and in the
1035  ! storage subroutines to keep track of how much data we're
1036  ! storing, and how much we need to "read" from the array to
1037  ! store the data in NetCDF4 efficiently and within bounds.
1038  !
1039  ! A quick example:
1040  ! -> If we have 2 elements, var_usage (variable memory usage)
1041  ! is initially 2, and rel_index (variable relative index,
1042  ! or our starting position) is initially 0.
1043  !
1044  ! -> We flush the buffer. Since we flushed our buffer,
1045  ! var_usage is reset to 0, and rel_index is now 2 since
1046  ! we stored 2 elements.
1047  !
1048  ! -> If we add 3 elements, we get a var_usage of 3 (for 3
1049  ! elements stored), and rel_index stays the same (2).
1050  !
1051  ! -> When we finally flush or write, this time we know to
1052  ! start at element number 3 (rel_index), and we know to
1053  ! write 3 elements from there (var_usage).
1054  !
1055  ! -> We now have a total of 5 elements! Indicies 1-2 were
1056  ! stored with the flush, and indicies 3-5 were stored
1057  ! afterwards - all thanks to buffer flushing!
1058  !
1059  ! Finally, if data flushing mode is enabled, the data_lock is
1060  ! not set to allow additional data to be written in the future.
1061  !
1062  ! However, if data flushing mode is not set, or it is disabled,
1063  ! we assume that we are writing only one more time (or once,
1064  ! depending on if buffer flushing was ever enabled or not).
1065  ! Therefore, we set the data_lock (data writing lock) to TRUE
1066  ! in this case, assuming data writing was successful.
1067  !
1068  ! If data writing has already been locked, this will error.
1069  !
1070  ! If data flushing mode is disabled, we will also check to see
1071  ! if each variable's data fills up the nchans dimension.
1072  !
1073  ! Depending on the strictness (strict_check), if the data is
1074  ! not filled to the nchans dimension, it could either result in
1075  ! an error (if strict_check is TRUE), or a warning (if
1076  ! strict_check is FALSE).
1077  !
1078  ! This is an internal subroutine, and is NOT meant to be called
1079  ! outside of nc_diag_write. Calling this subroutine in your
1080  ! program may result in unexpected behavior and/or data
1081  ! corruption!
1082  !
1083  ! Args:
1084  ! flush_data_only (logical, optional): whether to only flush
1085  ! the chaninfo data buffers or not. If we flush data,
1086  ! data locking will not be set.
1087  !
1088  ! Raises:
1089  ! If data writing has already been locked, and the data
1090  ! flushing argument is not set or is not TRUE, this will
1091  ! result in an error.
1092  !
1093  ! If the nchans dimension hasn't been defined yet, this will
1094  ! result in an error.
1095  !
1096  ! If strict checking (strict_check) is enabled, and a
1097  ! variable's data doesn't fill to the nchans dimension,
1098  ! this will result in an error.
1099  !
1100  ! If there is no file open (or the file is already closed),
1101  ! this will result in an error.
1102  !
1103  ! Other errors may result from invalid data storage, NetCDF
1104  ! errors, or even a bug. See the called subroutines'
1105  ! documentation for details.
1106  !
1107  subroutine nc_diag_chaninfo_write_data(flush_data_only)
1108  ! Optional internal flag to only flush data - if this is
1109  ! true, data flushing will be performed, and the data will
1110  ! NOT be locked.
1111  logical, intent(in), optional :: flush_data_only
1112 
1113  integer(i_byte) :: data_type
1114  integer(i_long) :: data_type_index
1115  character(len=100) :: data_name
1116 
1117  character(len=1000) :: nchan_empty_msg
1118 
1119  integer(i_llong) :: curdatindex, j
1120  integer(i_long) :: string_arr_maxlen
1121 
1122  character(len=:), allocatable :: string_arr(:)
1123 
1124 #ifdef ENABLE_ACTION_MSGS
1125  character(len=1000) :: action_str
1126 
1127  if (nclayer_enable_action) then
1128  if (present(flush_data_only)) then
1129  write(action_str, "(A, L, A)") "nc_diag_chaninfo_write_data(flush_data_only = ", flush_data_only, ")"
1130  else
1131  write(action_str, "(A)") "nc_diag_chaninfo_write_data(flush_data_only = (not specified))"
1132  end if
1133  call nclayer_actionm(trim(action_str))
1134  end if
1135 #endif
1136  ! Check to make sure a file is open / things are loaded!
1137  if (init_done .AND. allocated(diag_chaninfo_store)) then
1138  ! Check to see if we have any variables to write in the
1139  ! first place!
1140  if (diag_chaninfo_store%total > 0) then
1141  ! Check to make sure that we have nchans defined!
1142  if (diag_chaninfo_store%nchans /= -1) then
1143  ! Check if we can still write any data!
1144  if (.NOT. diag_chaninfo_store%data_lock) then
1145  ! Iterate through all of our variables!
1146  do curdatindex = 1, diag_chaninfo_store%total
1147  ! Fetch the variable's name and type!
1148  data_name = diag_chaninfo_store%names(curdatindex)
1149  data_type = diag_chaninfo_store%types(curdatindex)
1150 
1151  ! Figure out where our data is stored, given var_rel_pos
1152  ! and nchans... (see equation/discussion above for more
1153  ! details!)
1154  data_type_index = 1 + &
1155  ((diag_chaninfo_store%var_rel_pos(curdatindex) - 1) * diag_chaninfo_store%nchans)
1156 
1157  call nclayer_info("chaninfo: writing " // trim(data_name))
1158 
1159  ! Warn about low data filling... but only if we are finishing
1160  ! our data write (or writing once) - basically, we're NOT in
1161  ! flushing data mode!
1162  if ((.NOT. (present(flush_data_only) .AND. flush_data_only)) .AND. &
1163  ((diag_chaninfo_store%var_usage(curdatindex) + &
1164  diag_chaninfo_store%rel_indexes(curdatindex)) < diag_chaninfo_store%nchans)) then
1165  ! NOTE - I0 and TRIM are Fortran 95 specs
1166  write (nchan_empty_msg, "(A, A, A, I0, A, I0, A)") "Amount of data written in ", &
1167  trim(data_name), " (", &
1168  diag_chaninfo_store%var_usage(curdatindex) + &
1169  diag_chaninfo_store%rel_indexes(curdatindex), &
1170  ")" // char(10) // &
1171  " is less than nchans (", diag_chaninfo_store%nchans, ")!"
1172 
1173  ! If we are set to strict checking mode, error.
1174  ! Otherwise, just show a warning.
1175  if (diag_chaninfo_store%strict_check) then
1176  call nclayer_error(trim(nchan_empty_msg))
1177  else
1178  call nclayer_warning(trim(nchan_empty_msg))
1179  end if
1180  end if
1181 
1182 #ifdef _DEBUG_MEM_
1183  print *, "****** Processing ******"
1184  print *, "data_name:"
1185  print *, data_name
1186  print *, "data_type:"
1187  print *, data_type
1188  print *, "data_type_index:"
1189  print *, data_type_index
1190  print *, "diag_chaninfo_store%var_ids(curdatindex):"
1191  print *, diag_chaninfo_store%var_ids(curdatindex)
1192  print *, "diag_chaninfo_store%var_usage(curdatindex):"
1193  print *, diag_chaninfo_store%var_usage(curdatindex)
1194  print *, "Upper range (data_type_index + &"
1195  print *, " diag_chaninfo_store%var_usage(curdatindex) - 1):"
1196  print *, (data_type_index + &
1197  diag_chaninfo_store%var_usage(curdatindex) - 1)
1198 #endif
1199  ! Make sure we have variable data to write in the first place!
1200  !
1201  ! If we do, we essentially:
1202  !
1203  ! -> Find the right type to save to.
1204  !
1205  ! -> If we are NOT storing a string, we just store a subsection
1206  ! of our variable storage array at (1 + rel_index) in the
1207  ! NetCDF variable.
1208  !
1209  ! -> If we are storing a string, we create our own array to
1210  ! store all of our strings in to standardize the length
1211  ! (e.g. a 3, 4, and 5 character string is expanded to
1212  ! a 5, 5, and 5 character string array). This is needed
1213  ! to store all strings at once and match the NetCDF bounds.
1214  ! Once done, the array is sent through the NetCDF API for
1215  ! data storage. We deallocate the array once we're done!
1216  !
1217  if (diag_chaninfo_store%var_usage(curdatindex) > 0) then
1218  if (data_type == nlayer_byte) then
1219 #ifdef _DEBUG_MEM_
1220  print *, "Resulting data to be stored:"
1221  print *, diag_chaninfo_store%ci_byte(data_type_index:(data_type_index + &
1222  diag_chaninfo_store%var_usage(curdatindex) - 1))
1223 #endif
1224  call nclayer_check(nf90_put_var(ncid, diag_chaninfo_store%var_ids(curdatindex), &
1225  diag_chaninfo_store%ci_byte(data_type_index:(data_type_index + &
1226  diag_chaninfo_store%var_usage(curdatindex) - 1)), &
1227  start = (/ 1 + diag_chaninfo_store%rel_indexes(curdatindex) /) &
1228  ))
1229  else if (data_type == nlayer_short) then
1230  call nclayer_check(nf90_put_var(ncid, diag_chaninfo_store%var_ids(curdatindex), &
1231  diag_chaninfo_store%ci_short(data_type_index:(data_type_index + &
1232  diag_chaninfo_store%var_usage(curdatindex) - 1)), &
1233  start = (/ 1 + diag_chaninfo_store%rel_indexes(curdatindex) /) &
1234  ))
1235  else if (data_type == nlayer_long) then
1236 #ifdef _DEBUG_MEM_
1237  print *, "Resulting data to be stored:"
1238  print *, diag_chaninfo_store%ci_long(data_type_index:(data_type_index + &
1239  diag_chaninfo_store%var_usage(curdatindex) - 1))
1240  print *, "start index:"
1241  print *, 1 + diag_chaninfo_store%rel_indexes(curdatindex)
1242 #endif
1243  call nclayer_check(nf90_put_var(ncid, diag_chaninfo_store%var_ids(curdatindex), &
1244  diag_chaninfo_store%ci_long(data_type_index:(data_type_index + &
1245  diag_chaninfo_store%var_usage(curdatindex) - 1)), &
1246  start = (/ 1 + diag_chaninfo_store%rel_indexes(curdatindex) /) &
1247  ))
1248  else if (data_type == nlayer_float) then
1249  call nclayer_check(nf90_put_var(ncid, diag_chaninfo_store%var_ids(curdatindex), &
1250  diag_chaninfo_store%ci_rsingle(data_type_index:(data_type_index + &
1251  diag_chaninfo_store%var_usage(curdatindex) - 1)), &
1252  start = (/ 1 + diag_chaninfo_store%rel_indexes(curdatindex) /) &
1253  ))
1254  else if (data_type == nlayer_double) then
1255  call nclayer_check(nf90_put_var(ncid, diag_chaninfo_store%var_ids(curdatindex), &
1256  diag_chaninfo_store%ci_rdouble(data_type_index:(data_type_index + &
1257  diag_chaninfo_store%var_usage(curdatindex) - 1)), &
1258  start = (/ 1 + diag_chaninfo_store%rel_indexes(curdatindex) /) &
1259  ))
1260  else if (data_type == nlayer_string) then
1261  ! Storing to another variable may seem silly, but it's necessary
1262  ! to avoid "undefined variable" errors, thanks to the compiler's
1263  ! super optimization insanity...
1264  string_arr_maxlen = diag_chaninfo_store%max_str_lens(curdatindex)
1265  allocate(character(string_arr_maxlen) :: &
1266  string_arr(diag_chaninfo_store%var_usage(curdatindex)))
1267  if (enable_trim) then
1268  do j = data_type_index, data_type_index + &
1269  diag_chaninfo_store%var_usage(curdatindex) - 1
1270  string_arr(j - data_type_index + 1) = &
1271  trim(diag_chaninfo_store%ci_string(j))
1272  end do
1273 
1274 #ifdef _DEBUG_MEM_
1275  do j = 1, diag_chaninfo_store%var_usage(curdatindex)
1276  write (*, "(A, A, A)") "String: '", string_arr(j), "'"
1277  end do
1278 
1279  write (*, "(A, I0)") "string_arr_maxlen = ", string_arr_maxlen
1280  write (*, "(A, I0)") "diag_chaninfo_store%var_usage(curdatindex) = ", diag_chaninfo_store%var_usage(curdatindex)
1281 #endif
1282  else
1283  do j = data_type_index, data_type_index + &
1284  diag_chaninfo_store%var_usage(curdatindex) - 1
1285  string_arr(j - data_type_index + 1) = &
1286  diag_chaninfo_store%ci_string(j)
1287  end do
1288  end if
1289 
1290  call nclayer_check(nf90_put_var(ncid, diag_chaninfo_store%var_ids(curdatindex), &
1291  string_arr, &
1292  start = (/ 1, 1 + diag_chaninfo_store%rel_indexes(curdatindex) /), &
1293  count = (/ string_arr_maxlen, &
1294  diag_chaninfo_store%var_usage(curdatindex) /) ))
1295 
1296  deallocate(string_arr)
1297  else
1298  call nclayer_error("Critical error - unknown variable type!")
1299  end if
1300 
1301  ! Check for data flushing, and if so, update the relative indexes
1302  ! and set var_usage to 0.
1303  if (present(flush_data_only) .AND. flush_data_only) then
1304  diag_chaninfo_store%rel_indexes(curdatindex) = &
1305  diag_chaninfo_store%rel_indexes(curdatindex) + &
1306  diag_chaninfo_store%var_usage(curdatindex)
1307  diag_chaninfo_store%var_usage(curdatindex) = 0
1308 
1309 #ifdef _DEBUG_MEM_
1310  print *, "diag_chaninfo_store%rel_indexes(curdatindex) is now:"
1311  print *, diag_chaninfo_store%rel_indexes(curdatindex)
1312 #endif
1313  end if
1314  end if
1315  end do
1316 
1317  ! If we're flushing data, don't do anything...
1318  if (present(flush_data_only) .AND. flush_data_only) then
1319 #ifdef _DEBUG_MEM_
1320  print *, "In buffer flush mode!"
1321 #endif
1322  else
1323  ! Otherwise, lock data writing! Note that we do this,
1324  ! even if we have no data!
1325  diag_chaninfo_store%data_lock = .true.
1326 #ifdef _DEBUG_MEM_
1327  print *, "In data lock mode!"
1328 #endif
1329  end if
1330  else
1331  call nclayer_error("Can't write data - data have already been written and locked!")
1332  end if
1333  else
1334  call nclayer_error("Can't write data - number of chans not set yet!")
1335  end if
1336  end if
1337  else
1338  call nclayer_error("Can't write data - NetCDF4 layer not initialized yet!")
1339  end if
1340 
1341  end subroutine nc_diag_chaninfo_write_data
1342 
1343  ! Set the strict mode for chaninfo variables.
1344  !
1345  ! This sets the mode that determines how strict chaninfo's
1346  ! variable consistency checks will be.
1347  !
1348  ! During the final data write (nc_diag_chaninfo_write_data,
1349  ! without the buffering flag), chaninfo will check to see if all
1350  ! of the variables are filled, e.g. all variables have been
1351  ! stored up to nchans dimension.
1352  !
1353  ! If there are any variables that are not completely filled to
1354  ! the nchans dimension, one of the following may occur:
1355  !
1356  ! -> If strict mode is enabled, a consistency check error will
1357  ! occur and the program will exit.
1358  !
1359  ! -> If strict mode is disabled, this will only result in a
1360  ! consistency check warning. After the warning is
1361  ! displayed, normal operation will occur, including data
1362  ! writing. For values that are not in the variable (up to
1363  ! the nchans dimension), missing values will be placed.
1364  !
1365  ! By default, strict mode is disabled.
1366  !
1367  ! Since the strict mode is bound to the chaninfo type, it can
1368  ! only be set when a file is open and when diag_chaninfo_store
1369  ! is initialized. (It should be initialized if a file is open!)
1370  !
1371  ! If there isn't a file open / diag_chaninfo_store isn't
1372  ! initialized, an error will occur.
1373  !
1374  ! This is an internal subroutine, and is NOT meant to be called
1375  ! outside of nc_diag_write. Calling this subroutine in your
1376  ! program may result in unexpected behavior and/or data
1377  ! corruption!
1378  !
1379  ! Args:
1380  ! enable_strict (logical): boolean indicating whether to
1381  ! enable strict mode or not. If set to TRUE, strict mode
1382  ! will be enabled. Otherwise, it will be disabled.
1383  !
1384  ! Raises:
1385  ! If there is no file open (or the file is already closed),
1386  ! this will result in an error.
1387  !
1388  ! Although unlikely, other errors may indirectly occur.
1389  ! They may be general storage errors, or even a bug.
1390  ! See the called subroutines' documentation for details.
1391  !
1392  subroutine nc_diag_chaninfo_set_strict(enable_strict)
1393  logical, intent(in) :: enable_strict
1394 
1395  if (init_done .AND. allocated(diag_chaninfo_store)) then
1396  diag_chaninfo_store%strict_check = enable_strict
1397  else
1398  call nclayer_error("Can't set strictness level for chaninfo - NetCDF4 layer not initialized yet!")
1399  end if
1400  end subroutine nc_diag_chaninfo_set_strict
1401 
1402  ! Preallocate variable metadata storage (names, types, etc.).
1403  !
1404  ! This preallocates variable metadata storage for a given number
1405  ! of variables.
1406  !
1407  ! If properly defined, this can speed up chaninfo variable
1408  ! creation since reallocation will (hopefully) not be necessary
1409  ! for variable metadata storage, since it was preallocated here.
1410  !
1411  ! Variable metadata includes storing the variables' names,
1412  ! types, indicies, usage counts, etc. The metadata pre-allocated
1413  ! here is essentially the variable indexed arrays within our
1414  ! specific storage type!
1415  !
1416  ! Args:
1417  ! num_of_addl_vars (integer(i_llong)): the number of
1418  ! additional variables to preallocate metadata storage
1419  ! for.
1420  !
1421  ! Raises:
1422  ! If there is no file open (or the file is already closed),
1423  ! this will result in an error.
1424  !
1425  ! Other errors may result from invalid data storage, NetCDF
1426  ! errors, or even a bug. See the called subroutines'
1427  ! documentation for details.
1428  !
1429  subroutine nc_diag_chaninfo_prealloc_vars(num_of_addl_vars)
1430  integer(i_llong), intent(in) :: num_of_addl_vars
1431 #ifdef ENABLE_ACTION_MSGS
1432  character(len=1000) :: action_str
1433 
1434  if (nclayer_enable_action) then
1435  write(action_str, "(A, I0, A)") "nc_diag_chaninfo_prealloc_vars(num_of_addl_vars = ", num_of_addl_vars, ")"
1436  call nclayer_actionm(trim(action_str))
1437  end if
1438 #endif
1439  if (init_done .AND. allocated(diag_chaninfo_store)) then
1440  ! For all variable metadata fields:
1441  ! -> Check if the field is allocated.
1442  ! -> If not, allocate it with the default initial
1443  ! size, plus the number of additional variables
1444  ! specified in the argument.
1445  ! -> If it's allocated, check to see if the total
1446  ! number of variables exceeds our field's allocated
1447  ! size.
1448  ! -> If the size is exceeded, reallocate the field
1449  ! with the number of additional variables specified
1450  ! in the argument.
1451  !
1452  if (allocated(diag_chaninfo_store%names)) then
1453  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%names)) then
1454  call nc_diag_realloc(diag_chaninfo_store%names, num_of_addl_vars)
1455  end if
1456  else
1457  allocate(diag_chaninfo_store%names(nlayer_default_ent + num_of_addl_vars))
1458  end if
1459 
1460  if (allocated(diag_chaninfo_store%types)) then
1461  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%types)) then
1462  call nc_diag_realloc(diag_chaninfo_store%types, num_of_addl_vars)
1463  end if
1464  else
1465  allocate(diag_chaninfo_store%types(nlayer_default_ent + num_of_addl_vars))
1466  end if
1467 
1468  if (allocated(diag_chaninfo_store%var_rel_pos)) then
1469  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%var_rel_pos)) then
1470  call nc_diag_realloc(diag_chaninfo_store%var_rel_pos, num_of_addl_vars)
1471  end if
1472  else
1473  allocate(diag_chaninfo_store%var_rel_pos(nlayer_default_ent + num_of_addl_vars))
1474  diag_chaninfo_store%var_rel_pos = -1
1475  end if
1476 
1477  if (allocated(diag_chaninfo_store%var_usage)) then
1478  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%var_usage)) then
1479  call nc_diag_realloc(diag_chaninfo_store%var_usage, num_of_addl_vars)
1480  end if
1481  else
1482  allocate(diag_chaninfo_store%var_usage(nlayer_default_ent + num_of_addl_vars))
1483  diag_chaninfo_store%var_usage = 0
1484  end if
1485 
1486  if (allocated(diag_chaninfo_store%var_ids)) then
1487  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%var_ids)) then
1488  call nc_diag_realloc(diag_chaninfo_store%var_ids, num_of_addl_vars)
1489  end if
1490  else
1491  allocate(diag_chaninfo_store%var_ids(nlayer_default_ent + num_of_addl_vars))
1492  diag_chaninfo_store%var_ids = -1
1493  end if
1494 
1495  if (allocated(diag_chaninfo_store%max_str_lens)) then
1496  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%max_str_lens)) then
1497  call nc_diag_realloc(diag_chaninfo_store%max_str_lens, num_of_addl_vars)
1498  end if
1499  else
1500  allocate(diag_chaninfo_store%max_str_lens(nlayer_default_ent + num_of_addl_vars))
1501  diag_chaninfo_store%max_str_lens = -1
1502  end if
1503 
1504  if (allocated(diag_chaninfo_store%rel_indexes)) then
1505  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%rel_indexes)) then
1506  call nc_diag_realloc(diag_chaninfo_store%rel_indexes, num_of_addl_vars)
1507  end if
1508  else
1509  allocate(diag_chaninfo_store%rel_indexes(nlayer_default_ent + num_of_addl_vars))
1510  diag_chaninfo_store%rel_indexes = 0
1511  end if
1512  else
1513  call nclayer_error("NetCDF4 layer not initialized yet!")
1514  endif
1515  end subroutine nc_diag_chaninfo_prealloc_vars
1516 
1517  ! Preallocate actual variable data storage - the data itself.
1518  !
1519  ! This preallocates the variable data storage for a given
1520  ! variable type, and a given number of data elements or slots.
1521  !
1522  ! If properly defined, this can speed up chaninfo variable
1523  ! data insertion since reallocation will (hopefully) not be
1524  ! necessary for variable data storage, since it was preallocated
1525  ! here.
1526  !
1527  ! For example, if you have 10 float chaninfo variables, and
1528  ! nchans is 20, you can call:
1529  !
1530  ! nc_diag_chaninfo_prealloc_vars_storage(NLAYER_FLOAT, 200)
1531  !
1532  ! Assuming that no other float chaninfo variables get added,
1533  ! no reallocations should occur, therefore speeding up the
1534  ! variable data insertion process!
1535  !
1536  ! Note that this is a state-based subroutine call - by design,
1537  ! it preallocates the largest amount provided. For instance, if
1538  ! you attempted to preallocate 10 floats, then 9000 floats, then
1539  ! 5 floats, 20 floats will be preallocated.
1540  !
1541  ! Specifically, it looks like this:
1542  !
1543  ! -> Preallocate 10 floats - nothing allocated, so 10 floats
1544  ! allocated.
1545  !
1546  ! -> Preallocate 9000 floats - only 10 floats allocated, so
1547  ! reallocating to 9000.
1548  !
1549  ! -> Preallocate 20 floats - 9000 floats already allocated, so
1550  ! no need to do anything.
1551  !
1552  ! Args:
1553  ! nclayer_type (integer(i_byte)): the type of variable to
1554  ! preallocate data elements/slots for.
1555  ! num_of_addl_slots (integer(i_llong)): the number of
1556  ! additional variable data elements/slots to
1557  ! preallocate.
1558  !
1559  ! Raises:
1560  ! If the variable type is invalid, this will result in an
1561  ! error.
1562  !
1563  ! Other errors may result from invalid data storage, NetCDF
1564  ! errors, or even a bug. See the called subroutines'
1565  ! documentation for details.
1566  !
1567  subroutine nc_diag_chaninfo_prealloc_vars_storage(nclayer_type, num_of_addl_slots)
1568  integer(i_byte), intent(in) :: nclayer_type
1569  integer(i_llong), intent(in) :: num_of_addl_slots
1570 
1571 #ifdef ENABLE_ACTION_MSGS
1572  character(len=1000) :: action_str
1573 
1574  if (nclayer_enable_action) then
1575  write(action_str, "(A, I0, A, I0, A)") "nc_diag_chaninfo_prealloc_vars_storage(nclayer_type = ", nclayer_type, ", num_of_addl_slots = ", num_of_addl_slots, ")"
1576  call nclayer_actionm(trim(action_str))
1577  end if
1578 #endif
1579  ! Find the type specified, and attempt to pre-allocate.
1580  ! Note that FALSE is specified as an argument to ensure that
1581  ! the actual variable data storage usage count isn't
1582  ! incremented, since we're just preallocating here.
1583  !
1584  if (nclayer_type == nlayer_byte) then
1585  call nc_diag_chaninfo_resize_byte(num_of_addl_slots, .false.)
1586  else if (nclayer_type == nlayer_short) then
1587  call nc_diag_chaninfo_resize_short(num_of_addl_slots, .false.)
1588  else if (nclayer_type == nlayer_long) then
1589  call nc_diag_chaninfo_resize_long(num_of_addl_slots, .false.)
1590  else if (nclayer_type == nlayer_float) then
1591  call nc_diag_chaninfo_resize_rsingle(num_of_addl_slots, .false.)
1592  else if (nclayer_type == nlayer_double) then
1593  call nc_diag_chaninfo_resize_rdouble(num_of_addl_slots, .false.)
1594  else if (nclayer_type == nlayer_string) then
1595  call nc_diag_chaninfo_resize_string(num_of_addl_slots, .false.)
1596  else
1597  call nclayer_error("Invalid type specified for variable storage preallocation!")
1598  end if
1600 
1601  ! Expand variable metadata storage (names, types, etc.) for one
1602  ! single variable.
1603  !
1604  ! This ensures that there is enough variable metadata storage to
1605  ! add a single variable. If there isn't enough storage, it will
1606  ! reallocate as necessary. See this module's header for more
1607  ! information about how memory allocation works for variable
1608  ! metadata storage.
1609  !
1610  ! This is an internal subroutine, and is NOT meant to be called
1611  ! outside of nc_diag_write. Calling this subroutine in your
1612  ! program may result in unexpected behavior and/or data
1613  ! corruption!
1614  !
1615  ! Args:
1616  ! num_of_addl_vars (integer(i_llong)): the number of
1617  ! additional variables to preallocate metadata storage
1618  ! for.
1619  !
1620  ! Raises:
1621  ! If nchans has not been set yet, this will result in an
1622  ! error.
1623  !
1624  ! If there is no file open (or the file is already closed),
1625  ! this will result in an error.
1626  !
1627  ! Other errors may result from invalid data storage, NetCDF
1628  ! errors, or even a bug. See the called subroutines'
1629  ! documentation for details.
1630  !
1631  subroutine nc_diag_chaninfo_expand
1632  integer(i_llong) :: addl_fields
1633  ! Did we realloc at all?
1634  logical :: meta_realloc
1635  meta_realloc = .false.
1636 
1637  if (init_done .AND. allocated(diag_chaninfo_store)) then
1638  addl_fields = 1 + (nlayer_default_ent * (nlayer_multi_base ** diag_chaninfo_store%alloc_multi))
1639  if (diag_chaninfo_store%nchans /= -1) then
1640 
1641  ! For all variable metadata fields:
1642  ! -> Check if the field is allocated.
1643  ! -> If not, allocate it with the default initial
1644  ! size, and initialize it with blank values!
1645  ! -> If it's allocated, check to see if the total
1646  ! number of variables exceeds our field's
1647  ! allocated size.
1648  ! -> If the size is exceeded, reallocate the
1649  ! field, and indicate that a reallocation has
1650  ! occurred.
1651  !
1652  if (allocated(diag_chaninfo_store%names)) then
1653  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%names)) then
1654  call nc_diag_realloc(diag_chaninfo_store%names, addl_fields)
1655  meta_realloc = .true.
1656  end if
1657  else
1658  allocate(diag_chaninfo_store%names(nlayer_default_ent))
1659  end if
1660 
1661  if (allocated(diag_chaninfo_store%types)) then
1662  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%types)) then
1663  call nc_diag_realloc(diag_chaninfo_store%types, addl_fields)
1664  meta_realloc = .true.
1665  end if
1666  else
1667  allocate(diag_chaninfo_store%types(nlayer_default_ent))
1668  end if
1669 
1670  if (allocated(diag_chaninfo_store%var_rel_pos)) then
1671  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%var_rel_pos)) then
1672  call nc_diag_realloc(diag_chaninfo_store%var_rel_pos, addl_fields)
1673  meta_realloc = .true.
1674  end if
1675  else
1676  allocate(diag_chaninfo_store%var_rel_pos(nlayer_default_ent))
1677  diag_chaninfo_store%var_rel_pos = -1
1678  end if
1679 
1680  if (allocated(diag_chaninfo_store%var_usage)) then
1681  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%var_usage)) then
1682  call nc_diag_realloc(diag_chaninfo_store%var_usage, addl_fields)
1683  meta_realloc = .true.
1684  end if
1685  else
1686  allocate(diag_chaninfo_store%var_usage(nlayer_default_ent))
1687  diag_chaninfo_store%var_usage = 0
1688  end if
1689 
1690  if (allocated(diag_chaninfo_store%var_ids)) then
1691  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%var_ids)) then
1692  call nc_diag_realloc(diag_chaninfo_store%var_ids, addl_fields)
1693  meta_realloc = .true.
1694  end if
1695  else
1696  allocate(diag_chaninfo_store%var_ids(nlayer_default_ent))
1697  diag_chaninfo_store%var_ids = -1
1698  end if
1699 
1700  if (allocated(diag_chaninfo_store%max_str_lens)) then
1701  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%max_str_lens)) then
1702  call nc_diag_realloc(diag_chaninfo_store%max_str_lens, addl_fields)
1703  meta_realloc = .true.
1704  end if
1705  else
1706  allocate(diag_chaninfo_store%max_str_lens(nlayer_default_ent))
1707  diag_chaninfo_store%max_str_lens = -1
1708  end if
1709 
1710  if (allocated(diag_chaninfo_store%rel_indexes)) then
1711  if (diag_chaninfo_store%total >= size(diag_chaninfo_store%rel_indexes)) then
1712  call nc_diag_realloc(diag_chaninfo_store%rel_indexes, addl_fields)
1713  meta_realloc = .true.
1714  end if
1715  else
1716  allocate(diag_chaninfo_store%rel_indexes(nlayer_default_ent))
1717  diag_chaninfo_store%rel_indexes = 0
1718  end if
1719 
1720  ! If reallocation occurred, increment our multiplier
1721  ! to allocate more and speed things up in the
1722  ! future!
1723  if (meta_realloc) then
1724  diag_chaninfo_store%alloc_multi = diag_chaninfo_store%alloc_multi + 1
1725  end if
1726  else
1727  call nclayer_error("Number of chans not set yet!")
1728  end if
1729  else
1730  call nclayer_error("NetCDF4 layer not initialized yet!")
1731  end if
1732  end subroutine nc_diag_chaninfo_expand
1733 
1734  ! Add a single scalar byte integer to the given chaninfo
1735  ! variable.
1736  !
1737  ! This adds a single value to the specified chaninfo variable.
1738  !
1739  ! If the variable does not already exist, it will be created,
1740  ! and the value will be inserted as the variable's first
1741  ! element.
1742  !
1743  ! Otherwise, the value will be inserted to the next empty spot.
1744  !
1745  ! Values are inserted in the order of the calls made. As such,
1746  ! this subroutine is best designed to be used in a loop, where
1747  ! for every channel iteration, a value is added using this
1748  ! subroutine.
1749  !
1750  ! This is an internal subroutine, and is NOT meant to be called
1751  ! outside of nc_diag_write. Calling this subroutine in your
1752  ! program may result in unexpected behavior and/or data
1753  ! corruption! (You should use the generic nc_diag_chaninfo
1754  ! instead!)
1755  !
1756  ! Args:
1757  ! chaninfo_name (character(len=*)): the chaninfo variable
1758  ! to store to.
1759  ! chaninfo_value (integer(i_byte)): the value to store.
1760  !
1761  ! Raises:
1762  ! If the data has already been locked, this will result in
1763  ! an error.
1764  !
1765  ! If definitions have already been locked, and a new
1766  ! variable is being created, this will result in an error.
1767  !
1768  ! If the variable is already full (e.g. it has nchans number
1769  ! of elements), this will result in an error.
1770  !
1771  ! The following errors will trigger indirectly from other
1772  ! subroutines called here:
1773  !
1774  ! If nchans has not been set yet, this will result in an
1775  ! error.
1776  !
1777  ! If there is no file open (or the file is already closed),
1778  ! this will result in an error.
1779  !
1780  ! Other errors may result from invalid data storage, NetCDF
1781  ! errors, or even a bug. See the called subroutines'
1782  ! documentation for details.
1783  !
1784  subroutine nc_diag_chaninfo_byte(chaninfo_name, chaninfo_value)
1785  character(len=*), intent(in) :: chaninfo_name
1786  integer(i_byte), intent(in) :: chaninfo_value
1787 
1788  integer(i_long) :: i, var_index, var_rel_index, type_index
1789 
1790 #ifdef ENABLE_ACTION_MSGS
1791  character(len=1000) :: action_str
1792 
1793  if (nclayer_enable_action) then
1794  write(action_str, "(A, I0, A)") "nc_diag_chaninfo_byte(chaninfo_name = " // chaninfo_name // ", chaninfo_value = ", chaninfo_value, ")"
1795  call nclayer_actionm(trim(action_str))
1796  end if
1797 #endif
1798  ! Make sure that data hasn't been locked
1799  if (diag_chaninfo_store%data_lock) then
1800  call nclayer_error("Can't add new data - data have already been written and locked!")
1801  end if
1802 
1803  ! For byte, type index is 1
1804  type_index = 1
1805 
1806  ! Default to -1
1807  var_index = -1
1808 
1809  ! Attempt to match the variable name + fetch the variable
1810  ! index first!
1811  do i = 1, diag_chaninfo_store%total
1812  if (diag_chaninfo_store%names(i) == chaninfo_name) then
1813  var_rel_index = diag_chaninfo_store%var_rel_pos(i)
1814  var_index = i
1815  exit
1816  end if
1817  end do
1818 
1819  if (var_index == -1) then
1820  ! Entry does not exist yet...
1821 
1822  ! First, check to make sure we can still define new variables.
1823  if (diag_chaninfo_store%def_lock) then
1824  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1825  end if
1826 
1827  ! Expand variable metadata first!
1828  ! Make sure we have enough variable metadata storage
1829  ! (and if not, reallocate!)
1831 
1832  ! Add to the total!
1833  diag_chaninfo_store%total = diag_chaninfo_store%total + 1
1834 
1835  ! Store name and type!
1836  diag_chaninfo_store%names(diag_chaninfo_store%total) = chaninfo_name
1838 
1839  ! We just need to add one entry...
1840  ! Call resize subroutine to ensure we have enough space
1841  ! (and if not, realloc!)
1843 
1844  ! Now add a relative position... based on the next position!
1845 
1846  ! First, increment the number of variables stored for this type:
1847  diag_chaninfo_store%acount_v(type_index) = diag_chaninfo_store%acount_v(type_index) + 1
1848 
1849  ! Then, set the next variable's relative positioning,
1850  ! based on the number of variables stored for this type.
1851  diag_chaninfo_store%var_rel_pos(diag_chaninfo_store%total) = diag_chaninfo_store%acount_v(type_index)
1852 
1853  ! Initialize the amount of memory used to 1.
1854  diag_chaninfo_store%var_usage(diag_chaninfo_store%total) = 1
1855 
1856  ! Set var_index to the total
1857  var_index = diag_chaninfo_store%total
1858  else
1859  ! Variable already exists!
1860 
1861  ! Check to make sure we can fit more data!
1862  ! (# data < nchans)
1863  if (diag_chaninfo_store%var_usage(var_index) + &
1864  diag_chaninfo_store%rel_indexes(var_index) >= diag_chaninfo_store%nchans) then
1865  call nclayer_error("Can't add new data - data added is exceeding nchan! Data must fit within nchan constraint.")
1866  endif
1867 
1868  ! Increment current variable count
1869  diag_chaninfo_store%var_usage(var_index) = &
1870  diag_chaninfo_store%var_usage(var_index) + 1
1871  end if
1872 
1873  ! Now add the actual entry!
1874  diag_chaninfo_store%ci_byte(1 + &
1875  ((diag_chaninfo_store%var_rel_pos(var_index) - 1) &
1876  * diag_chaninfo_store%nchans) &
1877  + (diag_chaninfo_store%var_usage(var_index) - 1)) = chaninfo_value
1878  end subroutine nc_diag_chaninfo_byte
1879 
1880  ! Add a single scalar short integer to the given chaninfo
1881  ! variable.
1882  !
1883  ! This adds a single value to the specified chaninfo variable.
1884  !
1885  ! If the variable does not already exist, it will be created,
1886  ! and the value will be inserted as the variable's first
1887  ! element.
1888  !
1889  ! Otherwise, the value will be inserted to the next empty spot.
1890  !
1891  ! Values are inserted in the order of the calls made. As such,
1892  ! this subroutine is best designed to be used in a loop, where
1893  ! for every channel iteration, a value is added using this
1894  ! subroutine.
1895  !
1896  ! This is an internal subroutine, and is NOT meant to be called
1897  ! outside of nc_diag_write. Calling this subroutine in your
1898  ! program may result in unexpected behavior and/or data
1899  ! corruption! (You should use the generic nc_diag_chaninfo
1900  ! instead!)
1901  !
1902  ! Args:
1903  ! chaninfo_name (character(len=*)): the chaninfo variable
1904  ! to store to.
1905  ! chaninfo_value (integer(i_short)): the value to store.
1906  !
1907  ! Raises:
1908  ! If the data has already been locked, this will result in
1909  ! an error.
1910  !
1911  ! If definitions have already been locked, and a new
1912  ! variable is being created, this will result in an error.
1913  !
1914  ! If the variable is already full (e.g. it has nchans number
1915  ! of elements), this will result in an error.
1916  !
1917  ! The following errors will trigger indirectly from other
1918  ! subroutines called here:
1919  !
1920  ! If nchans has not been set yet, this will result in an
1921  ! error.
1922  !
1923  ! If there is no file open (or the file is already closed),
1924  ! this will result in an error.
1925  !
1926  ! Other errors may result from invalid data storage, NetCDF
1927  ! errors, or even a bug. See the called subroutines'
1928  ! documentation for details.
1929  !
1930  subroutine nc_diag_chaninfo_short(chaninfo_name, chaninfo_value)
1931  character(len=*), intent(in) :: chaninfo_name
1932  integer(i_short), intent(in) :: chaninfo_value
1933 
1934  integer(i_long) :: i, var_index, var_rel_index, type_index
1935 
1936 #ifdef ENABLE_ACTION_MSGS
1937  character(len=1000) :: action_str
1938 
1939  if (nclayer_enable_action) then
1940  write(action_str, "(A, I0, A)") "nc_diag_chaninfo_short(chaninfo_name = " // chaninfo_name // ", chaninfo_value = ", chaninfo_value, ")"
1941  call nclayer_actionm(trim(action_str))
1942  end if
1943 #endif
1944  ! Make sure that data hasn't been locked
1945  if (diag_chaninfo_store%data_lock) then
1946  call nclayer_error("Can't add new data - data have already been written and locked!")
1947  end if
1948 
1949  ! For short, type index is 2
1950  type_index = 2
1951 
1952  ! Default to -1
1953  var_index = -1
1954 
1955  ! Attempt to match the variable name + fetch the variable
1956  ! index first!
1957  do i = 1, diag_chaninfo_store%total
1958  if (diag_chaninfo_store%names(i) == chaninfo_name) then
1959  var_rel_index = diag_chaninfo_store%var_rel_pos(i)
1960  var_index = i
1961  exit
1962  end if
1963  end do
1964 
1965  if (var_index == -1) then
1966  ! Entry does not exist yet...
1967 
1968  ! First, check to make sure we can still define new variables.
1969  if (diag_chaninfo_store%def_lock) then
1970  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1971  end if
1972 
1973  ! Expand variable metadata first!
1974  ! Make sure we have enough variable metadata storage
1975  ! (and if not, reallocate!)
1977 
1978  ! Add to the total!
1979  diag_chaninfo_store%total = diag_chaninfo_store%total + 1
1980 
1981  ! Store name and type!
1982  diag_chaninfo_store%names(diag_chaninfo_store%total) = chaninfo_name
1984 
1985  ! We just need to add one entry...
1986  ! Call resize subroutine to ensure we have enough space
1987  ! (and if not, realloc!)
1989 
1990  ! Now add a relative position... based on the next position!
1991 
1992  ! First, increment the number of variables stored for this type:
1993  diag_chaninfo_store%acount_v(type_index) = diag_chaninfo_store%acount_v(type_index) + 1
1994 
1995  ! Then, set the next variable's relative positioning,
1996  ! based on the number of variables stored for this type.
1997  diag_chaninfo_store%var_rel_pos(diag_chaninfo_store%total) = diag_chaninfo_store%acount_v(type_index)
1998 
1999  ! Initialize the amount of memory used to 1.
2000  diag_chaninfo_store%var_usage(diag_chaninfo_store%total) = 1
2001 
2002  ! Set var_index to the total
2003  var_index = diag_chaninfo_store%total
2004  else
2005  ! Variable already exists!
2006 
2007  ! Check to make sure we can fit more data!
2008  ! (# data < nchans)
2009  if (diag_chaninfo_store%var_usage(var_index) + &
2010  diag_chaninfo_store%rel_indexes(var_index) >= diag_chaninfo_store%nchans) then
2011  call nclayer_error("Can't add new data - data added is exceeding nchan! Data must fit within nchan constraint.")
2012  endif
2013 
2014  ! Increment current variable count
2015  diag_chaninfo_store%var_usage(var_index) = &
2016  diag_chaninfo_store%var_usage(var_index) + 1
2017  end if
2018 
2019  ! Now add the actual entry!
2020  diag_chaninfo_store%ci_short(1 + &
2021  ((diag_chaninfo_store%var_rel_pos(var_index) - 1) &
2022  * diag_chaninfo_store%nchans) &
2023  + (diag_chaninfo_store%var_usage(var_index) - 1)) = chaninfo_value
2024  end subroutine nc_diag_chaninfo_short
2025 
2026  ! Add a single scalar long integer to the given chaninfo
2027  ! variable. (This is NOT a NetCDF "long", just a NetCDF "int".)
2028  !
2029  ! This adds a single value to the specified chaninfo variable.
2030  !
2031  ! If the variable does not already exist, it will be created,
2032  ! and the value will be inserted as the variable's first
2033  ! element.
2034  !
2035  ! Otherwise, the value will be inserted to the next empty spot.
2036  !
2037  ! Values are inserted in the order of the calls made. As such,
2038  ! this subroutine is best designed to be used in a loop, where
2039  ! for every channel iteration, a value is added using this
2040  ! subroutine.
2041  !
2042  ! This is an internal subroutine, and is NOT meant to be called
2043  ! outside of nc_diag_write. Calling this subroutine in your
2044  ! program may result in unexpected behavior and/or data
2045  ! corruption! (You should use the generic nc_diag_chaninfo
2046  ! instead!)
2047  !
2048  ! Args:
2049  ! chaninfo_name (character(len=*)): the chaninfo variable
2050  ! to store to.
2051  ! chaninfo_value (integer(i_long)): the value to store.
2052  !
2053  ! Raises:
2054  ! If the data has already been locked, this will result in
2055  ! an error.
2056  !
2057  ! If definitions have already been locked, and a new
2058  ! variable is being created, this will result in an error.
2059  !
2060  ! If the variable is already full (e.g. it has nchans number
2061  ! of elements), this will result in an error.
2062  !
2063  ! The following errors will trigger indirectly from other
2064  ! subroutines called here:
2065  !
2066  ! If nchans has not been set yet, this will result in an
2067  ! error.
2068  !
2069  ! If there is no file open (or the file is already closed),
2070  ! this will result in an error.
2071  !
2072  ! Other errors may result from invalid data storage, NetCDF
2073  ! errors, or even a bug. See the called subroutines'
2074  ! documentation for details.
2075  !
2076  subroutine nc_diag_chaninfo_long(chaninfo_name, chaninfo_value)
2077  character(len=*), intent(in) :: chaninfo_name
2078  integer(i_long), intent(in) :: chaninfo_value
2079 
2080  integer(i_long) :: i, var_index, var_rel_index, type_index
2081 
2082 #ifdef ENABLE_ACTION_MSGS
2083  character(len=1000) :: action_str
2084 
2085  if (nclayer_enable_action) then
2086  write(action_str, "(A, I0, A)") "nc_diag_chaninfo_long(chaninfo_name = " // chaninfo_name // ", chaninfo_value = ", chaninfo_value, ")"
2087  call nclayer_actionm(trim(action_str))
2088  end if
2089 #endif
2090  ! Make sure that data hasn't been locked
2091  if (diag_chaninfo_store%data_lock) then
2092  call nclayer_error("Can't add new data - data have already been written and locked!")
2093  end if
2094 
2095  ! For long, type index is 3
2096  type_index = 3
2097 
2098  ! Default to -1
2099  var_index = -1
2100 
2101  ! Attempt to match the variable name + fetch the variable
2102  ! index first!
2103  do i = 1, diag_chaninfo_store%total
2104  if (diag_chaninfo_store%names(i) == chaninfo_name) then
2105  var_rel_index = diag_chaninfo_store%var_rel_pos(i)
2106  var_index = i
2107  exit
2108  end if
2109  end do
2110 
2111 #ifdef _DEBUG_MEM_
2112  print *, " *** chaninfo_name"
2113  print *, chaninfo_name
2114  print *, " *** var_index is set to:"
2115  print *, var_index
2116 #endif
2117 
2118  if (var_index == -1) then
2119  ! Entry does not exist yet...
2120 
2121  ! First, check to make sure we can still define new variables.
2122  if (diag_chaninfo_store%def_lock) then
2123  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
2124  end if
2125 
2126  ! Expand variable metadata first!
2127  ! Make sure we have enough variable metadata storage
2128  ! (and if not, reallocate!)
2130 
2131  ! Add to the total!
2132  diag_chaninfo_store%total = diag_chaninfo_store%total + 1
2133 
2134  ! Store name and type!
2135  diag_chaninfo_store%names(diag_chaninfo_store%total) = chaninfo_name
2137 
2138  ! We just need to add one entry...
2139  ! Call resize subroutine to ensure we have enough space
2140  ! (and if not, realloc!)
2142 
2143  ! Now add a relative position... based on the next position!
2144 
2145  ! First, increment the number of variables stored for this type:
2146  diag_chaninfo_store%acount_v(type_index) = diag_chaninfo_store%acount_v(type_index) + 1
2147 
2148  ! Then, set the next variable's relative positioning,
2149  ! based on the number of variables stored for this type.
2150  diag_chaninfo_store%var_rel_pos(diag_chaninfo_store%total) = diag_chaninfo_store%acount_v(type_index)
2151 
2152  ! Initialize the amount of memory used to 1.
2153  diag_chaninfo_store%var_usage(diag_chaninfo_store%total) = 1
2154 
2155  ! Set var_index to the total
2156  var_index = diag_chaninfo_store%total
2157  else
2158  ! Variable already exists!
2159 
2160  ! Check to make sure we can fit more data!
2161  ! (# data < nchans)
2162  if (diag_chaninfo_store%var_usage(var_index) + &
2163  diag_chaninfo_store%rel_indexes(var_index) >= diag_chaninfo_store%nchans) then
2164 #ifdef _DEBUG_MEM_
2165  print *, "!!!! diag_chaninfo_store%var_usage(var_index)"
2166  print *, diag_chaninfo_store%var_usage(var_index)
2167 #endif
2168  call nclayer_error("Can't add new data - data added is exceeding nchan! Data must fit within nchan constraint.")
2169  endif
2170 
2171  ! Increment current variable count
2172  diag_chaninfo_store%var_usage(var_index) = &
2173  diag_chaninfo_store%var_usage(var_index) + 1
2174  end if
2175 
2176  ! Now add the actual entry!
2177 #ifdef _DEBUG_MEM_
2178  print *, "===================================="
2179  print *, "diag_chaninfo_store%total"
2180  print *, diag_chaninfo_store%total
2181  print *, "var_index"
2182  print *, var_index
2183  print *, "diag_chaninfo_store%var_rel_pos(var_index)"
2184  print *, diag_chaninfo_store%var_rel_pos(var_index)
2185  print *, "diag_chaninfo_store%nchans"
2186  print *, diag_chaninfo_store%nchans
2187  print *, "diag_chaninfo_store%var_usage(var_index)"
2188  print *, diag_chaninfo_store%var_usage(var_index)
2189  print *, "===================================="
2190 #endif
2191 
2192  diag_chaninfo_store%ci_long(1 + &
2193  ((diag_chaninfo_store%var_rel_pos(var_index) - 1) &
2194  * diag_chaninfo_store%nchans) &
2195  + (diag_chaninfo_store%var_usage(var_index) - 1)) = chaninfo_value
2196  end subroutine nc_diag_chaninfo_long
2197 
2198  ! Add a single scalar float to the given chaninfo variable.
2199  !
2200  ! This adds a single value to the specified chaninfo variable.
2201  !
2202  ! If the variable does not already exist, it will be created,
2203  ! and the value will be inserted as the variable's first
2204  ! element.
2205  !
2206  ! Otherwise, the value will be inserted to the next empty spot.
2207  !
2208  ! Values are inserted in the order of the calls made. As such,
2209  ! this subroutine is best designed to be used in a loop, where
2210  ! for every channel iteration, a value is added using this
2211  ! subroutine.
2212  !
2213  ! This is an internal subroutine, and is NOT meant to be called
2214  ! outside of nc_diag_write. Calling this subroutine in your
2215  ! program may result in unexpected behavior and/or data
2216  ! corruption! (You should use the generic nc_diag_chaninfo
2217  ! instead!)
2218  !
2219  ! Args:
2220  ! chaninfo_name (character(len=*)): the chaninfo variable
2221  ! to store to.
2222  ! chaninfo_value (real(r_single)): the value to store.
2223  !
2224  ! Raises:
2225  ! If the data has already been locked, this will result in
2226  ! an error.
2227  !
2228  ! If definitions have already been locked, and a new
2229  ! variable is being created, this will result in an error.
2230  !
2231  ! If the variable is already full (e.g. it has nchans number
2232  ! of elements), this will result in an error.
2233  !
2234  ! The following errors will trigger indirectly from other
2235  ! subroutines called here:
2236  !
2237  ! If nchans has not been set yet, this will result in an
2238  ! error.
2239  !
2240  ! If there is no file open (or the file is already closed),
2241  ! this will result in an error.
2242  !
2243  ! Other errors may result from invalid data storage, NetCDF
2244  ! errors, or even a bug. See the called subroutines'
2245  ! documentation for details.
2246  !
2247  subroutine nc_diag_chaninfo_rsingle(chaninfo_name, chaninfo_value)
2248  character(len=*), intent(in) :: chaninfo_name
2249  real(r_single), intent(in) :: chaninfo_value
2250 
2251  integer(i_long) :: i, var_index, var_rel_index, type_index
2252 
2253 #ifdef ENABLE_ACTION_MSGS
2254  character(len=1000) :: action_str
2255 
2256  if (nclayer_enable_action) then
2257  write(action_str, "(A, F0.5, A)") "nc_diag_chaninfo_rsingle(chaninfo_name = " // chaninfo_name // ", chaninfo_value = ", chaninfo_value, ")"
2258  call nclayer_actionm(trim(action_str))
2259  end if
2260 #endif
2261  ! Make sure that data hasn't been locked
2262  if (diag_chaninfo_store%data_lock) then
2263  call nclayer_error("Can't add new data - data have already been written and locked!")
2264  end if
2265 
2266  ! For rsingle, type index is 4
2267  type_index = 4
2268 
2269  ! Default to -1
2270  var_index = -1
2271 
2272  ! Attempt to match the variable name + fetch the variable
2273  ! index first!
2274  do i = 1, diag_chaninfo_store%total
2275  if (diag_chaninfo_store%names(i) == chaninfo_name) then
2276  var_rel_index = diag_chaninfo_store%var_rel_pos(i)
2277  var_index = i
2278  exit
2279  end if
2280  end do
2281 
2282 #ifdef _DEBUG_MEM_
2283  print *, " *** chaninfo_name"
2284  print *, chaninfo_name
2285  print *, " *** var_index is set to:"
2286  print *, var_index
2287 #endif
2288 
2289  if (var_index == -1) then
2290  ! Entry does not exist yet...
2291 
2292  ! First, check to make sure we can still define new variables.
2293  if (diag_chaninfo_store%def_lock) then
2294  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
2295  end if
2296 
2297  ! Expand variable metadata first!
2298  ! Make sure we have enough variable metadata storage
2299  ! (and if not, reallocate!)
2301 
2302  ! Add to the total!
2303  diag_chaninfo_store%total = diag_chaninfo_store%total + 1
2304 
2305  ! Store name and type!
2306  diag_chaninfo_store%names(diag_chaninfo_store%total) = chaninfo_name
2308 
2309  ! We just need to add one entry...
2310  ! Call resize subroutine to ensure we have enough space
2311  ! (and if not, realloc!)
2313 
2314  ! Now add a relative position... based on the next position!
2315 
2316  ! First, increment the number of variables stored for this type:
2317  diag_chaninfo_store%acount_v(type_index) = diag_chaninfo_store%acount_v(type_index) + 1
2318 
2319  ! Then, set the next variable's relative positioning,
2320  ! based on the number of variables stored for this type.
2321  diag_chaninfo_store%var_rel_pos(diag_chaninfo_store%total) = diag_chaninfo_store%acount_v(type_index)
2322 
2323  ! Initialize the amount of memory used to 1.
2324  diag_chaninfo_store%var_usage(diag_chaninfo_store%total) = 1
2325 
2326  ! Set var_index to the total
2327  var_index = diag_chaninfo_store%total
2328  else
2329  ! Variable already exists!
2330 
2331  ! Check to make sure we can fit more data!
2332  ! (# data < nchans)
2333  if (diag_chaninfo_store%var_usage(var_index) + &
2334  diag_chaninfo_store%rel_indexes(var_index) >= diag_chaninfo_store%nchans) then
2335  call nclayer_error("Can't add new data - data added is exceeding nchan! Data must fit within nchan constraint.")
2336  endif
2337 
2338  ! Increment current variable count
2339  diag_chaninfo_store%var_usage(var_index) = &
2340  diag_chaninfo_store%var_usage(var_index) + 1
2341  end if
2342 
2343 #ifdef _DEBUG_MEM_
2344  print *, "===================================="
2345  print *, "diag_chaninfo_store%total"
2346  print *, diag_chaninfo_store%total
2347  print *, "var_index"
2348  print *, var_index
2349  print *, "diag_chaninfo_store%var_rel_pos(var_index)"
2350  print *, diag_chaninfo_store%var_rel_pos(var_index)
2351  print *, "diag_chaninfo_store%nchans"
2352  print *, diag_chaninfo_store%nchans
2353  print *, "diag_chaninfo_store%var_usage(var_index)"
2354  print *, diag_chaninfo_store%var_usage(var_index)
2355  print *, "===================================="
2356 #endif
2357 
2358  ! Now add the actual entry!
2359  diag_chaninfo_store%ci_rsingle(1 + &
2360  ((diag_chaninfo_store%var_rel_pos(var_index) - 1) &
2361  * diag_chaninfo_store%nchans) &
2362  + (diag_chaninfo_store%var_usage(var_index) - 1)) = chaninfo_value
2363  end subroutine nc_diag_chaninfo_rsingle
2364 
2365  ! Add a single scalar double to the given chaninfo variable.
2366  !
2367  ! This adds a single value to the specified chaninfo variable.
2368  !
2369  ! If the variable does not already exist, it will be created,
2370  ! and the value will be inserted as the variable's first
2371  ! element.
2372  !
2373  ! Otherwise, the value will be inserted to the next empty spot.
2374  !
2375  ! Values are inserted in the order of the calls made. As such,
2376  ! this subroutine is best designed to be used in a loop, where
2377  ! for every channel iteration, a value is added using this
2378  ! subroutine.
2379  !
2380  ! This is an internal subroutine, and is NOT meant to be called
2381  ! outside of nc_diag_write. Calling this subroutine in your
2382  ! program may result in unexpected behavior and/or data
2383  ! corruption! (You should use the generic nc_diag_chaninfo
2384  ! instead!)
2385  !
2386  ! Args:
2387  ! chaninfo_name (character(len=*)): the chaninfo variable
2388  ! to store to.
2389  ! chaninfo_value (real(r_double)): the value to store.
2390  !
2391  ! Raises:
2392  ! If the data has already been locked, this will result in
2393  ! an error.
2394  !
2395  ! If definitions have already been locked, and a new
2396  ! variable is being created, this will result in an error.
2397  !
2398  ! If the variable is already full (e.g. it has nchans number
2399  ! of elements), this will result in an error.
2400  !
2401  ! The following errors will trigger indirectly from other
2402  ! subroutines called here:
2403  !
2404  ! If nchans has not been set yet, this will result in an
2405  ! error.
2406  !
2407  ! If there is no file open (or the file is already closed),
2408  ! this will result in an error.
2409  !
2410  ! Other errors may result from invalid data storage, NetCDF
2411  ! errors, or even a bug. See the called subroutines'
2412  ! documentation for details.
2413  !
2414  subroutine nc_diag_chaninfo_rdouble(chaninfo_name, chaninfo_value)
2415  character(len=*), intent(in) :: chaninfo_name
2416  real(r_double), intent(in) :: chaninfo_value
2417 
2418  integer(i_long) :: i, var_index, var_rel_index, type_index
2419 
2420 #ifdef ENABLE_ACTION_MSGS
2421  character(len=1000) :: action_str
2422 
2423  if (nclayer_enable_action) then
2424  write(action_str, "(A, F0.5, A)") "nc_diag_chaninfo_rdouble(chaninfo_name = " // chaninfo_name // ", chaninfo_value = ", chaninfo_value, ")"
2425  call nclayer_actionm(trim(action_str))
2426  end if
2427 #endif
2428  ! Make sure that data hasn't been locked
2429  if (diag_chaninfo_store%data_lock) then
2430  call nclayer_error("Can't add new data - data have already been written and locked!")
2431  end if
2432 
2433  ! For rdouble, type index is 5
2434  type_index = 5
2435 
2436  ! Default to -1
2437  var_index = -1
2438 
2439  ! Attempt to match the variable name + fetch the variable
2440  ! index first!
2441  do i = 1, diag_chaninfo_store%total
2442  if (diag_chaninfo_store%names(i) == chaninfo_name) then
2443  var_rel_index = diag_chaninfo_store%var_rel_pos(i)
2444  var_index = i
2445  exit
2446  end if
2447  end do
2448 
2449  if (var_index == -1) then
2450  ! Entry does not exist yet...
2451 
2452  ! First, check to make sure we can still define new variables.
2453  if (diag_chaninfo_store%def_lock) then
2454  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
2455  end if
2456 
2457  ! Expand variable metadata first!
2458  ! Make sure we have enough variable metadata storage
2459  ! (and if not, reallocate!)
2461 
2462  ! Add to the total!
2463  diag_chaninfo_store%total = diag_chaninfo_store%total + 1
2464 
2465  ! Store name and type!
2466  diag_chaninfo_store%names(diag_chaninfo_store%total) = chaninfo_name
2468 
2469  ! We just need to add one entry...
2470  ! Call resize subroutine to ensure we have enough space
2471  ! (and if not, realloc!)
2473 
2474  ! Now add a relative position... based on the next position!
2475 
2476  ! First, increment the number of variables stored for this type:
2477  diag_chaninfo_store%acount_v(type_index) = diag_chaninfo_store%acount_v(type_index) + 1
2478 
2479  ! Then, set the next variable's relative positioning,
2480  ! based on the number of variables stored for this type.
2481  diag_chaninfo_store%var_rel_pos(diag_chaninfo_store%total) = diag_chaninfo_store%acount_v(type_index)
2482 
2483  ! Initialize the amount of memory used to 1.
2484  diag_chaninfo_store%var_usage(diag_chaninfo_store%total) = 1
2485 
2486  ! Set var_index to the total
2487  var_index = diag_chaninfo_store%total
2488  else
2489  ! Variable already exists!
2490 
2491  ! Check to make sure we can fit more data!
2492  ! (# data < nchans)
2493  if (diag_chaninfo_store%var_usage(var_index) + &
2494  diag_chaninfo_store%rel_indexes(var_index) >= diag_chaninfo_store%nchans) then
2495  call nclayer_error("Can't add new data - data added is exceeding nchan! Data must fit within nchan constraint.")
2496  endif
2497 
2498  ! Increment current variable count
2499  diag_chaninfo_store%var_usage(var_index) = &
2500  diag_chaninfo_store%var_usage(var_index) + 1
2501  end if
2502 
2503  ! Now add the actual entry!
2504  diag_chaninfo_store%ci_rdouble(1 + &
2505  ((diag_chaninfo_store%var_rel_pos(var_index) - 1) &
2506  * diag_chaninfo_store%nchans) &
2507  + (diag_chaninfo_store%var_usage(var_index) - 1)) = chaninfo_value
2508  end subroutine nc_diag_chaninfo_rdouble
2509 
2510  ! Add a single scalar string to the given chaninfo variable.
2511  ! (This uses the NetCDF char type, stored internally as a 2D
2512  ! array of characters.)
2513  !
2514  ! This adds a single value to the specified chaninfo variable.
2515  !
2516  ! If the variable does not already exist, it will be created,
2517  ! and the value will be inserted as the variable's first
2518  ! element.
2519  !
2520  ! Otherwise, the value will be inserted to the next empty spot.
2521  !
2522  ! Values are inserted in the order of the calls made. As such,
2523  ! this subroutine is best designed to be used in a loop, where
2524  ! for every channel iteration, a value is added using this
2525  ! subroutine.
2526  !
2527  ! This is an internal subroutine, and is NOT meant to be called
2528  ! outside of nc_diag_write. Calling this subroutine in your
2529  ! program may result in unexpected behavior and/or data
2530  ! corruption! (You should use the generic nc_diag_chaninfo
2531  ! instead!)
2532  !
2533  ! Args:
2534  ! chaninfo_name (character(len=*)): the chaninfo variable
2535  ! to store to.
2536  ! chaninfo_value (character(len=*)): the value to store.
2537  !
2538  ! Raises:
2539  ! If the data has already been locked, this will result in
2540  ! an error.
2541  !
2542  ! If definitions have already been locked, and a new
2543  ! variable is being created, this will result in an error.
2544  !
2545  ! If the variable is already full (e.g. it has nchans number
2546  ! of elements), this will result in an error.
2547  !
2548  ! The following errors will trigger indirectly from other
2549  ! subroutines called here:
2550  !
2551  ! If nchans has not been set yet, this will result in an
2552  ! error.
2553  !
2554  ! If there is no file open (or the file is already closed),
2555  ! this will result in an error.
2556  !
2557  ! Other errors may result from invalid data storage, NetCDF
2558  ! errors, or even a bug. See the called subroutines'
2559  ! documentation for details.
2560  !
2561  subroutine nc_diag_chaninfo_string(chaninfo_name, chaninfo_value)
2562  character(len=*), intent(in) :: chaninfo_name
2563  character(len=*), intent(in) :: chaninfo_value
2564 
2565  integer(i_long) :: i, var_index, var_rel_index, type_index
2566 
2567 #ifdef ENABLE_ACTION_MSGS
2568  character(len=1000) :: action_str
2569 
2570  if (nclayer_enable_action) then
2571  write(action_str, "(A)") "nc_diag_chaninfo_string(chaninfo_name = " // chaninfo_name // ", chaninfo_value = " // trim(chaninfo_value) // ")"
2572  call nclayer_actionm(trim(action_str))
2573  end if
2574 #endif
2575  ! Make sure that data hasn't been locked
2576  if (diag_chaninfo_store%data_lock) then
2577  call nclayer_error("Can't add new data - data have already been written and locked!")
2578  end if
2579 
2580  ! For string, type index is 6
2581  type_index = 6
2582 
2583  ! Default to -1
2584  var_index = -1
2585 
2586  ! Attempt to match the variable name + fetch the variable
2587  ! index first!
2588  do i = 1, diag_chaninfo_store%total
2589  if (diag_chaninfo_store%names(i) == chaninfo_name) then
2590  var_rel_index = diag_chaninfo_store%var_rel_pos(i)
2591  var_index = i
2592  exit
2593  end if
2594  end do
2595 
2596  if (var_index == -1) then
2597  ! Entry does not exist yet...
2598 
2599  ! First, check to make sure we can still define new variables.
2600  if (diag_chaninfo_store%def_lock) then
2601  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
2602  end if
2603 
2604  ! Expand variable metadata first!
2605  ! Make sure we have enough variable metadata storage
2606  ! (and if not, reallocate!)
2608 
2609  ! Add to the total!
2610  diag_chaninfo_store%total = diag_chaninfo_store%total + 1
2611 
2612  ! Store name and type!
2613  diag_chaninfo_store%names(diag_chaninfo_store%total) = chaninfo_name
2615 
2616  ! We just need to add one entry...
2617  ! Call resize subroutine to ensure we have enough space
2618  ! (and if not, realloc!)
2620 
2621  ! Now add a relative position... based on the next position!
2622 
2623  ! First, increment the number of variables stored for this type:
2624  diag_chaninfo_store%acount_v(type_index) = diag_chaninfo_store%acount_v(type_index) + 1
2625 
2626  ! Then, set the next variable's relative positioning,
2627  ! based on the number of variables stored for this type.
2628  diag_chaninfo_store%var_rel_pos(diag_chaninfo_store%total) = diag_chaninfo_store%acount_v(type_index)
2629 
2630  ! Initialize the amount of memory used to 1.
2631  diag_chaninfo_store%var_usage(diag_chaninfo_store%total) = 1
2632 
2633  ! Set var_index to the total
2634  var_index = diag_chaninfo_store%total
2635  else
2636  ! Variable already exists!
2637 
2638  ! Check to make sure we can fit more data!
2639  ! (# data < nchans)
2640  if (diag_chaninfo_store%var_usage(var_index) + &
2641  diag_chaninfo_store%rel_indexes(var_index) >= diag_chaninfo_store%nchans) then
2642  call nclayer_error("Can't add new data - data added is exceeding nchan! Data must fit within nchan constraint.")
2643  endif
2644 
2645  ! Check max string length
2646  if ((diag_chaninfo_store%def_lock) .AND. &
2647  (len_trim(chaninfo_value) > diag_chaninfo_store%max_str_lens(var_index))) &
2648  call nclayer_error("Cannot expand variable string length after locking variable definitions!")
2649 
2650  ! Increment current variable count
2651  diag_chaninfo_store%var_usage(var_index) = &
2652  diag_chaninfo_store%var_usage(var_index) + 1
2653  end if
2654 
2655  ! If trim isn't enabled, set our maximum string length here!
2656  if (.NOT. enable_trim) then
2657  if (diag_chaninfo_store%max_str_lens(var_index) == -1) then
2658  diag_chaninfo_store%max_str_lens(var_index) = len(chaninfo_value)
2659  else
2660  ! Validate that our non-first value isn't different from
2661  ! the initial string length
2662  if (diag_chaninfo_store%max_str_lens(var_index) /= len(chaninfo_value)) &
2663  call nclayer_error("Cannot change string size when trimming is disabled!")
2664  end if
2665  end if
2666 
2667  ! Now add the actual entry!
2668  diag_chaninfo_store%ci_string(1 + &
2669  ((diag_chaninfo_store%var_rel_pos(var_index) - 1) &
2670  * diag_chaninfo_store%nchans) &
2671  + (diag_chaninfo_store%var_usage(var_index) - 1)) = chaninfo_value
2672  end subroutine nc_diag_chaninfo_string
2673 end module ncdw_chaninfo
subroutine nc_diag_chaninfo_resize_string(addl_num_entries, update_acount_in)
type(diag_chaninfo), allocatable diag_chaninfo_store
Definition: ncdw_state.f90:16
integer(i_long), parameter nlayer_multi_base
Definition: ncdw_types.F90:27
subroutine nc_diag_chaninfo_resize_short(addl_num_entries, update_acount_in)
integer, parameter, public i_byte
Definition: ncd_kinds.F90:45
subroutine nc_diag_chaninfo_rdouble(chaninfo_name, chaninfo_value)
subroutine nc_diag_chaninfo_write_def(internal)
integer(i_byte), parameter nlayer_short
Definition: ncdw_types.F90:11
subroutine nc_diag_chaninfo_short(chaninfo_name, chaninfo_value)
integer(i_long), parameter nlayer_fill_long
Definition: ncdw_types.F90:31
integer, parameter, public i_long
Definition: ncd_kinds.F90:47
integer(i_byte), parameter nlayer_fill_byte
Definition: ncdw_types.F90:29
logical init_done
Definition: ncdw_state.f90:9
integer(i_byte), parameter nlayer_double
Definition: ncdw_types.F90:14
logical enable_trim
Definition: ncdw_state.f90:12
subroutine nc_diag_chaninfo_expand
subroutine nc_diag_chaninfo_resize_byte(addl_num_entries, update_acount_in)
subroutine nc_diag_chaninfo_resize_long(addl_num_entries, update_acount_in)
integer function max_len_string_array(str_arr, arr_length)
integer(i_byte), parameter nlayer_string
Definition: ncdw_types.F90:15
integer(i_long) ncid
Definition: ncdw_state.f90:8
subroutine nc_diag_chaninfo_dim_set(nchans)
integer(i_long), parameter nlayer_compression
Definition: ncdw_types.F90:21
integer(i_short), parameter nlayer_fill_short
Definition: ncdw_types.F90:30
logical append_only
Definition: ncdw_state.f90:10
subroutine nc_diag_chaninfo_set_strict(enable_strict)
subroutine nc_diag_chaninfo_resize_rsingle(addl_num_entries, update_acount_in)
subroutine nc_diag_chaninfo_write_data(flush_data_only)
integer, parameter, public i_short
Definition: ncd_kinds.F90:46
subroutine nc_diag_chaninfo_prealloc_vars_storage(nclayer_type, num_of_addl_slots)
integer(i_byte), parameter nlayer_byte
Definition: ncdw_types.F90:10
character, parameter nlayer_fill_char
Definition: ncdw_types.F90:34
subroutine nc_diag_varattr_add_var(var_name, var_type, var_id)
integer(i_short), parameter nlayer_default_ent
Definition: ncdw_types.F90:18
subroutine nc_diag_chaninfo_resize_rdouble(addl_num_entries, update_acount_in)
subroutine nc_diag_chaninfo_rsingle(chaninfo_name, chaninfo_value)
subroutine nc_diag_chaninfo_allocmulti(multiplier)
subroutine nc_diag_chaninfo_load_def
real(r_single), parameter nlayer_fill_float
Definition: ncdw_types.F90:32
integer, parameter, public r_double
Definition: ncd_kinds.F90:80
subroutine nc_diag_chaninfo_prealloc_vars(num_of_addl_vars)
integer, parameter, public r_single
Definition: ncd_kinds.F90:79
integer, parameter, public i_llong
Definition: ncd_kinds.F90:49
integer(i_byte), parameter nlayer_float
Definition: ncdw_types.F90:13
integer(i_byte), parameter nlayer_long
Definition: ncdw_types.F90:12
subroutine nc_diag_chaninfo_string(chaninfo_name, chaninfo_value)
real(r_double), parameter nlayer_fill_double
Definition: ncdw_types.F90:33
subroutine nc_diag_chaninfo_long(chaninfo_name, chaninfo_value)
subroutine nc_diag_chaninfo_byte(chaninfo_name, chaninfo_value)