FV3 Bundle
ncdw_data2d.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 ! data2d module - ncdw_data2d
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 numeric
43  ! 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
53  !
54  ! The derived type used, diag_chaninfo, has these variables to help
55  ! 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 data
59  ! in. However, don't mistaken "throw in" with "disorganized" -
60  ! chaninfo uses a very structured format for these variables!
61  ! Keep reading to find out how we structure it...
62  !
63  ! -> nchans - the number of channels to use. Remember that chaninfo
64  ! variables have dimensions 1 x nchans - basically, we need to
65  ! store nchans values. We'll need this a LOT to do consistency
66  ! checks, and to keep track of everything!
67  !
68  ! -> names - all of the chaninfo variable names! We'll be using
69  ! this array to store and lookup chaninfo variables, as well as
70  ! storing them!
71  !
72  ! -> types - all of the chaninfo variable types! These are byte
73  ! integers that get compared to our NLAYER_* type constants
74  ! (see: ncdw_types.F90).
75  !
76  ! -> var_usage - the amount of entries we've stored in our chaninfo
77  ! variable! For instance, if we called
78  ! nc_diag_chaninfo("myvar", 1) three times, for that particular
79  ! var_usage(i), we would have an entry of 3.
80  !
81  ! -> var_rel_pos - the star of the show! This is an abbreviation
82  ! of "variable relative positioning". Recall that we store
83  ! our variable data in ci_* specific type arrays. We know
84  ! the nchans amount, and we know the type. This variable stores
85  ! the "block" that our data is in within the type array.
86  !
87  ! Therefore, we can use the equation to find our starting
88  ! position: 1 + [(REL_VAL - 1) * nchans]
89  !
90  ! For instance, if var_rel_pos(1) for variable names(1) = "bla"
91  ! is set to 2, nchans is 10, and the type is NLAYER_FLOAT, we
92  ! can deduce that in ci_rsingle, our data can be found starting
93  ! at 1 + (1 * 10) = 11. This makes sense, as seen with our mini
94  ! diagram below:
95  !
96  ! ci_rsingle:
97  ! / ci_rsingle index \
98  ! 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20
99  ! [ x, x, x, x, x, x, x, x, x, x, y, y, y, y, y, y, y, y, y, y ]
100  ! \ ci_rsingle array /
101  !
102  ! Indeed, our second block does start at index 11!
103  ! As a bonus, since our data is in blocks, things can be super
104  ! fast since we're just cutting our big array into small ones!
105  !
106  ! -> acount_v: Finally, we have dynamic allocation. We have in our
107  ! type a variable called acount_v. This tells us how many
108  ! variables are stored in each type. Using the same equation
109  ! above, and combining with var_usage, we can figure out where
110  ! to put our data!
111  !
112  ! Assume var_usage(i) = 2, block starts at index 11 with the
113  ! equation above.
114  !
115  ! Again, with our fun little diagram:
116  !
117  ! ci_rsingle:
118  ! / ci_rsingle index \
119  ! 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20
120  ! [ x, x, x, x, x, x, x, x, x, x, y, y, Y, y, y, y, y, y, y, y ]
121  ! [ BLOCK 1 SEEK = 1->10->11 ][var_u=2|---block 2 area 11->20]
122  ! \ ci_rsingle array /
123  !
124  ! The capital Y marks the place we store our data!
125  !
126  ! For the non-data variables (e.g. variable names, types, etc.),
127  ! they are indexed by variable insertion order. This allows for
128  ! easy lookup by looking up the variable name, and using the
129  ! resulting index for fetching other information.
130  !
131  ! Example:
132  ! names: [ 'asdf', 'ghjk', 'zxcv' ]
133  ! types: [ BYTE, FLOAT, BYTE ]
134  ! var_rel_pos: [ 1, 1, 2 ]
135  !
136  ! Lookup: "ghjk", result index = 2
137  !
138  ! Therefore, the "ghjk" variable type is types(2) = FLOAT, and
139  ! the var_rel_pos for "ghjk" variable is var_rel_pos(2) = 1.
140  !
141  ! These variables are allocated and reallocated, as needed.
142  !
143  ! For the variable metadata fields (variable names, types,
144  ! relative indicies, etc.), these are reallocated incrementally
145  ! when a new variable is added.
146  !
147  ! For the data storage fields, these are reallocated incrementally
148  ! when new data is added.
149  !
150  ! Initial allocation and subsequent reallocation is done by
151  ! chunks. Allocating one element and/or reallocating and adding
152  ! just one element is inefficient, since it's likely that much
153  ! more data (and variables) will be added. Thus, allocation and
154  ! reallocation is done by (re-)allocating exponentially increasing
155  ! chunk sizes. See nc_diag_chaninfo_allocmulti help for more
156  ! details.
157  !
158  ! Because (re-)allocation is done in chunks, we keep a count of
159  ! how much of the memory we're using so that we know when it's
160  ! time to (re-)allocate. Once we need to (re-)allocate, we
161  ! perform it, and then update our total memory counter to keep
162  ! track of the memory already allocated.
163  !
164  ! With all of these variables (and a few more state variables),
165  ! we can reliably store our chaninfo data quickly and
166  ! efficiently!
167  !
168 
169  use ncd_kinds, only: i_byte, i_short, i_long, i_llong, r_single, &
170  r_double
171  use ncdw_state, only: init_done, append_only, ncid, &
172  enable_trim, &
180  use ncdw_strarrutils, only: &
181 #ifdef _debug_mem_
182  string_array_dump, &
183 #endif
184  max_len_string_array, max_len_notrim_string_array
187 
193  use ncdw_realloc, only: nc_diag_realloc
194 
195  use netcdf, only: nf90_inquire, nf90_inquire_variable, &
196  nf90_inquire_dimension, nf90_def_dim, nf90_def_var, &
197  nf90_def_var_deflate, nf90_def_var_chunking, nf90_put_var, &
198  nf90_byte, nf90_short, nf90_int, nf90_float, nf90_double, &
199  nf90_char, nf90_max_name, nf90_chunked
200 
201  use ncdw_climsg, only: &
202 #ifdef enable_action_msgs
203  nclayer_enable_action, nclayer_actionm, &
204 #endif
205 #ifdef _DEBUG_MEM_
206  nclayer_debug, &
207 #endif
208  nclayer_error, nclayer_warning, nclayer_info, nclayer_check
209 
210  implicit none
211 
212  interface nc_diag_data2d
213  module procedure nc_diag_data2d_byte, &
217  end interface nc_diag_data2d
218 
219  contains
220  subroutine nc_diag_data2d_allocmulti(multiplier)
221  integer(i_long), intent(in) :: multiplier
222  if (init_done) then
223  ! # of times we needed to realloc simple data2d
224  ! also the multiplier factor for allocation (2^x)
225  diag_data2d_store%alloc_s_multi = multiplier
226 
227  ! # of times we needed to realloc data2d data storage
228  ! also the multiplier factor for allocation (2^x)
229  diag_data2d_store%alloc_m_multi = multiplier
230 
231  ! # of times we needed to realloc data2d INDEX data storage
232  ! also the multiplier factor for allocation (2^x)
233  diag_data2d_store%alloc_mi_multi = multiplier
234  end if
235  end subroutine nc_diag_data2d_allocmulti
236 
237  function nc_diag_data2d_max_len_var(var_index) result(max_len)
238  integer(i_llong), intent(in) :: var_index
239 
240  integer :: i, max_len
241 
242  character(len=1000) :: data_uneven_msg
243 
244  max_len = -1
245 
246  do i = 1, diag_data2d_store%stor_i_arr(var_index)%icount
247  ! Only show a message if strict checking is enabled.
248  ! Otherwise, show the message later in data writing.
249  if (diag_data2d_store%strict_check .AND. &
250  (diag_data2d_store%stor_i_arr(var_index)%length_arr(i) /= max_len) .AND. &
251  (max_len /= -1)) then
252  ! Show message!
253  ! NOTE - I0 and TRIM are Fortran 95 specs
254  write (data_uneven_msg, "(A, I0, A, I0, A)") "Amount of data written in " // &
255  trim(diag_data2d_store%names(var_index)) // " (", &
256  diag_data2d_store%stor_i_arr(var_index)%length_arr(i), &
257  ")" // char(10) // &
258  " does not match the variable length" // &
259  " (", max_len, ")!"
260 
261  ! Probably not needed, since this only triggers on a
262  ! strict check... but just in case...
263  if (diag_data2d_store%strict_check) then
264  call nclayer_error(trim(data_uneven_msg))
265  else
266  call nclayer_warning(trim(data_uneven_msg))
267  end if
268  end if
269 
270  if (diag_data2d_store%stor_i_arr(var_index)%length_arr(i) > max_len) &
271  max_len = diag_data2d_store%stor_i_arr(var_index)%length_arr(i)
272  end do
273  end function nc_diag_data2d_max_len_var
274 
275  subroutine nc_diag_data2d_load_def
276  integer(i_long) :: ndims, nvars, var_index, type_index
277  integer(i_long) :: rel_index, i, nobs_size
278 
279  character(len=NF90_MAX_NAME) :: tmp_var_name
280  integer(i_long) :: tmp_var_type, tmp_var_ndims
281 
282  integer(i_long), dimension(:), allocatable :: tmp_var_dimids, tmp_var_dim_sizes
283  character(len=NF90_MAX_NAME) , allocatable :: tmp_var_dim_names(:)
284 
285  logical :: is_data2d_var
286 
287  ! Get top level info about the file!
288  call nclayer_check(nf90_inquire(ncid, ndimensions = ndims, &
289  nvariables = nvars))
290 
291  ! Now search for variables that use data2d storage!
292  ! Loop through each variable!
293  do var_index = 1, nvars
294  ! Grab number of dimensions and attributes first
295  call nclayer_check(nf90_inquire_variable(ncid, var_index, name = tmp_var_name, ndims = tmp_var_ndims))
296 
297  ! Allocate temporary variable dimids storage!
298  allocate(tmp_var_dimids(tmp_var_ndims))
299  allocate(tmp_var_dim_names(tmp_var_ndims))
300  allocate(tmp_var_dim_sizes(tmp_var_ndims))
301 
302  ! Grab the actual dimension IDs and attributes
303 
304  call nclayer_check(nf90_inquire_variable(ncid, var_index, dimids = tmp_var_dimids, &
305  xtype = tmp_var_type))
306 
307  if ((tmp_var_ndims == 2) .OR. &
308  ((tmp_var_ndims == 3) .AND. (tmp_var_type == nf90_char))) then
309  is_data2d_var = .false.
310 
311  do i = 1, tmp_var_ndims
312  call nclayer_check(nf90_inquire_dimension(ncid, tmp_var_dimids(i), tmp_var_dim_names(i), &
313  tmp_var_dim_sizes(i)))
314 
315  if (tmp_var_dim_names(i) == "nobs") then
316  nobs_size = tmp_var_dim_sizes(i)
317  if (tmp_var_type /= nf90_char) then
318  is_data2d_var = .true.
319  else if (tmp_var_type == nf90_char) then
320  if (index(tmp_var_dim_names(1), "_str_dim") /= 0) &
321  is_data2d_var = .true.
322  end if
323  end if
324  end do
325 
326  if (is_data2d_var) then
327  ! Expand things first!
329 
330  ! Add to the total!
331  diag_data2d_store%total = diag_data2d_store%total + 1
332 
333  ! Store name and type!
334  diag_data2d_store%names(diag_data2d_store%total) = trim(tmp_var_name)
335 
336  ! The relative index is the total nobs
337  rel_index = nobs_size
338 
339  if (tmp_var_type == nf90_byte) then
341  type_index = 1
342  else if (tmp_var_type == nf90_short) then
344  type_index = 2
345  else if (tmp_var_type == nf90_int) then
347  type_index = 3
348  else if (tmp_var_type == nf90_float) then
350  type_index = 4
351  else if (tmp_var_type == nf90_double) then
353  type_index = 5
354  else if (tmp_var_type == nf90_char) then
355  diag_data2d_store%max_str_lens(diag_data2d_store%total) = tmp_var_dim_sizes(1)
357  type_index = 6
358  else
359  call nclayer_error("NetCDF4 type invalid!")
360  end if
361 
362  if (tmp_var_type == nf90_char) then
363  diag_data2d_store%max_lens(diag_data2d_store%total) = tmp_var_dim_sizes(2)
364  else
365  diag_data2d_store%max_lens(diag_data2d_store%total) = tmp_var_dim_sizes(1)
366  end if
367 
368  print *, trim(tmp_var_name), "rel index", rel_index
369 
370  ! Now add a relative position... based on the next position!
371 
372  ! Set relative index!
373  diag_data2d_store%rel_indexes(diag_data2d_store%total) = rel_index
374 
375  ! Set variable ID! Note that var_index here is the actual variable ID.
376  diag_data2d_store%var_ids(diag_data2d_store%total) = var_index
377  end if
378  end if
379 
380  ! Deallocate
381  deallocate(tmp_var_dimids)
382  deallocate(tmp_var_dim_names)
383  deallocate(tmp_var_dim_sizes)
384  end do
385 
386  diag_data2d_store%def_lock = .true.
387  end subroutine nc_diag_data2d_load_def
388 
389  subroutine nc_diag_data2d_write_def(internal)
390  logical, intent(in), optional :: internal
391 
392  integer(i_byte) :: data_type
393  character(len=100) :: data2d_name
394 
395  integer(i_llong) :: curdatindex, j
396  integer(i_long) :: nc_data_type
397  integer(i_long) :: tmp_dim_id
398  character(len=120) :: data_dim_name
399  character(len=120) :: data_dim_str_name
400  integer(i_long) :: max_len
401  integer(i_long) :: max_str_len, msl_tmp
402 
403  character(len=:), allocatable :: string_arr(:)
404 
405 #ifdef ENABLE_ACTION_MSGS
406  character(len=1000) :: action_str
407 
408  if (nclayer_enable_action) then
409  if (present(internal)) then
410  write(action_str, "(A, L, A)") "nc_diag_data2d_write_def(internal = ", internal, ")"
411  else
412  write(action_str, "(A)") "nc_diag_data2d_write_def(internal = (not specified))"
413  end if
414  call nclayer_actionm(trim(action_str))
415  end if
416 #endif
417 
418  if (init_done) then
419  if (.NOT. diag_data2d_store%def_lock) then
420  ! Use global nobs ID!
421  ! Call subroutine to ensure the nobs dim is created already...
423 
424  do curdatindex = 1, diag_data2d_store%total
425  data2d_name = diag_data2d_store%names(curdatindex)
426  data_type = diag_data2d_store%types(curdatindex)
427 
428  call nclayer_info("data2d: defining " // trim(data2d_name))
429 
430  if (data_type == nlayer_byte) nc_data_type = nf90_byte
431  if (data_type == nlayer_short) nc_data_type = nf90_short
432  if (data_type == nlayer_long) nc_data_type = nf90_int
433  if (data_type == nlayer_float) nc_data_type = nf90_float
434  if (data_type == nlayer_double) nc_data_type = nf90_double
435  if (data_type == nlayer_string) nc_data_type = nf90_char
436 
437 #ifdef _DEBUG_MEM_
438  print *, "data2d part 1"
439 #endif
440 
441  ! We need to create a new dimension...
442  write (data_dim_name, "(A, A)") trim(data2d_name), "_arr_dim"
443 
444  ! Find the maximum array length of this variable!
445  max_len = nc_diag_data2d_max_len_var(curdatindex)
446 
447  ! Create this maximum array length dimension for this variable
448  if (.NOT. append_only) &
449  call nclayer_check(nf90_def_dim(ncid, data_dim_name, max_len, diag_data2d_store%var_dim_ids(curdatindex)))
450 
451  ! Store maximum length
452  diag_data2d_store%max_lens(curdatindex) = max_len;
453 
454  if (data_type == nlayer_string) then
455  max_str_len = 0
456  write (data_dim_name, "(A, A)") trim(data2d_name), "_maxstrlen"
457 
458  ! If trimming is enabled, we haven't found our max_str_len yet.
459  ! Go find it!
460  if (enable_trim) then
461  ! Dimension is # of chars by # of obs (unlimited)
462  do j = 1, diag_data2d_store%stor_i_arr(curdatindex)%icount
463  allocate(character(10000) :: string_arr(diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j)))
464  string_arr = &
465  diag_data2d_store%m_string(diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) &
466  : diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) + &
467  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j))
468 
469 #ifdef _DEBUG_MEM_
470  write(*, "(A, I0)") "DEBUG DATA2D: tmp array size is: ", diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j)
471 #endif
472 
473  ! Now we can calculate the length!
474  msl_tmp = max_len_string_array(string_arr, &
475  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j))
476 
477  if (msl_tmp > max_str_len) max_str_len = msl_tmp
478 
479 #ifdef _DEBUG_MEM_
480  write (*, "(A, A, A, I0, A, I0)") "DEBUG DATA2D DEF WRITE: at data2d_name ", trim(data2d_name), ", msl_tmp computes to ", msl_tmp, ", max_str_len computes to ", max_str_len
481  print *, "DEBUG DATA2D DEF WRITE: string array dump follows:"
482  call string_array_dump(string_arr)
483 #endif
484 
485  ! Deallocate right after we're done!
486  deallocate(string_arr)
487  end do
488 #ifdef _DEBUG_MEM_
489  write (*, "(A, A, A, I0, A, I0)") "DEBUG DATA2D DEF WRITE: ** at data2d_name ", trim(data2d_name), ", FINAL max_str_len computes to ", max_str_len, ", max_len computes to ", max_len
490 #endif
491 
492  ! Save the max string len
493  diag_data2d_store%max_str_lens(curdatindex) = max_str_len
494  end if
495 
496  ! Create dimension needed!
497  write (data_dim_str_name, "(A, A)") trim(data2d_name), "_str_dim"
498  if (.NOT. append_only) &
499  call nclayer_check(nf90_def_dim(ncid, data_dim_str_name, &
500  diag_data2d_store%max_str_lens(curdatindex), tmp_dim_id))
501 
502 #ifdef _DEBUG_MEM_
503  print *, "Defining char var type..."
504 #endif
505 
506  if (.NOT. append_only) &
507  call nclayer_check(nf90_def_var(ncid, data2d_name, nc_data_type, &
508  (/ tmp_dim_id, diag_data2d_store%var_dim_ids(curdatindex), &
509  diag_varattr_store%nobs_dim_id /), &
510  diag_data2d_store%var_ids(curdatindex)))
511 
512 #ifdef _DEBUG_MEM_
513  write (*, "(A, A, A, I0, A, I0)") "DEBUG DATA2D DEF WRITE: ** at data2d_name ", trim(data2d_name), ", result VID is ", diag_data2d_store%var_ids(curdatindex)
514  write (*, "(A, I0, A, I0)") "DEBUG DATA2D DEF WRITE: ** result dim is unlim x max_len = ", max_len, " x max_str_len = ", diag_data2d_store%max_str_lens(curdatindex)
515  print *, "data2d part 2"
516 #endif
517 
518 #ifdef _DEBUG_MEM_
519  print *, "Done defining char var type..."
520 #endif
521  else
522 #ifdef _DEBUG_MEM_
523  print *, "Definition for variable " // trim(data2d_name) // ":"
524  print *, diag_data2d_store%max_lens(curdatindex), "x unlimited (NetCDF order)"
525 #endif
526  if (.NOT. append_only) &
527  call nclayer_check(nf90_def_var(ncid, data2d_name, nc_data_type, &
528  (/ diag_data2d_store%var_dim_ids(curdatindex), diag_varattr_store%nobs_dim_id /), &
529  diag_data2d_store%var_ids(curdatindex)))
530  end if
531 
532  call nc_diag_varattr_add_var(diag_data2d_store%names(curdatindex), &
533  diag_data2d_store%types(curdatindex), &
534  diag_data2d_store%var_ids(curdatindex))
535 
536  ! Enable compression
537  ! Args: ncid, varid, enable_shuffle (yes), enable_deflate (yes), deflate_level
538 #ifdef _DEBUG_MEM_
539  print *, "Defining compression 1 (chunking)..."
540 #endif
541 
542  if (.NOT. append_only) then
543  if (data_type == nlayer_string) then
544  call nclayer_check(nf90_def_var_chunking(ncid, diag_data2d_store%var_ids(curdatindex), &
545  nf90_chunked, (/ diag_data2d_store%max_str_lens(curdatindex), &
546  diag_data2d_store%max_lens(curdatindex), nlayer_chunking /)))
547  else
548  call nclayer_check(nf90_def_var_chunking(ncid, diag_data2d_store%var_ids(curdatindex), &
549  nf90_chunked, (/ diag_data2d_store%max_lens(curdatindex), nlayer_chunking /)))
550  end if
551  end if
552 #ifdef _DEBUG_MEM_
553  print *, "Defining compression 2 (gzip)..."
554 #endif
555  if (.NOT. append_only) &
556  call nclayer_check(nf90_def_var_deflate(ncid, diag_data2d_store%var_ids(curdatindex), &
557  1, 1, nlayer_compression))
558 
559 #ifdef _DEBUG_MEM_
560  print *, "Done defining compression..."
561 #endif
562 
563  ! Lock the definitions!
564  diag_data2d_store%def_lock = .true.
565  end do
566  else
567  if(.NOT. present(internal)) &
568  call nclayer_error("Can't write definitions - definitions have already been written and locked!")
569  end if
570  end if
571  end subroutine nc_diag_data2d_write_def
572 
573  subroutine nc_diag_data2d_write_data(flush_data_only)
574  ! Optional internal flag to only flush data - if this is
575  ! true, data flushing will be performed, and the data will
576  ! NOT be locked.
577  logical, intent(in), optional :: flush_data_only
578 
579  integer(i_byte) :: data_type
580  character(len=100) :: data2d_name
581 
582  ! For some strange reason, curdatindex needs to be
583  ! initialized here to 1, otherwise a runtime error of using
584  ! an undefined variable occurs... even though it's set
585  ! by the DO loop...
586  integer(i_long) :: curdatindex = 1, j
587 #ifdef _DEBUG_MEM_
588  ! Index counter for inner loop (intermediate) array debug
589  integer(i_long) :: i
590 #endif
591 
592  integer(i_byte), dimension(:, :), allocatable :: byte_arr
593  integer(i_short),dimension(:, :), allocatable :: short_arr
594  integer(i_long), dimension(:, :), allocatable :: long_arr
595  real(r_single), dimension(:, :), allocatable :: rsingle_arr
596  real(r_double), dimension(:, :), allocatable :: rdouble_arr
597  character(len=:),dimension(:, :), allocatable :: string_arr
598 
599  integer(i_long) :: max_str_len
600 
601  integer(i_llong) :: data_length_counter
602  character(len=100) :: counter_data_name
603  integer(i_llong) :: current_length_count
604  character(len=1000) :: data_uneven_msg
605 
606 #ifdef ENABLE_ACTION_MSGS
607  character(len=1000) :: action_str
608 
609  if (nclayer_enable_action) then
610  if (present(flush_data_only)) then
611  write(action_str, "(A, L, A)") "nc_diag_data2d_write_data(flush_data_only = ", flush_data_only, ")"
612  else
613  write(action_str, "(A)") "nc_diag_data2d_write_data(flush_data_only = (not specified))"
614  end if
615  call nclayer_actionm(trim(action_str))
616  end if
617 #endif
618 
619  ! Initialization MUST occur here, not in decl...
620  ! Otherwise, it'll initialize once, and never again...
621  !
622  ! This will cause scary issues in the future, where closing
623  ! and opening a new file shows strange errors about a file
624  ! opened in the past...
625  max_str_len = -1
626  data_length_counter = -1
627  current_length_count = -1
628 
629  if (init_done .AND. allocated(diag_data2d_store)) then
630  if (.NOT. diag_data2d_store%data_lock) then
631  do curdatindex = 1, diag_data2d_store%total
632 #ifdef _DEBUG_MEM_
633  print *, curdatindex
634 #endif
635  data2d_name = diag_data2d_store%names(curdatindex)
636  data_type = diag_data2d_store%types(curdatindex)
637 
638  call nclayer_info("data2d: writing " // trim(data2d_name))
639 
640  ! Warn about data inconsistencies
641  if (.NOT. (present(flush_data_only) .AND. flush_data_only)) then
642  current_length_count = diag_data2d_store%stor_i_arr(curdatindex)%icount + &
643  diag_data2d_store%rel_indexes(curdatindex)
644 
645  if (data_length_counter == -1) then
646  data_length_counter = current_length_count
647  counter_data_name = data2d_name
648  else
649  if (data_length_counter /= current_length_count) then
650  ! Show message!
651  ! NOTE - I0 and TRIM are Fortran 95 specs
652  write (data_uneven_msg, "(A, I0, A, I0, A)") "Amount of data written in " // &
653  trim(data2d_name) // " (", &
654  current_length_count, &
655  ")" // char(10) // &
656  " differs from variable " // trim(counter_data_name) // &
657  " (", data_length_counter, ")!"
658 
659  if (diag_data2d_store%strict_check) then
660  call nclayer_error(trim(data_uneven_msg))
661  else
662  call nclayer_warning(trim(data_uneven_msg))
663  end if
664  end if
665  end if
666  end if
667 
668  ! Make sure we have data to write in the first place!
669  if (diag_data2d_store%stor_i_arr(curdatindex)%icount > 0) then
670  ! MAJOR GOTCHA:
671  ! Fortran is weird... and by weird, we mean Fortran's indexing
672  ! system! Fortran uses a column-major system, which means that
673  ! we MUST allocate and store in a column-major format! Each
674  ! column needs to store a single array of data. Before, with
675  ! single dimensions, this didn't matter since the data itself
676  ! was automatically stored into a column. With 2D data,
677  ! we MUST be aware of the reversed dimensions!
678  ! (NetCDF4 respects the Fortran way, and takes in a "row" of
679  ! data via columns!)
680 
681  if (data_type == nlayer_byte) then
682  allocate(byte_arr(diag_data2d_store%max_lens(curdatindex), &
683  diag_data2d_store%stor_i_arr(curdatindex)%icount))
684 
685  byte_arr = nlayer_fill_byte
686 
687  do j = 1, diag_data2d_store%stor_i_arr(curdatindex)%icount
688  ! Just in case our definition checks failed...
689  if (diag_data2d_store%max_lens(curdatindex) /= &
690  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j)) then
691  ! Show message!
692  ! NOTE - I0 and TRIM are Fortran 95 specs
693  write (data_uneven_msg, "(A, I0, A, I0, A)") "Amount of data written in " // &
694  trim(data2d_name) // " (", &
695  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), &
696  ")" // char(10) // &
697  " does not match the variable length" // &
698  " (", diag_data2d_store%max_lens(curdatindex), ")!"
699 
700  if (diag_data2d_store%strict_check) then
701  call nclayer_error(trim(data_uneven_msg))
702  else
703  call nclayer_warning(trim(data_uneven_msg))
704  end if
705  end if
706 
707  byte_arr(1 : diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), j) = &
708  diag_data2d_store%m_byte( &
709  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) : &
710  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) + &
711  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j) - 1)
712  end do
713 
714  call nclayer_check(nf90_put_var(&
715  ncid, diag_data2d_store%var_ids(curdatindex), &
716  byte_arr, &
717  (/ 1, 1 + diag_data2d_store%rel_indexes(curdatindex) /), &
718  (/ diag_data2d_store%max_lens(curdatindex), &
719  diag_data2d_store%stor_i_arr(curdatindex)%icount /) &
720  ))
721 
722  deallocate(byte_arr)
723  else if (data_type == nlayer_short) then
724  allocate(short_arr(diag_data2d_store%max_lens(curdatindex), &
725  diag_data2d_store%stor_i_arr(curdatindex)%icount))
726 
727  short_arr = nlayer_fill_short
728 
729  do j = 1, diag_data2d_store%stor_i_arr(curdatindex)%icount
730  ! Just in case our definition checks failed...
731  if (diag_data2d_store%max_lens(curdatindex) /= &
732  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j)) then
733  ! Show message!
734  ! NOTE - I0 and TRIM are Fortran 95 specs
735  write (data_uneven_msg, "(A, I0, A, I0, A)") "Amount of data written in " // &
736  trim(data2d_name) // " (", &
737  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), &
738  ")" // char(10) // &
739  " does not match the variable length" // &
740  " (", diag_data2d_store%max_lens(curdatindex), ")!"
741 
742  if (diag_data2d_store%strict_check) then
743  call nclayer_error(trim(data_uneven_msg))
744  else
745  call nclayer_warning(trim(data_uneven_msg))
746  end if
747  end if
748 
749  short_arr(1 : diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), j) = &
750  diag_data2d_store%m_short( &
751  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) : &
752  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) + &
753  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j) - 1)
754  end do
755 
756  call nclayer_check(nf90_put_var(&
757  ncid, diag_data2d_store%var_ids(curdatindex), &
758  short_arr, &
759  (/ 1, 1 + diag_data2d_store%rel_indexes(curdatindex) /), &
760  (/ diag_data2d_store%max_lens(curdatindex), &
761  diag_data2d_store%stor_i_arr(curdatindex)%icount /) &
762  ))
763 
764  deallocate(short_arr)
765  else if (data_type == nlayer_long) then
766  !allocate(long_arr(diag_data2d_store%stor_i_arr(curdatindex)%icount, &
767  ! diag_data2d_store%max_lens(curdatindex)))
768 
769  allocate(long_arr(diag_data2d_store%max_lens(curdatindex), &
770  diag_data2d_store%stor_i_arr(curdatindex)%icount))
771 
772 #ifdef _DEBUG_MEM_
773  write (*, "(A, I0)") "NLAYER_FILL_LONG = ", nlayer_fill_long
774 #endif
775 
776  long_arr = nlayer_fill_long
777 
778 #ifdef _DEBUG_MEM_
779  write (*, "(A)") "************ DEBUG: INITIAL var array for " // trim(data2d_name)
780  do j = 1, diag_data2d_store%stor_i_arr(curdatindex)%icount
781  print *, long_arr(:, j)
782  end do
783 #endif
784 
785  do j = 1, diag_data2d_store%stor_i_arr(curdatindex)%icount
786  ! Just in case our definition checks failed...
787  if (diag_data2d_store%max_lens(curdatindex) /= &
788  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j)) then
789  ! Show message!
790  ! NOTE - I0 and TRIM are Fortran 95 specs
791  write (data_uneven_msg, "(A, I0, A, I0, A)") "Amount of data written in " // &
792  trim(data2d_name) // " (", &
793  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), &
794  ")" // char(10) // &
795  " does not match the variable length" // &
796  " (", diag_data2d_store%max_lens(curdatindex), ")!"
797 
798  if (diag_data2d_store%strict_check) then
799  call nclayer_error(trim(data_uneven_msg))
800  else
801  call nclayer_warning(trim(data_uneven_msg))
802  end if
803  end if
804 
805 #ifdef _DEBUG_MEM_
806  write (*, "(A, I0, A)") "Adding to long_arr, index ", j, ":"
807  print *, diag_data2d_store%m_long( &
808  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) : &
809  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) + &
810  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j) - 1)
811  write (*, "(A, I0)") " -> length of subarr: ", diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j)
812 #endif
813 
814  long_arr(1 : diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), j) = &
815  diag_data2d_store%m_long( &
816  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) : &
817  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) + &
818  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j) - 1)
819 
820 #ifdef _DEBUG_MEM_
821  write (*, "(A)") "************ DEBUG: INTERMEDIATE var array for " // trim(data2d_name)
822  do i = 1, diag_data2d_store%stor_i_arr(curdatindex)%icount
823  print *, long_arr(:, i)
824  end do
825 #endif
826  end do
827 
828 #ifdef _DEBUG_MEM_
829  write (*, "(A, I0, A, I0, A, I0, A, I0, A)") &
830  "Writing long with start = (", 1, ", ", &
831  1 + diag_data2d_store%rel_indexes(curdatindex), &
832  "), count = (", diag_data2d_store%stor_i_arr(curdatindex)%icount, &
833  ", ", 1, ")"
834 
835  write (*, "(A, I0, A, I0)") "************ DEBUG: dim for " // trim(data2d_name) // ": ", &
836  diag_data2d_store%stor_i_arr(curdatindex)%icount, " by ", &
837  diag_data2d_store%max_lens(curdatindex)
838  write (*, "(A)") "************ DEBUG: var array for " // trim(data2d_name)
839 
840  do j = 1, diag_data2d_store%stor_i_arr(curdatindex)%icount
841  print *, long_arr(:, j)
842  end do
843 #endif
844 
845  call nclayer_check(nf90_put_var(&
846  ncid, diag_data2d_store%var_ids(curdatindex), &
847  long_arr, &
848  (/ 1, 1 + diag_data2d_store%rel_indexes(curdatindex) /), &
849  (/ diag_data2d_store%max_lens(curdatindex), &
850  diag_data2d_store%stor_i_arr(curdatindex)%icount /) &
851  ))
852 
853  deallocate(long_arr)
854  else if (data_type == nlayer_float) then
855  allocate(rsingle_arr(diag_data2d_store%max_lens(curdatindex), &
856  diag_data2d_store%stor_i_arr(curdatindex)%icount))
857 
858  rsingle_arr = nlayer_fill_float
859 
860  do j = 1, diag_data2d_store%stor_i_arr(curdatindex)%icount
861  ! Just in case our definition checks failed...
862  if (diag_data2d_store%max_lens(curdatindex) /= &
863  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j)) then
864  ! Show message!
865  ! NOTE - I0 and TRIM are Fortran 95 specs
866  write (data_uneven_msg, "(A, I0, A, I0, A)") "Amount of data written in " // &
867  trim(data2d_name) // " (", &
868  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), &
869  ")" // char(10) // &
870  " does not match the variable length" // &
871  " (", diag_data2d_store%max_lens(curdatindex), ")!"
872 
873  if (diag_data2d_store%strict_check) then
874  call nclayer_error(trim(data_uneven_msg))
875  else
876  call nclayer_warning(trim(data_uneven_msg))
877  end if
878  end if
879 
880  rsingle_arr(1 : diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), j) = &
881  diag_data2d_store%m_rsingle( &
882  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) : &
883  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) + &
884  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j) - 1)
885  end do
886 
887  !print *, "end queue / start put"
888  call nclayer_check(nf90_put_var(&
889  ncid, diag_data2d_store%var_ids(curdatindex), &
890  rsingle_arr, &
891  (/ 1, 1 + diag_data2d_store%rel_indexes(curdatindex) /), &
892  (/ diag_data2d_store%max_lens(curdatindex), &
893  diag_data2d_store%stor_i_arr(curdatindex)%icount /) &
894  ))
895  !call nclayer_check(nf90_sync(ncid))
896  deallocate(rsingle_arr)
897  !print *, "end put"
898 
899  else if (data_type == nlayer_double) then
900  allocate(rdouble_arr(diag_data2d_store%max_lens(curdatindex), &
901  diag_data2d_store%stor_i_arr(curdatindex)%icount))
902 
903  rdouble_arr = nlayer_fill_double
904 
905  do j = 1, diag_data2d_store%stor_i_arr(curdatindex)%icount
906  ! Just in case our definition checks failed...
907  if (diag_data2d_store%max_lens(curdatindex) /= &
908  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j)) then
909  ! Show message!
910  ! NOTE - I0 and TRIM are Fortran 95 specs
911  write (data_uneven_msg, "(A, I0, A, I0, A)") "Amount of data written in " // &
912  trim(data2d_name) // " (", &
913  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), &
914  ")" // char(10) // &
915  " does not match the variable length" // &
916  " (", diag_data2d_store%max_lens(curdatindex), ")!"
917 
918  if (diag_data2d_store%strict_check) then
919  call nclayer_error(trim(data_uneven_msg))
920  else
921  call nclayer_warning(trim(data_uneven_msg))
922  end if
923  end if
924 
925  rdouble_arr(1 : diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), j) = &
926  diag_data2d_store%m_rdouble( &
927  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) : &
928  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) + &
929  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j) - 1)
930  end do
931 
932  call nclayer_check(nf90_put_var(&
933  ncid, diag_data2d_store%var_ids(curdatindex), &
934  rdouble_arr, &
935  (/ 1, 1 + diag_data2d_store%rel_indexes(curdatindex) /), &
936  (/ diag_data2d_store%max_lens(curdatindex), &
937  diag_data2d_store%stor_i_arr(curdatindex)%icount /) &
938  ))
939  deallocate(rdouble_arr)
940  else if (data_type == nlayer_string) then
941  ! We need to seperate everything because the Intel Fortran compiler loves
942  ! to optimize... and then assume that I'll try to use an unallocated variable,
943  ! even with checks.
944  if (allocated(diag_data2d_store%max_str_lens)) then
945  max_str_len = diag_data2d_store%max_str_lens(curdatindex)
946  else
947  call nclayer_error("BUG: diag_data2d_store%max_str_lens not allocated yet!")
948  end if
949 
950  allocate(character(max_str_len) :: &
951  string_arr(diag_data2d_store%max_lens(curdatindex), &
952  diag_data2d_store%stor_i_arr(curdatindex)%icount &
953  ))
954 
955  string_arr = nlayer_fill_char
956 
957  do j = 1, diag_data2d_store%stor_i_arr(curdatindex)%icount
958  ! Just in case our definition checks failed...
959  if (diag_data2d_store%max_lens(curdatindex) /= &
960  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j)) then
961  ! Show message!
962  ! NOTE - I0 and TRIM are Fortran 95 specs
963  write (data_uneven_msg, "(A, I0, A, I0, A)") "Amount of data written in " // &
964  trim(data2d_name) // " (", &
965  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), &
966  ")" // char(10) // &
967  " does not match the variable length" // &
968  " (", diag_data2d_store%max_lens(curdatindex), ")!"
969 
970  if (diag_data2d_store%strict_check) then
971  call nclayer_error(trim(data_uneven_msg))
972  else
973  call nclayer_warning(trim(data_uneven_msg))
974  end if
975  end if
976 
977  string_arr(1 : diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j), j) = &
978  diag_data2d_store%m_string( &
979  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) : &
980  diag_data2d_store%stor_i_arr(curdatindex)%index_arr(j) + &
981  diag_data2d_store%stor_i_arr(curdatindex)%length_arr(j) - 1)
982  end do
983 
984  if (allocated(diag_data2d_store%max_str_lens)) then
985  call nclayer_check(nf90_put_var(&
986  ncid, diag_data2d_store%var_ids(curdatindex), &
987  string_arr, &
988  (/ 1, 1, 1 + diag_data2d_store%rel_indexes(curdatindex) /), &
989  (/ diag_data2d_store%max_str_lens(curdatindex), &
990  diag_data2d_store%max_lens(curdatindex), &
991  diag_data2d_store%stor_i_arr(curdatindex)%icount /) &
992  ))
993  else
994  call nclayer_error("BUG: diag_data2d_store%max_str_lens not allocated yet!")
995  end if
996 
997  deallocate(string_arr)
998  end if
999 
1000  ! Check for data flushing, and if so, update the relative indexes
1001  ! and set icount to 0.
1002  if (present(flush_data_only) .AND. flush_data_only) then
1003  diag_data2d_store%rel_indexes(curdatindex) = &
1004  diag_data2d_store%rel_indexes(curdatindex) + &
1005  diag_data2d_store%stor_i_arr(curdatindex)%icount
1006  diag_data2d_store%stor_i_arr(curdatindex)%icount = 0
1007 
1008 #ifdef _DEBUG_MEM_
1009  print *, "diag_data2d_store%rel_indexes(curdatindex) is now:"
1010  print *, diag_data2d_store%rel_indexes(curdatindex)
1011 #endif
1012  end if
1013 
1014  end if
1015  end do
1016 
1017  if (present(flush_data_only) .AND. flush_data_only) then
1018 #ifdef _DEBUG_MEM_
1019  print *, "In buffer flush mode!"
1020 #endif
1021 
1022  ! We need to reset all array counts to zero!
1023  diag_data2d_store%acount = 0
1024  else
1025  ! Lock data writing
1026  diag_data2d_store%data_lock = .true.
1027 #ifdef _DEBUG_MEM_
1028  print *, "In data lock mode!"
1029 #endif
1030  end if
1031  else
1032  call nclayer_error("Can't write data - data have already been written and locked!")
1033  end if
1034  else
1035  call nclayer_error("Can't write data - NetCDF4 layer not initialized yet!")
1036  end if
1037 
1038 #ifdef _DEBUG_MEM_
1039  print *, "All done writing data2d data"
1040 #endif
1041  end subroutine nc_diag_data2d_write_data
1042 
1043  ! Set strict checking
1044  subroutine nc_diag_data2d_set_strict(enable_strict)
1045  logical, intent(in) :: enable_strict
1046 
1047  if (init_done .AND. allocated(diag_data2d_store)) then
1048  diag_data2d_store%strict_check = enable_strict
1049  else
1050  call nclayer_error("Can't set strictness level for data2d - NetCDF4 layer not initialized yet!")
1051  end if
1052  end subroutine nc_diag_data2d_set_strict
1053 
1054  ! Preallocate variable name/type/etc. storage.
1055  subroutine nc_diag_data2d_prealloc_vars(num_of_addl_vars)
1056  integer(i_llong), intent(in) :: num_of_addl_vars
1057 #ifdef ENABLE_ACTION_MSGS
1058  character(len=1000) :: action_str
1059 
1060  if (nclayer_enable_action) then
1061  write(action_str, "(A, I0, A)") "nc_diag_data2d_prealloc_vars(num_of_addl_vars = ", num_of_addl_vars, ")"
1062  call nclayer_actionm(trim(action_str))
1063  end if
1064 #endif
1065  if (init_done .AND. allocated(diag_data2d_store)) then
1066  if (allocated(diag_data2d_store%names)) then
1067  if (diag_data2d_store%total >= size(diag_data2d_store%names)) then
1068  call nc_diag_realloc(diag_data2d_store%names, num_of_addl_vars)
1069  end if
1070  else
1071  allocate(diag_data2d_store%names(nlayer_default_ent + num_of_addl_vars))
1072  end if
1073 
1074  if (allocated(diag_data2d_store%types)) then
1075  if (diag_data2d_store%total >= size(diag_data2d_store%types)) then
1076  call nc_diag_realloc(diag_data2d_store%types, num_of_addl_vars)
1077  end if
1078  else
1079  allocate(diag_data2d_store%types(nlayer_default_ent + num_of_addl_vars))
1080  end if
1081 
1082  if (allocated(diag_data2d_store%stor_i_arr)) then
1083  if (diag_data2d_store%total >= size(diag_data2d_store%stor_i_arr)) then
1084  call nc_diag_data2d_resize_iarr_type(num_of_addl_vars)
1085  end if
1086  else
1087  allocate(diag_data2d_store%stor_i_arr(nlayer_default_ent + num_of_addl_vars))
1088  end if
1089 
1090  if (allocated(diag_data2d_store%var_ids)) then
1091  if (diag_data2d_store%total >= size(diag_data2d_store%var_ids)) then
1092  call nc_diag_realloc(diag_data2d_store%var_ids, num_of_addl_vars)
1093  end if
1094  else
1095  allocate(diag_data2d_store%var_ids(nlayer_default_ent + num_of_addl_vars))
1096  diag_data2d_store%var_ids = -1
1097  end if
1098 
1099  if (allocated(diag_data2d_store%var_dim_ids)) then
1100  if (diag_data2d_store%total >= size(diag_data2d_store%var_dim_ids)) then
1101  call nc_diag_realloc(diag_data2d_store%var_dim_ids, num_of_addl_vars)
1102  end if
1103  else
1104  allocate(diag_data2d_store%var_dim_ids(nlayer_default_ent + num_of_addl_vars))
1105  diag_data2d_store%var_dim_ids = -1
1106  end if
1107 
1108  if (allocated(diag_data2d_store%alloc_sia_multi)) then
1109  if (diag_data2d_store%total >= size(diag_data2d_store%alloc_sia_multi)) then
1110  call nc_diag_realloc(diag_data2d_store%alloc_sia_multi, num_of_addl_vars)
1111  end if
1112  else
1113  allocate(diag_data2d_store%alloc_sia_multi(nlayer_default_ent + num_of_addl_vars))
1114  diag_data2d_store%alloc_sia_multi = 0
1115  end if
1116 
1117  if (allocated(diag_data2d_store%max_str_lens)) then
1118  if (diag_data2d_store%total >= size(diag_data2d_store%max_str_lens)) then
1119  call nc_diag_realloc(diag_data2d_store%max_str_lens, num_of_addl_vars)
1120  end if
1121  else
1122  allocate(diag_data2d_store%max_str_lens(nlayer_default_ent + num_of_addl_vars))
1123  diag_data2d_store%max_str_lens = -1
1124  end if
1125 
1126  if (allocated(diag_data2d_store%rel_indexes)) then
1127  if (diag_data2d_store%total >= size(diag_data2d_store%rel_indexes)) then
1128  call nc_diag_realloc(diag_data2d_store%rel_indexes, num_of_addl_vars)
1129  end if
1130  else
1131  allocate(diag_data2d_store%rel_indexes(nlayer_default_ent + num_of_addl_vars))
1132  diag_data2d_store%rel_indexes = 0
1133  end if
1134 
1135  if (allocated(diag_data2d_store%max_lens)) then
1136  if (diag_data2d_store%total >= size(diag_data2d_store%max_lens)) then
1137  call nc_diag_realloc(diag_data2d_store%max_lens, num_of_addl_vars)
1138  end if
1139  else
1140  allocate(diag_data2d_store%max_lens(nlayer_default_ent + num_of_addl_vars))
1141  diag_data2d_store%max_lens = 0
1142  end if
1143 
1144  diag_data2d_store%prealloc_total = diag_data2d_store%prealloc_total + num_of_addl_vars
1145  else
1146  call nclayer_error("NetCDF4 layer not initialized yet!")
1147  endif
1148  end subroutine nc_diag_data2d_prealloc_vars
1149 
1150  ! Preallocate actual variable data storage
1151  subroutine nc_diag_data2d_prealloc_vars_storage(nclayer_type, num_of_addl_slots)
1152  integer(i_byte), intent(in) :: nclayer_type
1153  integer(i_llong), intent(in) :: num_of_addl_slots
1154 
1155 #ifdef ENABLE_ACTION_MSGS
1156  character(len=1000) :: action_str
1157 
1158  if (nclayer_enable_action) then
1159  write(action_str, "(A, I0, A, I0, A)") "nc_diag_data2d_prealloc_vars_storage(nclayer_type = ", nclayer_type, ", num_of_addl_slots = ", num_of_addl_slots, ")"
1160  call nclayer_actionm(trim(action_str))
1161  end if
1162 #endif
1163 
1164  if (nclayer_type == nlayer_byte) then
1165  call nc_diag_data2d_resize_byte(num_of_addl_slots, .false.)
1166  else if (nclayer_type == nlayer_short) then
1167  call nc_diag_data2d_resize_short(num_of_addl_slots, .false.)
1168  else if (nclayer_type == nlayer_long) then
1169  call nc_diag_data2d_resize_long(num_of_addl_slots, .false.)
1170  else if (nclayer_type == nlayer_float) then
1171  call nc_diag_data2d_resize_rsingle(num_of_addl_slots, .false.)
1172  else if (nclayer_type == nlayer_double) then
1173  call nc_diag_data2d_resize_rdouble(num_of_addl_slots, .false.)
1174  else if (nclayer_type == nlayer_string) then
1175  call nc_diag_data2d_resize_string(num_of_addl_slots, .false.)
1176  else
1177  call nclayer_error("Invalid type specified for variable storage preallocation!")
1178  end if
1180 
1181  ! Preallocate index storage
1182  subroutine nc_diag_data2d_prealloc_vars_storage_all(num_of_addl_slots)
1183  integer(i_llong), intent(in) :: num_of_addl_slots
1184  integer(i_long) :: i
1185 
1186 #ifdef ENABLE_ACTION_MSGS
1187  character(len=1000) :: action_str
1188 
1189  if (nclayer_enable_action) then
1190  write(action_str, "(A, I0, A)") "nc_diag_data2d_prealloc_vars_storage_all(num_of_addl_slots = ", num_of_addl_slots, ")"
1191  call nclayer_actionm(trim(action_str))
1192  end if
1193 #endif
1194 
1195  do i = 1, diag_data2d_store%prealloc_total
1196  call nc_diag_data2d_resize_iarr(i, num_of_addl_slots, .false.)
1197  end do
1199 
1200  subroutine nc_diag_data2d_expand
1201  integer(i_llong) :: addl_fields
1202 
1203  ! Did we realloc at all?
1204  logical :: meta_realloc
1205 
1206  meta_realloc = .false.
1207 
1208  if (init_done .AND. allocated(diag_data2d_store)) then
1209  addl_fields = 1 + (nlayer_default_ent * (nlayer_multi_base ** diag_data2d_store%alloc_s_multi))
1210 
1211 #ifdef _DEBUG_MEM_
1212  call nclayer_debug("INITIAL value of diag_data2d_store%alloc_s_multi:")
1213  print *, diag_data2d_store%alloc_s_multi
1214 #endif
1215 
1216  if (allocated(diag_data2d_store%names)) then
1217  if (diag_data2d_store%total >= size(diag_data2d_store%names)) then
1218 #ifdef _DEBUG_MEM_
1219  call nclayer_debug("Reallocating diag_data2d_store%names...")
1220  print *, (nlayer_multi_base ** diag_data2d_store%alloc_s_multi)
1221  print *, addl_fields
1222 #endif
1223  call nc_diag_realloc(diag_data2d_store%names, addl_fields)
1224 #ifdef _DEBUG_MEM_
1225  call nclayer_debug("Reallocated diag_data2d_store%names. Size:")
1226  print *, size(diag_data2d_store%names)
1227 #endif
1228  meta_realloc = .true.
1229  end if
1230  else
1231 #ifdef _DEBUG_MEM_
1232  call nclayer_debug("Allocating diag_data2d_store%names for first time...")
1233  print *, nlayer_default_ent
1234 #endif
1235 
1236  allocate(diag_data2d_store%names(nlayer_default_ent))
1237 
1238 #ifdef _DEBUG_MEM_
1239  call nclayer_debug("Allocated diag_data2d_store%names. Size:")
1240  print *, size(diag_data2d_store%names)
1241 #endif
1242  end if
1243 
1244  if (allocated(diag_data2d_store%types)) then
1245  if (diag_data2d_store%total >= size(diag_data2d_store%types)) then
1246 #ifdef _DEBUG_MEM_
1247  call nclayer_debug("Reallocating diag_data2d_store%types...")
1248  print *, (nlayer_multi_base ** diag_data2d_store%alloc_s_multi)
1249  print *, addl_fields
1250 #endif
1251  call nc_diag_realloc(diag_data2d_store%types, addl_fields)
1252  meta_realloc = .true.
1253  end if
1254  else
1255  allocate(diag_data2d_store%types(nlayer_default_ent))
1256  end if
1257 
1258  if (allocated(diag_data2d_store%stor_i_arr)) then
1259  if (diag_data2d_store%total >= size(diag_data2d_store%stor_i_arr)) then
1260 #ifdef _DEBUG_MEM_
1261  call nclayer_debug("Reallocating diag_data2d_store%stor_i_arr...")
1262  print *, (nlayer_multi_base ** diag_data2d_store%alloc_s_multi)
1263  print *, (1 + (nlayer_default_ent * (nlayer_multi_base ** diag_data2d_store%alloc_s_multi)))
1264 #endif
1265  call nc_diag_data2d_resize_iarr_type(addl_fields)
1266 
1267  meta_realloc = .true.
1268  end if
1269  else
1270  allocate(diag_data2d_store%stor_i_arr(nlayer_default_ent))
1271  end if
1272 
1273  if (allocated(diag_data2d_store%var_ids)) then
1274  if (diag_data2d_store%total >= size(diag_data2d_store%var_ids)) then
1275  call nc_diag_realloc(diag_data2d_store%var_ids, addl_fields)
1276  meta_realloc = .true.
1277  end if
1278  else
1279  allocate(diag_data2d_store%var_ids(nlayer_default_ent))
1280  diag_data2d_store%var_ids = -1
1281  end if
1282 
1283  if (allocated(diag_data2d_store%var_dim_ids)) then
1284  if (diag_data2d_store%total >= size(diag_data2d_store%var_dim_ids)) then
1285  call nc_diag_realloc(diag_data2d_store%var_dim_ids, addl_fields)
1286  meta_realloc = .true.
1287  end if
1288  else
1289  allocate(diag_data2d_store%var_dim_ids(nlayer_default_ent))
1290  diag_data2d_store%var_dim_ids = -1
1291  end if
1292 
1293  if (allocated(diag_data2d_store%alloc_sia_multi)) then
1294  if (diag_data2d_store%total >= size(diag_data2d_store%alloc_sia_multi)) then
1295  call nc_diag_realloc(diag_data2d_store%alloc_sia_multi, addl_fields)
1296  meta_realloc = .true.
1297  end if
1298  else
1299  allocate(diag_data2d_store%alloc_sia_multi(nlayer_default_ent))
1300  diag_data2d_store%alloc_sia_multi = 0
1301  end if
1302 
1303  if (allocated(diag_data2d_store%max_str_lens)) then
1304  if (diag_data2d_store%total >= size(diag_data2d_store%max_str_lens)) then
1305  call nc_diag_realloc(diag_data2d_store%max_str_lens, addl_fields)
1306  meta_realloc = .true.
1307  end if
1308  else
1309  allocate(diag_data2d_store%max_str_lens(nlayer_default_ent))
1310  diag_data2d_store%max_str_lens = -1
1311  end if
1312 
1313  if (allocated(diag_data2d_store%rel_indexes)) then
1314  if (diag_data2d_store%total >= size(diag_data2d_store%rel_indexes)) then
1315  call nc_diag_realloc(diag_data2d_store%rel_indexes, addl_fields)
1316  meta_realloc = .true.
1317  end if
1318  else
1319  allocate(diag_data2d_store%rel_indexes(nlayer_default_ent))
1320  diag_data2d_store%rel_indexes = 0
1321  end if
1322 
1323  if (allocated(diag_data2d_store%max_lens)) then
1324  if (diag_data2d_store%total >= size(diag_data2d_store%max_lens)) then
1325  call nc_diag_realloc(diag_data2d_store%max_lens, addl_fields)
1326  meta_realloc = .true.
1327  end if
1328  else
1329  allocate(diag_data2d_store%max_lens(nlayer_default_ent))
1330  diag_data2d_store%max_lens = 0
1331  end if
1332 
1333  if (meta_realloc) then
1334  diag_data2d_store%alloc_s_multi = diag_data2d_store%alloc_s_multi + 1
1335 #ifdef _DEBUG_MEM_
1336  print *, "Incrementing alloc_s_multi... new value:"
1337  print *, diag_data2d_store%alloc_s_multi
1338 #endif
1339  endif
1340  else
1341  call nclayer_error("NetCDF4 layer not initialized yet!")
1342  endif
1343 
1344  end subroutine nc_diag_data2d_expand
1345 
1346  function nc_diag_data2d_lookup_var(data2d_name) result(ind)
1347  character(len=*), intent(in) :: data2d_name
1348  integer :: i, ind
1349 
1350  ind = -1
1351 
1352  if (init_done .AND. allocated(diag_data2d_store)) then
1353  do i = 1, diag_data2d_store%total
1354  if (diag_data2d_store%names(i) == data2d_name) then
1355  ind = i
1356  exit
1357  end if
1358  end do
1359  end if
1360  end function nc_diag_data2d_lookup_var
1361 
1362  ! nc_diag_data2d - input integer(i_byte)
1363  ! Corresponding NetCDF4 type: byte
1364  subroutine nc_diag_data2d_byte(data2d_name, data2d_value)
1365  character(len=*), intent(in) :: data2d_name
1366  integer(i_byte), dimension(:), intent(in) :: data2d_value
1367 
1368  integer(i_long) :: var_index
1369  integer(i_llong) :: input_size
1370 
1371 #ifdef ENABLE_ACTION_MSGS
1372  character(len=1000) :: action_str
1373  integer(i_llong) :: data_value_size
1374 
1375  if (nclayer_enable_action) then
1376  data_value_size = size(data2d_value)
1377  write(action_str, "(A, I0, A, I0, A, I0, A, I0, A)") &
1378  "nc_diag_data2d_byte(data2d_name = " // data2d_name // &
1379  ", data2d_value = array with length of ", &
1380  data_value_size, &
1381  " [", &
1382  data2d_value(1), &
1383  " ... ", &
1384  data2d_value(data_value_size), &
1385  "]"
1386  call nclayer_actionm(trim(action_str))
1387  end if
1388 #endif
1389 
1390  if (diag_data2d_store%data_lock) then
1391  call nclayer_error("Can't add new data - data have already been written and locked!")
1392  end if
1393 
1394  var_index = nc_diag_data2d_lookup_var(data2d_name)
1395 
1396  if (var_index == -1) then
1397  ! First, check to make sure we can still define new variables.
1398  if (diag_data2d_store%def_lock) then
1399  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1400  end if
1401 
1403 
1404  diag_data2d_store%total = diag_data2d_store%total + 1
1405 
1406  diag_data2d_store%names(diag_data2d_store%total) = data2d_name
1408 
1409  var_index = diag_data2d_store%total
1410  end if
1411 
1412  ! Get input size and do size checks!
1413  input_size = size(data2d_value)
1414 
1415  if ((diag_data2d_store%def_lock) .AND. &
1416  (size(data2d_value) > diag_data2d_store%max_lens(var_index))) then
1417  call nclayer_error("Cannot expand variable size after locking variable definitions!")
1418  end if
1419 
1420  ! We just need to add one entry...
1421  call nc_diag_data2d_resize_iarr(var_index, 1_i_llong)
1422  call nc_diag_data2d_resize_byte(input_size)
1423 
1424  ! Now add the actual entry!
1425  diag_data2d_store%m_byte(diag_data2d_store%acount(1) - input_size + 1:diag_data2d_store%acount(1)) = &
1426  data2d_value
1427  diag_data2d_store%stor_i_arr(var_index)%index_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1428  diag_data2d_store%acount(1) - input_size + 1
1429  diag_data2d_store%stor_i_arr(var_index)%length_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1430  input_size
1431  end subroutine nc_diag_data2d_byte
1432 
1433  ! nc_diag_data2d - input integer(i_short)
1434  ! Corresponding NetCDF4 type: short
1435  subroutine nc_diag_data2d_short(data2d_name, data2d_value)
1436  character(len=*), intent(in) :: data2d_name
1437  integer(i_short), dimension(:), intent(in) :: data2d_value
1438 
1439  integer(i_long) :: var_index
1440  integer(i_llong) :: input_size
1441 
1442 #ifdef ENABLE_ACTION_MSGS
1443  character(len=1000) :: action_str
1444  integer(i_llong) :: data_value_size
1445 
1446  if (nclayer_enable_action) then
1447  data_value_size = size(data2d_value)
1448  write(action_str, "(A, I0, A, I0, A, I0, A, I0, A)") &
1449  "nc_diag_data2d_short(data2d_name = " // data2d_name // &
1450  ", data2d_value = array with length of ", &
1451  data_value_size, &
1452  " [", &
1453  data2d_value(1), &
1454  " ... ", &
1455  data2d_value(data_value_size), &
1456  "]"
1457  call nclayer_actionm(trim(action_str))
1458  end if
1459 #endif
1460 
1461  if (diag_data2d_store%data_lock) then
1462  call nclayer_error("Can't add new data - data have already been written and locked!")
1463  end if
1464 
1465  var_index = nc_diag_data2d_lookup_var(data2d_name)
1466 
1467  if (var_index == -1) then
1468  ! First, check to make sure we can still define new variables.
1469  if (diag_data2d_store%def_lock) then
1470  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1471  end if
1472 
1474 
1475  diag_data2d_store%total = diag_data2d_store%total + 1
1476 
1477  diag_data2d_store%names(diag_data2d_store%total) = data2d_name
1479 
1480  var_index = diag_data2d_store%total
1481  end if
1482 
1483  ! Get input size and do size checks!
1484  input_size = size(data2d_value)
1485 
1486  if ((diag_data2d_store%def_lock) .AND. &
1487  (size(data2d_value) > diag_data2d_store%max_lens(var_index))) then
1488  call nclayer_error("Cannot expand variable size after locking variable definitions!")
1489  end if
1490 
1491  ! We just need to add one entry...
1492  call nc_diag_data2d_resize_iarr(var_index, 1_i_llong)
1493  call nc_diag_data2d_resize_short(input_size)
1494 
1495  ! Now add the actual entry!
1496  diag_data2d_store%m_short(diag_data2d_store%acount(2) - input_size + 1:diag_data2d_store%acount(2)) = &
1497  data2d_value
1498  diag_data2d_store%stor_i_arr(var_index)%index_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1499  diag_data2d_store%acount(2) - input_size + 1
1500  diag_data2d_store%stor_i_arr(var_index)%length_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1501  input_size
1502  end subroutine nc_diag_data2d_short
1503 
1504  ! nc_diag_data2d - input integer(i_long)
1505  ! Corresponding NetCDF4 type: int (old: long)
1506  subroutine nc_diag_data2d_long(data2d_name, data2d_value)
1507  character(len=*), intent(in) :: data2d_name
1508  integer(i_long), dimension(:), intent(in) :: data2d_value
1509 
1510  integer(i_long) :: var_index
1511  integer(i_llong) :: input_size
1512 
1513 #ifdef ENABLE_ACTION_MSGS
1514  character(len=1000) :: action_str
1515  integer(i_llong) :: data_value_size
1516 
1517  if (nclayer_enable_action) then
1518  data_value_size = size(data2d_value)
1519  write(action_str, "(A, I0, A, I0, A, I0, A, I0, A)") &
1520  "nc_diag_data2d_long(data2d_name = " // data2d_name // &
1521  ", data2d_value = array with length of ", &
1522  data_value_size, &
1523  " [", &
1524  data2d_value(1), &
1525  " ... ", &
1526  data2d_value(data_value_size), &
1527  "]"
1528  call nclayer_actionm(trim(action_str))
1529  end if
1530 #endif
1531 
1532  if (diag_data2d_store%data_lock) then
1533  call nclayer_error("Can't add new data - data have already been written and locked!")
1534  end if
1535 
1536  var_index = nc_diag_data2d_lookup_var(data2d_name)
1537 
1538  if (var_index == -1) then
1539  ! First, check to make sure we can still define new variables.
1540  if (diag_data2d_store%def_lock) then
1541  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1542  end if
1543 
1545 
1546  diag_data2d_store%total = diag_data2d_store%total + 1
1547 
1548  diag_data2d_store%names(diag_data2d_store%total) = data2d_name
1550 
1551  var_index = diag_data2d_store%total
1552  end if
1553 
1554 #ifdef _DEBUG_MEM_
1555  call nclayer_debug("Current total:")
1556  print *, diag_data2d_store%total
1557 #endif
1558 
1559  ! Get input size and do size checks!
1560  input_size = size(data2d_value)
1561 
1562  if ((diag_data2d_store%def_lock) .AND. &
1563  (size(data2d_value) > diag_data2d_store%max_lens(var_index))) then
1564  call nclayer_error("Cannot expand variable size after locking variable definitions!")
1565  end if
1566 
1567  ! We just need to add one entry...
1568  call nc_diag_data2d_resize_iarr(var_index, 1_i_llong)
1569  call nc_diag_data2d_resize_long(input_size)
1570 
1571  ! Now add the actual entry!
1572  diag_data2d_store%m_long(diag_data2d_store%acount(3) - input_size + 1:diag_data2d_store%acount(3)) = &
1573  data2d_value
1574  diag_data2d_store%stor_i_arr(var_index)%index_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1575  diag_data2d_store%acount(3) - input_size + 1
1576  diag_data2d_store%stor_i_arr(var_index)%length_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1577  input_size
1578  end subroutine nc_diag_data2d_long
1579 
1580  ! nc_diag_data2d - input real(r_single)
1581  ! Corresponding NetCDF4 type: float (or real)
1582  subroutine nc_diag_data2d_rsingle(data2d_name, data2d_value)
1583  character(len=*), intent(in) :: data2d_name
1584  real(r_single), dimension(:), intent(in) :: data2d_value
1585 
1586  integer(i_long) :: var_index
1587  integer(i_llong) :: input_size
1588 
1589 #ifdef ENABLE_ACTION_MSGS
1590  character(len=1000) :: action_str
1591  integer(i_llong) :: data_value_size
1592 
1593  if (nclayer_enable_action) then
1594  data_value_size = size(data2d_value)
1595  write(action_str, "(A, I0, A, F0.5, A, F0.5, A)") &
1596  "nc_diag_data2d_rsingle(data2d_name = " // data2d_name // &
1597  ", data2d_value = array with length of ", &
1598  data_value_size, &
1599  " [", &
1600  data2d_value(1), &
1601  " ... ", &
1602  data2d_value(data_value_size), &
1603  "]"
1604  call nclayer_actionm(trim(action_str))
1605  end if
1606 #endif
1607 
1608  if (diag_data2d_store%data_lock) then
1609  call nclayer_error("Can't add new data - data have already been written and locked!")
1610  end if
1611 
1612  var_index = nc_diag_data2d_lookup_var(data2d_name)
1613 
1614  if (var_index == -1) then
1615  ! First, check to make sure we can still define new variables.
1616  if (diag_data2d_store%def_lock) then
1617  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1618  end if
1619 #ifdef _DEBUG_MEM_
1620  write (*, "(A, A, A, F)") "NEW data2d: ", data2d_name, " | First value: ", data2d_value
1621 #endif
1623 
1624  diag_data2d_store%total = diag_data2d_store%total + 1
1625 
1626  diag_data2d_store%names(diag_data2d_store%total) = data2d_name
1628 
1629  var_index = diag_data2d_store%total
1630  end if
1631 
1632  ! Get input size and do size checks!
1633  input_size = size(data2d_value)
1634 
1635  if ((diag_data2d_store%def_lock) .AND. &
1636  (size(data2d_value) > diag_data2d_store%max_lens(var_index))) then
1637  call nclayer_error("Cannot expand variable size after locking variable definitions!")
1638  end if
1639 
1640  ! We just need to add one entry...
1641  call nc_diag_data2d_resize_iarr(var_index, 1_i_llong)
1642  call nc_diag_data2d_resize_rsingle(input_size)
1643 
1644  ! Now add the actual entry!
1645  diag_data2d_store%m_rsingle(diag_data2d_store%acount(4) - input_size + 1:diag_data2d_store%acount(4)) = &
1646  data2d_value
1647  diag_data2d_store%stor_i_arr(var_index)%index_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1648  diag_data2d_store%acount(4) - input_size + 1
1649  diag_data2d_store%stor_i_arr(var_index)%length_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1650  input_size
1651  end subroutine nc_diag_data2d_rsingle
1652 
1653  ! nc_diag_data2d - input real(r_double)
1654  ! Corresponding NetCDF4 type: double
1655  subroutine nc_diag_data2d_rdouble(data2d_name, data2d_value)
1656  character(len=*), intent(in) :: data2d_name
1657  real(r_double), dimension(:), intent(in) :: data2d_value
1658 
1659  integer(i_long) :: var_index
1660  integer(i_llong) :: input_size
1661 
1662 #ifdef ENABLE_ACTION_MSGS
1663  character(len=1000) :: action_str
1664  integer(i_llong) :: data_value_size
1665 
1666  if (nclayer_enable_action) then
1667  data_value_size = size(data2d_value)
1668  write(action_str, "(A, I0, A, F0.5, A, F0.5, A)") &
1669  "nc_diag_data2d_rdouble(data2d_name = " // data2d_name // &
1670  ", data2d_value = array with length of ", &
1671  data_value_size, &
1672  " [", &
1673  data2d_value(1), &
1674  " ... ", &
1675  data2d_value(data_value_size), &
1676  "]"
1677  call nclayer_actionm(trim(action_str))
1678  end if
1679 #endif
1680 
1681  if (diag_data2d_store%data_lock) then
1682  call nclayer_error("Can't add new data - data have already been written and locked!")
1683  end if
1684 
1685  var_index = nc_diag_data2d_lookup_var(data2d_name)
1686 
1687  if (var_index == -1) then
1688  ! First, check to make sure we can still define new variables.
1689  if (diag_data2d_store%def_lock) then
1690  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1691  end if
1692 
1694 
1695  diag_data2d_store%total = diag_data2d_store%total + 1
1696 
1697  diag_data2d_store%names(diag_data2d_store%total) = data2d_name
1699 
1700  var_index = diag_data2d_store%total
1701  end if
1702 
1703  ! Get input size and do size checks!
1704  input_size = size(data2d_value)
1705 
1706  if ((diag_data2d_store%def_lock) .AND. &
1707  (size(data2d_value) > diag_data2d_store%max_lens(var_index))) then
1708  call nclayer_error("Cannot expand variable size after locking variable definitions!")
1709  end if
1710 
1711  ! We just need to add one entry...
1712  call nc_diag_data2d_resize_iarr(var_index, 1_i_llong)
1713  call nc_diag_data2d_resize_rdouble(input_size)
1714 
1715  ! Now add the actual entry!
1716  diag_data2d_store%m_rdouble(diag_data2d_store%acount(5) - input_size + 1:diag_data2d_store%acount(5)) = &
1717  data2d_value
1718  diag_data2d_store%stor_i_arr(var_index)%index_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1719  diag_data2d_store%acount(5) - input_size + 1
1720  diag_data2d_store%stor_i_arr(var_index)%length_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1721  input_size
1722  end subroutine nc_diag_data2d_rdouble
1723 
1724  ! nc_diag_data2d - input character(len=*)
1725  ! Corresponding NetCDF4 type: string? char?
1726  subroutine nc_diag_data2d_string(data2d_name, data2d_value)
1727  character(len=*), intent(in) :: data2d_name
1728  character(len=*), dimension(:), intent(in) :: data2d_value
1729 
1730  integer(i_long) :: var_index
1731  integer(i_long) :: max_str_len
1732  integer(i_llong) :: input_size
1733 
1734 #ifdef ENABLE_ACTION_MSGS
1735  character(len=1000) :: action_str
1736  integer(i_llong) :: data_value_size
1737 
1738  if (nclayer_enable_action) then
1739  data_value_size = size(data2d_value)
1740  write(action_str, "(A, I0, A, A)") &
1741  "nc_diag_data2d_string(data2d_name = " // data2d_name // &
1742  ", data2d_value = array with length of ", &
1743  data_value_size, &
1744  " [" // &
1745  trim(data2d_value(1)) // &
1746  " ... " // &
1747  trim(data2d_value(data_value_size)) // &
1748  "]"
1749  call nclayer_actionm(trim(action_str))
1750  end if
1751 #endif
1752 
1753  if (diag_data2d_store%data_lock) then
1754  call nclayer_error("Can't add new data - data have already been written and locked!")
1755  end if
1756 
1757  var_index = nc_diag_data2d_lookup_var(data2d_name)
1758 
1759  if (var_index == -1) then
1760  ! First, check to make sure we can still define new variables.
1761  if (diag_data2d_store%def_lock) then
1762  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1763  end if
1764 
1766 
1767  diag_data2d_store%total = diag_data2d_store%total + 1
1768 
1769  diag_data2d_store%names(diag_data2d_store%total) = data2d_name
1771 
1772  var_index = diag_data2d_store%total
1773  else
1774  ! Check max string length
1775 #ifdef _DEBUG_MEM_
1776  print *, "len_trim(data2d_value) = ", len_trim(data2d_value)
1777  print *, "diag_data2d_store%max_str_lens(var_index) = ", diag_data2d_store%max_str_lens(var_index)
1778 #endif
1779  end if
1780 
1781  ! Get input size and do size checks!
1782  input_size = size(data2d_value)
1783 
1784  if (diag_data2d_store%def_lock) then
1785  if (input_size > diag_data2d_store%max_lens(var_index)) &
1786  call nclayer_error("Cannot expand variable size after locking variable definitions!")
1787 
1788  ! Check max string length
1789  max_str_len = max_len_string_array(data2d_value, &
1790  int(input_size))
1791 
1792 #ifdef _DEBUG_MEM_
1793  print *, "max_str_len: ", max_str_len
1794  print *, "diag_data2d_store%max_str_lens(var_index): ", diag_data2d_store%max_str_lens(var_index)
1795 #endif
1796 
1797  if (max_str_len > diag_data2d_store%max_str_lens(var_index)) &
1798  call nclayer_error("Cannot expand variable string length after locking variable definitions!")
1799  end if
1800 
1801  ! We just need to add one entry...
1802  call nc_diag_data2d_resize_iarr(var_index, 1_i_llong)
1803  call nc_diag_data2d_resize_string(input_size)
1804 
1805  ! If trim isn't enabled, set our maximum string length here!
1806  if (.NOT. enable_trim) then
1807  if (diag_data2d_store%max_str_lens(var_index) == -1) then
1808  diag_data2d_store%max_str_lens(var_index) = len(data2d_value(1))
1809  else
1810  ! Validate that our non-first value isn't different from
1811  ! the initial string length
1812  if (max_len_notrim_string_array(data2d_value, int(input_size)) /= &
1813  diag_data2d_store%max_str_lens(var_index)) &
1814  call nclayer_error("Cannot change string size when trimming is disabled!")
1815  end if
1816  end if
1817 
1818  ! Now add the actual entry!
1819  diag_data2d_store%m_string(diag_data2d_store%acount(6) - input_size + 1:diag_data2d_store%acount(6)) = &
1820  data2d_value
1821  diag_data2d_store%stor_i_arr(var_index)%index_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1822  diag_data2d_store%acount(6) - input_size + 1
1823  diag_data2d_store%stor_i_arr(var_index)%length_arr(diag_data2d_store%stor_i_arr(var_index)%icount) = &
1824  input_size
1825  end subroutine nc_diag_data2d_string
1826 end module ncdw_data2d
type(diag_data2d), allocatable diag_data2d_store
Definition: ncdw_state.f90:18
subroutine nc_diag_data2d_load_def
subroutine nc_diag_data2d_allocmulti(multiplier)
integer(i_long), parameter nlayer_multi_base
Definition: ncdw_types.F90:27
integer, parameter, public i_byte
Definition: ncd_kinds.F90:45
integer(i_byte), parameter nlayer_short
Definition: ncdw_types.F90:11
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
subroutine nc_diag_data2d_rsingle(data2d_name, data2d_value)
logical init_done
Definition: ncdw_state.f90:9
subroutine nc_diag_data2d_resize_iarr_type(addl_num_entries)
integer(i_byte), parameter nlayer_double
Definition: ncdw_types.F90:14
subroutine nc_diag_data2d_resize_rdouble(addl_num_entries, update_acount_in)
logical enable_trim
Definition: ncdw_state.f90:12
subroutine nc_diag_data2d_resize_rsingle(addl_num_entries, update_acount_in)
type(diag_varattr), allocatable diag_varattr_store
Definition: ncdw_state.f90:19
integer(i_byte), parameter nlayer_string
Definition: ncdw_types.F90:15
subroutine nc_diag_data2d_write_data(flush_data_only)
subroutine nc_diag_data2d_long(data2d_name, data2d_value)
integer(i_long) ncid
Definition: ncdw_state.f90:8
subroutine nc_diag_data2d_resize_byte(addl_num_entries, update_acount_in)
integer function nc_diag_data2d_lookup_var(data2d_name)
subroutine nc_diag_data2d_prealloc_vars(num_of_addl_vars)
subroutine nc_diag_data2d_resize_long(addl_num_entries, update_acount_in)
subroutine nc_diag_data2d_prealloc_vars_storage_all(num_of_addl_slots)
integer(i_long), parameter nlayer_compression
Definition: ncdw_types.F90:21
integer(i_short), parameter nlayer_fill_short
Definition: ncdw_types.F90:30
integer function nc_diag_data2d_max_len_var(var_index)
subroutine nc_diag_data2d_resize_iarr(iarr_index, addl_num_entries, update_icount_in)
logical append_only
Definition: ncdw_state.f90:10
subroutine nc_diag_data2d_resize_short(addl_num_entries, update_acount_in)
subroutine nc_diag_data2d_write_def(internal)
subroutine nc_diag_data2d_expand
subroutine nc_diag_data2d_prealloc_vars_storage(nclayer_type, num_of_addl_slots)
integer, parameter, public i_short
Definition: ncd_kinds.F90:46
subroutine nc_diag_data2d_rdouble(data2d_name, data2d_value)
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_data2d_short(data2d_name, data2d_value)
subroutine nc_diag_data2d_set_strict(enable_strict)
subroutine nc_diag_data2d_string(data2d_name, data2d_value)
subroutine nc_diag_data2d_resize_string(addl_num_entries, update_acount_in)
real(r_single), parameter nlayer_fill_float
Definition: ncdw_types.F90:32
integer, parameter, public r_double
Definition: ncd_kinds.F90:80
integer(i_long), parameter nlayer_chunking
Definition: ncdw_types.F90:24
integer, parameter, public r_single
Definition: ncd_kinds.F90:79
subroutine nc_diag_data2d_byte(data2d_name, data2d_value)
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
real(r_double), parameter nlayer_fill_double
Definition: ncdw_types.F90:33
subroutine nc_diag_varattr_make_nobs_dim