FV3 Bundle
ncdw_metadata.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 ! metadata module - ncdw_metadata
17 !
19  ! Module that provides metadata variable storage support.
20  !
21  ! This module has all of the subroutines needed to store metadata
22  ! data. It includes the metadata storing subroutine
23  ! (nc_diag_chaninfo), subroutines for controlling chaninfo data
24  ! (loading definitions, saving definitions, saving data, etc.),
25  ! and preallocation subroutines.
26  !
27  ! Background:
28  ! metadata is an unlimited storage variable, with dimensions of
29  ! 1 x nobs, where nobs is an unlimited dimension. With unlimited
30  ! dimensions in this variable type, an unlimited amount of
31  ! metadata data can be stored in metadata variables.
32  !
33  ! Unlike chaninfo, we can NOT make any assumptions, since the
34  ! dimensions are now unlimited instead of fixed. This time, the
35  ! variables will have to be stored differently!
36  !
37  ! At the time of development, there were two ideas of approaching
38  ! this new type:
39  !
40  ! -> Same variable metadata storage, but now storing data inside
41  ! a derived type array instead of in a giant variable data
42  ! storage. The derived type array is now filled with the
43  ! various type arrays. Only one type array is allocated and
44  ! filled so that the array itself has the complete data, and
45  ! the array can be written out directly to NetCDF.
46  !
47  ! -> Same variable metadata storage, same variable data storage,
48  ! but with an addition of a derived type containing an array
49  ! of indicies referring to the location where the variable's
50  ! data is stored.
51  !
52  ! In the end, the array of indicies option was chosen. This was
53  ! due to these reasons:
54  !
55  ! -> Although writing the data would be rather quick (since the
56  ! data is already in a vector), several factors would make
57  ! the costs outweight this benefit. In particular...
58  !
59  ! -> Writing to the array would require more time, since it has
60  ! to seek to the allocatable array, then seek to the position,
61  ! and then write. This is due to the many other non-allocated
62  ! types in the derived type.
63  !
64  ! -> Reallocation would occur more often, since the arrays are
65  ! allocated by variable, and not allocated by type. Instead of
66  ! reallocating 6 times, it would reallocate (# of variables)
67  ! times, assuming all variables are appended to equally.
68  !
69  ! -> More counters (specifically, (# of variables) amount of
70  ! counters) will have to be used to keep track of the total
71  ! and the amount of data used for the allocatable arrays,
72  ! making reallocation even more costly.
73  !
74  ! -> With regards to the indicies option, appending and writing
75  ! times are equal. They essentially boil down to store index,
76  ! store value vs read index, read value.
77  !
78  ! -> The indicies array is stored in a derived type, but since
79  ! it is the sole element within the derived type, the array
80  ! access is much quicker.
81  !
82  ! -> Finally, the indicies array still uses the 6 type variable
83  ! data array storage, which just uses 6*2 counters and only
84  ! a maximum of 6 reallocations, which is much more efficient.
85  !
86  ! That said, we can therefore apply this method to our metadata
87  ! data storage!
88  !
89  ! Like with chaninfo, we support the following types:
90  ! i_byte, i_short, i_long, r_single, r_double, character(len=*)
91  !
92  ! Again, we store everything within a derived type, diag_metadata:
93  !
94  ! -> m_* - these arrays store the variable data for the types
95  ! listed above. This time, we organize the metadata using
96  ! array indicies for each variable, as mentionned above. More
97  ! details about the storage method to follow...
98  !
99  ! -> names - all of the metadata variable names! We'll be using
100  ! this array to store and lookup metadata variables, as well as
101  ! storing them!
102  !
103  ! -> types - all of the metadata variable types! These are byte
104  ! integers that get compared to our NLAYER_* type constants
105  ! (see: ncdw_types.F90).
106  !
107  !
108  !
109  ! -> stor_i_arr - for metadata, this is the star of the show! This
110  ! is an abbreviation for "storage index array". This is the
111  ! implementation of our indicies array idea. Since this is a
112  ! array of "diag_md_iarr" derived types, let's peek at the
113  ! derived type itself:
114  !
115  ! -> index_arr - the array of indicies. This stores all of the
116  ! variable data storage indicies, which indicate where our
117  ! data is stored within the variable type-specific data
118  ! storage.
119  !
120  ! -> icount - the number of indicies stored within this derived
121  ! type.
122  !
123  ! -> isize - the current indicies array size. Used for
124  ! reallocation when adding more elements.
125  !
126 
127  use ncd_kinds, only: i_byte, i_short, i_long, i_llong, r_single, &
128  r_double
129  use ncdw_state, only: init_done, ncid, append_only, &
130  enable_trim, &
137 
138  use ncdw_realloc, only: nc_diag_realloc
139  use ncdw_mresize, only: &
145 
148 
149  use ncdw_climsg, only: &
150 #ifdef enable_action_msgs
151  nclayer_enable_action, nclayer_actionm, &
152 #endif
153 #ifdef _DEBUG_MEM_
154  nclayer_debug, &
155 #endif
156  nclayer_error, nclayer_warning, nclayer_info, nclayer_check
157 
158  use netcdf, only: nf90_inquire, nf90_inquire_variable, &
159  nf90_inquire_dimension, nf90_def_dim, nf90_def_var, &
160  nf90_put_var, nf90_def_var_chunking, nf90_def_var_deflate, &
161  nf90_byte, nf90_short, nf90_int, nf90_float, nf90_double, &
162  nf90_char, nf90_max_name, nf90_chunked
163 
164  implicit none
165 
167  module procedure nc_diag_metadata_byte, &
171  end interface nc_diag_metadata
172 
173  contains
174  subroutine nc_diag_metadata_allocmulti(multiplier)
175  integer(i_long), intent(in) :: multiplier
176  if (init_done) then
177  ! # of times we needed to realloc simple metadata
178  ! also the multiplier factor for allocation (2^x)
179  diag_metadata_store%alloc_s_multi = multiplier
180 
181  ! # of times we needed to realloc metadata data storage
182  ! also the multiplier factor for allocation (2^x)
183  diag_metadata_store%alloc_m_multi = multiplier
184 
185  ! # of times we needed to realloc metadata INDEX data storage
186  ! also the multiplier factor for allocation (2^x)
187  diag_metadata_store%alloc_mi_multi = multiplier
188  end if
189  end subroutine nc_diag_metadata_allocmulti
190 
191  subroutine nc_diag_metadata_load_def
192  integer(i_long) :: ndims, nvars, var_index, type_index
193  integer(i_long) :: rel_index, i, nobs_size
194 
195  character(len=NF90_MAX_NAME) :: tmp_var_name
196  integer(i_long) :: tmp_var_type, tmp_var_ndims
197 
198  integer(i_long), dimension(:), allocatable :: tmp_var_dimids, tmp_var_dim_sizes
199  character(len=NF90_MAX_NAME) , allocatable :: tmp_var_dim_names(:)
200 
201  logical :: is_metadata_var
202 
203  ! Get top level info about the file!
204  call nclayer_check(nf90_inquire(ncid, ndimensions = ndims, &
205  nvariables = nvars))
206 
207  ! Now search for variables that use metadata storage!
208  ! Loop through each variable!
209  do var_index = 1, nvars
210  ! Grab number of dimensions and attributes first
211  call nclayer_check(nf90_inquire_variable(ncid, var_index, name = tmp_var_name, ndims = tmp_var_ndims))
212 
213  ! Allocate temporary variable dimids storage!
214  allocate(tmp_var_dimids(tmp_var_ndims))
215  allocate(tmp_var_dim_names(tmp_var_ndims))
216  allocate(tmp_var_dim_sizes(tmp_var_ndims))
217 
218  ! Grab the actual dimension IDs and attributes
219 
220  call nclayer_check(nf90_inquire_variable(ncid, var_index, dimids = tmp_var_dimids, &
221  xtype = tmp_var_type))
222 
223  if ((tmp_var_ndims == 1) .OR. &
224  ((tmp_var_ndims == 2) .AND. (tmp_var_type == nf90_char))) then
225  is_metadata_var = .false.
226 
227  do i = 1, tmp_var_ndims
228  call nclayer_check(nf90_inquire_dimension(ncid, tmp_var_dimids(i), tmp_var_dim_names(i), &
229  tmp_var_dim_sizes(i)))
230 
231  if (tmp_var_dim_names(i) == "nobs") then
232  nobs_size = tmp_var_dim_sizes(i)
233  if (tmp_var_type /= nf90_char) then
234  is_metadata_var = .true.
235  else if (tmp_var_type == nf90_char) then
236  if (index(tmp_var_dim_names(1), "_maxstrlen") /= 0) &
237  is_metadata_var = .true.
238  end if
239  end if
240  end do
241 
242  if (is_metadata_var) then
243  ! Expand things first!
245 
246  ! Add to the total!
247  diag_metadata_store%total = diag_metadata_store%total + 1
248 
249  ! Store name and type!
250  diag_metadata_store%names(diag_metadata_store%total) = trim(tmp_var_name)
251 
252  ! The relative index is the total nobs
253  rel_index = nobs_size
254 
255  if (tmp_var_type == nf90_byte) then
257  !call nc_diag_metadata_resize_byte(int8(diag_metadata_store%nchans), .FALSE.)
258  type_index = 1
259  else if (tmp_var_type == nf90_short) then
261  !call nc_diag_metadata_resize_short(int8(diag_metadata_store%nchans), .FALSE.)
262  type_index = 2
263  else if (tmp_var_type == nf90_int) then
265  !call nc_diag_metadata_resize_long(int8(diag_metadata_store%nchans), .FALSE.)
266  type_index = 3
267  else if (tmp_var_type == nf90_float) then
269  !call nc_diag_metadata_resize_rsingle(int8(diag_metadata_store%nchans), .FALSE.)
270  type_index = 4
271  else if (tmp_var_type == nf90_double) then
273  !call nc_diag_metadata_resize_rdouble(int8(diag_metadata_store%nchans), .FALSE.)
274  type_index = 5
275  else if (tmp_var_type == nf90_char) then
277  diag_metadata_store%max_str_lens(diag_metadata_store%total) = tmp_var_dim_sizes(1)
278  !call nc_diag_metadata_resize_string(int8(diag_metadata_store%nchans), .FALSE.)
279  type_index = 6
280  else
281  call nclayer_error("NetCDF4 type invalid!")
282  end if
283 
284 ! print *, trim(tmp_var_name), "rel index", rel_index
285 
286  ! Now add a relative position... based on the next position!
287 
288  ! Set relative index!
289  diag_metadata_store%rel_indexes(diag_metadata_store%total) = rel_index
290 
291  ! Set variable ID! Note that var_index here is the actual variable ID.
292  diag_metadata_store%var_ids(diag_metadata_store%total) = var_index
293 
294 ! print *, var_index
295 ! print *, diag_metadata_store%var_ids(diag_metadata_store%total)
296  end if
297  end if
298 
299  ! Deallocate
300  deallocate(tmp_var_dimids)
301  deallocate(tmp_var_dim_names)
302  deallocate(tmp_var_dim_sizes)
303  end do
304 
305  diag_metadata_store%def_lock = .true.
306  end subroutine nc_diag_metadata_load_def
307 
308  subroutine nc_diag_metadata_write_def(internal)
309  logical, intent(in), optional :: internal
310 
311  integer(i_byte) :: data_type
312  character(len=100) :: data_name
313 
314  integer(i_llong) :: curdatindex, j
315  integer(i_long) :: nc_data_type
316  integer(i_long) :: tmp_dim_id
317  character(len=120) :: data_dim_name
318 
319  character(len=:), allocatable :: string_arr(:)
320 
321 #ifdef ENABLE_ACTION_MSGS
322  character(len=1000) :: action_str
323 
324  if (nclayer_enable_action) then
325  if (present(internal)) then
326  write(action_str, "(A, L, A)") "nc_diag_metadata_write_def(internal = ", internal, ")"
327  else
328  write(action_str, "(A)") "nc_diag_metadata_write_def(internal = (not specified))"
329  end if
330  call nclayer_actionm(trim(action_str))
331  end if
332 #endif
333 
334  if (init_done) then
335  if (.NOT. diag_metadata_store%def_lock) then
336  ! Use global nobs ID!
337  ! Call subroutine to ensure the nobs dim is created already...
339 
340  do curdatindex = 1, diag_metadata_store%total
341  data_name = diag_metadata_store%names(curdatindex)
342  data_type = diag_metadata_store%types(curdatindex)
343 
344  call nclayer_info("metadata: defining " // trim(data_name))
345 
346  if (data_type == nlayer_byte) nc_data_type = nf90_byte
347  if (data_type == nlayer_short) nc_data_type = nf90_short
348  if (data_type == nlayer_long) nc_data_type = nf90_int
349  if (data_type == nlayer_float) nc_data_type = nf90_float
350  if (data_type == nlayer_double) nc_data_type = nf90_double
351  if (data_type == nlayer_string) nc_data_type = nf90_char
352 
353 #ifdef _DEBUG_MEM_
354  print *, "metadata part 1"
355 #endif
356 
357  if (data_type == nlayer_string) then
358  write (data_dim_name, "(A, A)") trim(data_name), "_maxstrlen"
359 
360  ! If trimming is enabled, we haven't found our max_str_len yet.
361  ! Go find it!
362  if (enable_trim) then
363  ! Dimension is # of chars by # of obs (unlimited)
364  allocate(character(10000) :: string_arr(diag_metadata_store%stor_i_arr(curdatindex)%icount))
365  do j = 1, diag_metadata_store%stor_i_arr(curdatindex)%icount
366  string_arr(j) = diag_metadata_store%m_string(diag_metadata_store%stor_i_arr(curdatindex)%index_arr(j))
367  end do
368 
369  ! Save the max string len
370  diag_metadata_store%max_str_lens(curdatindex) = max_len_string_array(string_arr, &
371  diag_metadata_store%stor_i_arr(curdatindex)%icount)
372 
373  deallocate(string_arr)
374  end if
375 
376  if (.NOT. append_only) &
377  call nclayer_check(nf90_def_dim(ncid, data_dim_name, &
378  diag_metadata_store%max_str_lens(curdatindex), tmp_dim_id))
379 
380 #ifdef _DEBUG_MEM_
381  print *, "Defining char var type..."
382 #endif
383 
384  if (.NOT. append_only) &
385  call nclayer_check(nf90_def_var(ncid, data_name, nc_data_type, &
386  (/ tmp_dim_id, diag_varattr_store%nobs_dim_id /), &
387  diag_metadata_store%var_ids(curdatindex)))
388 
389 #ifdef _DEBUG_MEM_
390  print *, "Done defining char var type..."
391 #endif
392  else
393  if (.NOT. append_only) &
394  call nclayer_check(nf90_def_var(ncid, data_name, nc_data_type, diag_varattr_store%nobs_dim_id, &
395  diag_metadata_store%var_ids(curdatindex)))
396  end if
397 
398 #ifdef _DEBUG_MEM_
399  print *, "metadata part 2"
400 #endif
401 
402  call nc_diag_varattr_add_var(diag_metadata_store%names(curdatindex), &
403  diag_metadata_store%types(curdatindex), &
404  diag_metadata_store%var_ids(curdatindex))
405 
406  ! Enable compression
407  ! Args: ncid, varid, enable_shuffle (yes), enable_deflate (yes), deflate_level
408 #ifdef _DEBUG_MEM_
409  print *, "Defining compression 1 (chunking)..."
410 #endif
411 
412  if (.NOT. append_only) then
413  if (data_type == nlayer_string) then
414  call nclayer_check(nf90_def_var_chunking(ncid, diag_metadata_store%var_ids(curdatindex), &
415  nf90_chunked, (/ diag_metadata_store%max_str_lens(curdatindex), nlayer_chunking /)))
416  else
417  call nclayer_check(nf90_def_var_chunking(ncid, diag_metadata_store%var_ids(curdatindex), &
418  nf90_chunked, (/ nlayer_chunking /)))
419  end if
420 
421 #ifdef _DEBUG_MEM_
422  print *, "Defining compression 2 (gzip)..."
423 #endif
424  call nclayer_check(nf90_def_var_deflate(ncid, diag_metadata_store%var_ids(curdatindex), &
425  1, 1, nlayer_compression))
426 
427 #ifdef _DEBUG_MEM_
428  print *, "Done defining compression..."
429 #endif
430  end if
431 
432  ! Lock the definitions!
433  diag_metadata_store%def_lock = .true.
434  end do
435  else
436  if(.NOT. present(internal)) &
437  call nclayer_error("Can't write definitions - definitions have already been written and locked!")
438  end if
439  end if
440  end subroutine nc_diag_metadata_write_def
441 
442  subroutine nc_diag_metadata_write_data(flush_data_only)
443  ! Optional internal flag to only flush data - if this is
444  ! true, data flushing will be performed, and the data will
445  ! NOT be locked.
446  logical, intent(in), optional :: flush_data_only
447 
448  integer(i_byte) :: data_type
449  character(len=100) :: data_name
450 
451  integer(i_long) :: curdatindex, j
452 
453  integer(i_byte), dimension(:), allocatable :: byte_arr
454  integer(i_short),dimension(:), allocatable :: short_arr
455  integer(i_long), dimension(:), allocatable :: long_arr
456  real(r_single), dimension(:), allocatable :: rsingle_arr
457  real(r_double), dimension(:), allocatable :: rdouble_arr
458  character(len=:), allocatable :: string_arr(:)
459 
460  integer(i_llong) :: string_arr_maxlen
461 
462  integer(i_llong) :: data_length_counter
463  character(len=100) :: counter_data_name
464  integer(i_llong) :: current_length_count
465  character(len=1000) :: data_uneven_msg
466 
467 #ifdef ENABLE_ACTION_MSGS
468  character(len=1000) :: action_str
469 
470  if (nclayer_enable_action) then
471  if (present(flush_data_only)) then
472  write(action_str, "(A, L, A)") "nc_diag_metadata_write_data(flush_data_only = ", flush_data_only, ")"
473  else
474  write(action_str, "(A)") "nc_diag_metadata_write_data(flush_data_only = (not specified))"
475  end if
476  call nclayer_actionm(trim(action_str))
477  end if
478 #endif
479  ! Initialization MUST occur here, not in decl...
480  ! Otherwise, it'll initialize once, and never again...
481  !
482  ! This will cause scary issues in the future, where closing
483  ! and opening a new file shows strange errors about a file
484  ! opened in the past...
485  data_length_counter = -1
486  current_length_count = -1
487 
488  if (init_done .AND. allocated(diag_metadata_store)) then
489  if (.NOT. diag_metadata_store%data_lock) then
490  do curdatindex = 1, diag_metadata_store%total
491 #ifdef _DEBUG_MEM_
492  print *, curdatindex
493 #endif
494  data_name = diag_metadata_store%names(curdatindex)
495  data_type = diag_metadata_store%types(curdatindex)
496 
497  call nclayer_info("metadata: writing " // trim(data_name))
498 
499  ! Warn about data inconsistencies
500  if (.NOT. (present(flush_data_only) .AND. flush_data_only)) then
501  current_length_count = diag_metadata_store%stor_i_arr(curdatindex)%icount + &
502  diag_metadata_store%rel_indexes(curdatindex)
503 
504  if (data_length_counter == -1) then
505  data_length_counter = current_length_count
506  counter_data_name = data_name
507  else
508  if (data_length_counter /= current_length_count) then
509  ! Show message!
510  ! NOTE - I0 and TRIM are Fortran 95 specs
511  write (data_uneven_msg, "(A, I0, A, I0, A)") "Amount of data written in " // &
512  trim(data_name) // " (", &
513  current_length_count, &
514  ")" // char(10) // &
515  " differs from variable " // trim(counter_data_name) // &
516  " (", data_length_counter, ")!"
517 
518  if (diag_metadata_store%strict_check) then
519  call nclayer_error(trim(data_uneven_msg))
520  else
521  call nclayer_warning(trim(data_uneven_msg))
522  end if
523  end if
524  end if
525  end if
526 
527  ! Make sure we have data to write in the first place!
528  if (diag_metadata_store%stor_i_arr(curdatindex)%icount > 0) then
529  if (data_type == nlayer_byte) then
530  allocate(byte_arr(diag_metadata_store%stor_i_arr(curdatindex)%icount))
531  do j = 1, diag_metadata_store%stor_i_arr(curdatindex)%icount
532  byte_arr(j) = diag_metadata_store%m_byte(diag_metadata_store%stor_i_arr(curdatindex)%index_arr(j))
533  end do
534  call nclayer_check(nf90_put_var(&
535  ncid, diag_metadata_store%var_ids(curdatindex), &
536  byte_arr, &
537  (/ 1 + diag_metadata_store%rel_indexes(curdatindex) /) &
538  ))
539 
540  deallocate(byte_arr)
541  else if (data_type == nlayer_short) then
542  allocate(short_arr(diag_metadata_store%stor_i_arr(curdatindex)%icount))
543  do j = 1, diag_metadata_store%stor_i_arr(curdatindex)%icount
544  short_arr(j) = diag_metadata_store%m_short(diag_metadata_store%stor_i_arr(curdatindex)%index_arr(j))
545  end do
546  call nclayer_check(nf90_put_var(&
547  ncid, diag_metadata_store%var_ids(curdatindex), &
548  short_arr, &
549  (/ 1 + diag_metadata_store%rel_indexes(curdatindex) /) &
550  ))
551 
552  deallocate(short_arr)
553  else if (data_type == nlayer_long) then
554  allocate(long_arr(diag_metadata_store%stor_i_arr(curdatindex)%icount))
555  do j = 1, diag_metadata_store%stor_i_arr(curdatindex)%icount
556  long_arr(j) = diag_metadata_store%m_long(diag_metadata_store%stor_i_arr(curdatindex)%index_arr(j))
557  end do
558 
559  call nclayer_check(nf90_put_var(&
560  ncid, diag_metadata_store%var_ids(curdatindex), &
561  long_arr, &
562  (/ 1 + diag_metadata_store%rel_indexes(curdatindex) /) &
563  ))
564 
565  deallocate(long_arr)
566  else if (data_type == nlayer_float) then
567  allocate(rsingle_arr(diag_metadata_store%stor_i_arr(curdatindex)%icount))
568  do j = 1, diag_metadata_store%stor_i_arr(curdatindex)%icount
569  rsingle_arr(j) = diag_metadata_store%m_rsingle(diag_metadata_store%stor_i_arr(curdatindex)%index_arr(j))
570  end do
571 
572  call nclayer_check(nf90_put_var(&
573  ncid, diag_metadata_store%var_ids(curdatindex), &
574  rsingle_arr, &
575  (/ 1 + diag_metadata_store%rel_indexes(curdatindex) /) &
576  ))
577 
578  deallocate(rsingle_arr)
579  else if (data_type == nlayer_double) then
580  allocate(rdouble_arr(diag_metadata_store%stor_i_arr(curdatindex)%icount))
581  do j = 1, diag_metadata_store%stor_i_arr(curdatindex)%icount
582  rdouble_arr(j) = diag_metadata_store%m_rdouble(diag_metadata_store%stor_i_arr(curdatindex)%index_arr(j))
583  end do
584 
585  call nclayer_check(nf90_put_var(&
586  ncid, diag_metadata_store%var_ids(curdatindex), &
587  rdouble_arr, &
588  (/ 1 + diag_metadata_store%rel_indexes(curdatindex) /) &
589  ))
590  deallocate(rdouble_arr)
591  else if (data_type == nlayer_string) then
592  ! Only get maximum if we haven't already done that in the define step!
593  if (diag_metadata_store%max_str_lens(curdatindex) == -1) then
594  allocate(character(10000) :: string_arr(diag_metadata_store%stor_i_arr(curdatindex)%icount))
595  do j = 1, diag_metadata_store%stor_i_arr(curdatindex)%icount
596  string_arr(j) = diag_metadata_store%m_string(diag_metadata_store%stor_i_arr(curdatindex)%index_arr(j))
597  end do
598 
599  string_arr_maxlen = max_len_string_array(string_arr, &
600  diag_metadata_store%stor_i_arr(curdatindex)%icount)
601 
602  deallocate(string_arr)
603  else
604  string_arr_maxlen = diag_metadata_store%max_str_lens(curdatindex)
605  end if
606 
607  allocate(character(string_arr_maxlen) :: string_arr(diag_metadata_store%stor_i_arr(curdatindex)%icount))
608  do j = 1, diag_metadata_store%stor_i_arr(curdatindex)%icount
609  string_arr(j) = diag_metadata_store%m_string(diag_metadata_store%stor_i_arr(curdatindex)%index_arr(j))
610  end do
611 
612  call nclayer_check(nf90_put_var(&
613  ncid, diag_metadata_store%var_ids(curdatindex), &
614  string_arr, &
615  (/ 1, 1 + diag_metadata_store%rel_indexes(curdatindex) /) &
616  ))
617  deallocate(string_arr)
618  end if
619 
620  ! Check for data flushing, and if so, update the relative indexes
621  ! and set icount to 0.
622  if (present(flush_data_only) .AND. flush_data_only) then
623  diag_metadata_store%rel_indexes(curdatindex) = &
624  diag_metadata_store%rel_indexes(curdatindex) + &
625  diag_metadata_store%stor_i_arr(curdatindex)%icount
626  diag_metadata_store%stor_i_arr(curdatindex)%icount = 0
627 
628 #ifdef _DEBUG_MEM_
629  print *, "diag_metadata_store%rel_indexes(curdatindex) is now:"
630  print *, diag_metadata_store%rel_indexes(curdatindex)
631 #endif
632  end if
633 
634  end if
635  end do
636 
637  if (present(flush_data_only) .AND. flush_data_only) then
638 #ifdef _DEBUG_MEM_
639  print *, "In buffer flush mode!"
640 #endif
641 
642  ! We need to reset all array counts to zero!
643  diag_metadata_store%acount = 0
644  else
645  ! Lock data writing
646  diag_metadata_store%data_lock = .true.
647 #ifdef _DEBUG_MEM_
648  print *, "In data lock mode!"
649 #endif
650  end if
651  else
652  call nclayer_error("Can't write data - data have already been written and locked!")
653  end if
654  else
655  call nclayer_error("Can't write data - NetCDF4 layer not initialized yet!")
656  end if
657 
658 #ifdef _DEBUG_MEM_
659  print *, "All done writing metadata data"
660 #endif
661  end subroutine nc_diag_metadata_write_data
662 
663  subroutine nc_diag_metadata_set_strict(enable_strict)
664  logical, intent(in) :: enable_strict
665 
666  if (init_done .AND. allocated(diag_metadata_store)) then
667  diag_metadata_store%strict_check = enable_strict
668  else
669  call nclayer_error("Can't set strictness level for metadata - NetCDF4 layer not initialized yet!")
670  end if
671  end subroutine nc_diag_metadata_set_strict
672 
673  ! Preallocate variable name/type/etc. storage.
674  subroutine nc_diag_metadata_prealloc_vars(num_of_addl_vars)
675  integer(i_llong), intent(in) :: num_of_addl_vars
676 #ifdef ENABLE_ACTION_MSGS
677  character(len=1000) :: action_str
678 
679  if (nclayer_enable_action) then
680  write(action_str, "(A, I0, A)") "nc_diag_metadata_prealloc_vars(num_of_addl_vars = ", num_of_addl_vars, ")"
681  call nclayer_actionm(trim(action_str))
682  end if
683 #endif
684  if (init_done .AND. allocated(diag_metadata_store)) then
685  if (allocated(diag_metadata_store%names)) then
686  if (diag_metadata_store%total >= size(diag_metadata_store%names)) then
687  call nc_diag_realloc(diag_metadata_store%names, num_of_addl_vars)
688  end if
689  else
690  allocate(diag_metadata_store%names(nlayer_default_ent + num_of_addl_vars))
691  end if
692 
693  if (allocated(diag_metadata_store%types)) then
694  if (diag_metadata_store%total >= size(diag_metadata_store%types)) then
695  call nc_diag_realloc(diag_metadata_store%types, num_of_addl_vars)
696  end if
697  else
698  allocate(diag_metadata_store%types(nlayer_default_ent + num_of_addl_vars))
699  end if
700 
701  if (allocated(diag_metadata_store%stor_i_arr)) then
702  if (diag_metadata_store%total >= size(diag_metadata_store%stor_i_arr)) then
703  call nc_diag_metadata_resize_iarr_type(num_of_addl_vars)
704  end if
705  else
706  allocate(diag_metadata_store%stor_i_arr(nlayer_default_ent + num_of_addl_vars))
707  end if
708 
709  if (allocated(diag_metadata_store%var_ids)) then
710  if (diag_metadata_store%total >= size(diag_metadata_store%var_ids)) then
711  call nc_diag_realloc(diag_metadata_store%var_ids, num_of_addl_vars)
712  end if
713  else
714  allocate(diag_metadata_store%var_ids(nlayer_default_ent + num_of_addl_vars))
715  diag_metadata_store%var_ids = -1
716  end if
717 
718  if (allocated(diag_metadata_store%alloc_sia_multi)) then
719  if (diag_metadata_store%total >= size(diag_metadata_store%alloc_sia_multi)) then
720  call nc_diag_realloc(diag_metadata_store%alloc_sia_multi, num_of_addl_vars)
721  end if
722  else
723  allocate(diag_metadata_store%alloc_sia_multi(nlayer_default_ent + num_of_addl_vars))
724  diag_metadata_store%alloc_sia_multi = 0
725  end if
726 
727  if (allocated(diag_metadata_store%max_str_lens)) then
728  if (diag_metadata_store%total >= size(diag_metadata_store%max_str_lens)) then
729  call nc_diag_realloc(diag_metadata_store%max_str_lens, num_of_addl_vars)
730  end if
731  else
732  allocate(diag_metadata_store%max_str_lens(nlayer_default_ent + num_of_addl_vars))
733  diag_metadata_store%max_str_lens = -1
734  end if
735 
736  if (allocated(diag_metadata_store%rel_indexes)) then
737  if (diag_metadata_store%total >= size(diag_metadata_store%rel_indexes)) then
738  call nc_diag_realloc(diag_metadata_store%rel_indexes, num_of_addl_vars)
739  end if
740  else
741  allocate(diag_metadata_store%rel_indexes(nlayer_default_ent + num_of_addl_vars))
742  diag_metadata_store%rel_indexes = 0
743  end if
744 
745  diag_metadata_store%prealloc_total = diag_metadata_store%prealloc_total + num_of_addl_vars
746  else
747  call nclayer_error("NetCDF4 layer not initialized yet!")
748  endif
749  end subroutine nc_diag_metadata_prealloc_vars
750 
751  ! Preallocate actual variable data storage
752  subroutine nc_diag_metadata_prealloc_vars_storage(nclayer_type, num_of_addl_slots)
753  integer(i_byte), intent(in) :: nclayer_type
754  integer(i_llong), intent(in) :: num_of_addl_slots
755 
756 #ifdef ENABLE_ACTION_MSGS
757  character(len=1000) :: action_str
758 
759  if (nclayer_enable_action) then
760  write(action_str, "(A, I0, A, I0, A)") "nc_diag_metadata_prealloc_vars_storage(nclayer_type = ", nclayer_type, ", num_of_addl_slots = ", num_of_addl_slots, ")"
761  call nclayer_actionm(trim(action_str))
762  end if
763 #endif
764 
765  if (nclayer_type == nlayer_byte) then
766  call nc_diag_metadata_resize_byte(num_of_addl_slots, .false.)
767  else if (nclayer_type == nlayer_short) then
768  call nc_diag_metadata_resize_short(num_of_addl_slots, .false.)
769  else if (nclayer_type == nlayer_long) then
770  call nc_diag_metadata_resize_long(num_of_addl_slots, .false.)
771  else if (nclayer_type == nlayer_float) then
772  call nc_diag_metadata_resize_rsingle(num_of_addl_slots, .false.)
773  else if (nclayer_type == nlayer_double) then
774  call nc_diag_metadata_resize_rdouble(num_of_addl_slots, .false.)
775  else if (nclayer_type == nlayer_string) then
776  call nc_diag_metadata_resize_string(num_of_addl_slots, .false.)
777  else
778  call nclayer_error("Invalid type specified for variable storage preallocation!")
779  end if
781 
782  ! Preallocate index storage
783  subroutine nc_diag_metadata_prealloc_vars_storage_all(num_of_addl_slots)
784  integer(i_llong), intent(in) :: num_of_addl_slots
785  integer(i_long) :: i
786 
787 #ifdef ENABLE_ACTION_MSGS
788  character(len=1000) :: action_str
789 
790  if (nclayer_enable_action) then
791  write(action_str, "(A, I0, A)") "nc_diag_metadata_prealloc_vars_storage_all(num_of_addl_slots = ", num_of_addl_slots, ")"
792  call nclayer_actionm(trim(action_str))
793  end if
794 #endif
795 
796  do i = 1, diag_metadata_store%prealloc_total
797  call nc_diag_metadata_resize_iarr(i, num_of_addl_slots, .false.)
798  end do
800 
801  subroutine nc_diag_metadata_expand
802  integer(i_llong) :: addl_fields
803 
804  ! Did we realloc at all?
805  logical :: meta_realloc
806 
807  meta_realloc = .false.
808 
809  if (init_done .AND. allocated(diag_metadata_store)) then
810  addl_fields = 1 + (nlayer_default_ent * (nlayer_multi_base ** diag_metadata_store%alloc_s_multi))
811 
812 #ifdef _DEBUG_MEM_
813  call nclayer_debug("INITIAL value of diag_metadata_store%alloc_s_multi:")
814  print *, diag_metadata_store%alloc_s_multi
815 #endif
816 
817  if (allocated(diag_metadata_store%names)) then
818  if (diag_metadata_store%total >= size(diag_metadata_store%names)) then
819 #ifdef _DEBUG_MEM_
820  call nclayer_debug("Reallocating diag_metadata_store%names...")
821  print *, (nlayer_multi_base ** diag_metadata_store%alloc_s_multi)
822  print *, addl_fields
823 #endif
824  call nc_diag_realloc(diag_metadata_store%names, addl_fields)
825 #ifdef _DEBUG_MEM_
826  call nclayer_debug("Reallocated diag_metadata_store%names. Size:")
827  print *, size(diag_metadata_store%names)
828 #endif
829  meta_realloc = .true.
830  end if
831  else
832 #ifdef _DEBUG_MEM_
833  call nclayer_debug("Allocating diag_metadata_store%names for first time...")
834  print *, nlayer_default_ent
835 #endif
836 
837  allocate(diag_metadata_store%names(nlayer_default_ent))
838 
839 #ifdef _DEBUG_MEM_
840  call nclayer_debug("Allocated diag_metadata_store%names. Size:")
841  print *, size(diag_metadata_store%names)
842 #endif
843  end if
844 
845  if (allocated(diag_metadata_store%types)) then
846  if (diag_metadata_store%total >= size(diag_metadata_store%types)) then
847 #ifdef _DEBUG_MEM_
848  call nclayer_debug("Reallocating diag_metadata_store%types...")
849  print *, (nlayer_multi_base ** diag_metadata_store%alloc_s_multi)
850  print *, addl_fields
851 #endif
852  call nc_diag_realloc(diag_metadata_store%types, addl_fields)
853  meta_realloc = .true.
854  end if
855  else
856  allocate(diag_metadata_store%types(nlayer_default_ent))
857  end if
858 
859  if (allocated(diag_metadata_store%stor_i_arr)) then
860  if (diag_metadata_store%total >= size(diag_metadata_store%stor_i_arr)) then
861 #ifdef _DEBUG_MEM_
862  call nclayer_debug("Reallocating diag_metadata_store%stor_i_arr...")
863  print *, (nlayer_multi_base ** diag_metadata_store%alloc_s_multi)
864  print *, (1 + (nlayer_default_ent * (nlayer_multi_base ** diag_metadata_store%alloc_s_multi)))
865 #endif
866  call nc_diag_metadata_resize_iarr_type(addl_fields)
867 
868  meta_realloc = .true.
869  end if
870  else
871  allocate(diag_metadata_store%stor_i_arr(nlayer_default_ent))
872  end if
873 
874  if (allocated(diag_metadata_store%var_ids)) then
875  if (diag_metadata_store%total >= size(diag_metadata_store%var_ids)) then
876  call nc_diag_realloc(diag_metadata_store%var_ids, addl_fields)
877  meta_realloc = .true.
878  end if
879  else
880  allocate(diag_metadata_store%var_ids(nlayer_default_ent))
881  diag_metadata_store%var_ids = -1
882  end if
883 
884  if (allocated(diag_metadata_store%alloc_sia_multi)) then
885  if (diag_metadata_store%total >= size(diag_metadata_store%alloc_sia_multi)) then
886  call nc_diag_realloc(diag_metadata_store%alloc_sia_multi, addl_fields)
887  meta_realloc = .true.
888  end if
889  else
890  allocate(diag_metadata_store%alloc_sia_multi(nlayer_default_ent))
891  diag_metadata_store%alloc_sia_multi = 0
892  end if
893 
894  if (allocated(diag_metadata_store%max_str_lens)) then
895  if (diag_metadata_store%total >= size(diag_metadata_store%max_str_lens)) then
896  call nc_diag_realloc(diag_metadata_store%max_str_lens, addl_fields)
897  meta_realloc = .true.
898  end if
899  else
900  allocate(diag_metadata_store%max_str_lens(nlayer_default_ent))
901  diag_metadata_store%max_str_lens = -1
902  end if
903 
904  if (allocated(diag_metadata_store%rel_indexes)) then
905  if (diag_metadata_store%total >= size(diag_metadata_store%rel_indexes)) then
906  call nc_diag_realloc(diag_metadata_store%rel_indexes, addl_fields)
907  meta_realloc = .true.
908  end if
909  else
910  allocate(diag_metadata_store%rel_indexes(nlayer_default_ent))
911  diag_metadata_store%rel_indexes = 0
912  end if
913 
914  if (meta_realloc) then
915  diag_metadata_store%alloc_s_multi = diag_metadata_store%alloc_s_multi + 1
916 #ifdef _DEBUG_MEM_
917  print *, "Incrementing alloc_s_multi... new value:"
918  print *, diag_metadata_store%alloc_s_multi
919 #endif
920  endif
921  else
922  call nclayer_error("NetCDF4 layer not initialized yet!")
923  endif
924 
925  end subroutine nc_diag_metadata_expand
926 
927  function nc_diag_metadata_lookup_var(metadata_name) result(ind)
928  character(len=*), intent(in) :: metadata_name
929  integer :: i, ind
930 
931  ind = -1
932 
933  if (init_done .AND. allocated(diag_metadata_store)) then
934  do i = 1, diag_metadata_store%total
935  if (diag_metadata_store%names(i) == metadata_name) then
936  ind = i
937  exit
938  end if
939  end do
940  end if
941  end function nc_diag_metadata_lookup_var
942 
943  ! nc_diag_metadata - input integer(i_byte)
944  ! Corresponding NetCDF4 type: byte
945  subroutine nc_diag_metadata_byte(metadata_name, metadata_value)
946  character(len=*), intent(in) :: metadata_name
947  integer(i_byte), intent(in) :: metadata_value
948 
949  integer(i_long) :: var_index
950 
951 #ifdef ENABLE_ACTION_MSGS
952  character(len=1000) :: action_str
953 
954  if (nclayer_enable_action) then
955  write(action_str, "(A, I0, A)") "nc_diag_metadata_byte(metadata_name = " // metadata_name // ", metadata_value = ", metadata_value, ")"
956  call nclayer_actionm(trim(action_str))
957  end if
958 #endif
959 
960  if (diag_metadata_store%data_lock) then
961  call nclayer_error("Can't add new data - data have already been written and locked!")
962  end if
963 
964  var_index = nc_diag_metadata_lookup_var(metadata_name)
965 
966  if (var_index == -1) then
967  ! First, check to make sure we can still define new variables.
968  if (diag_metadata_store%def_lock) then
969  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
970  end if
971 
973 
974  diag_metadata_store%total = diag_metadata_store%total + 1
975 
976  diag_metadata_store%names(diag_metadata_store%total) = metadata_name
978 
979  var_index = diag_metadata_store%total
980  end if
981 
982  ! We just need to add one entry...
983  call nc_diag_metadata_resize_iarr(var_index, 1_i_llong)
984  call nc_diag_metadata_resize_byte(1_i_llong)
985 
986  ! Now add the actual entry!
987  diag_metadata_store%m_byte(diag_metadata_store%acount(1)) = metadata_value
988  diag_metadata_store%stor_i_arr(var_index)%index_arr(diag_metadata_store%stor_i_arr(var_index)%icount) = &
989  diag_metadata_store%acount(1)
990  end subroutine nc_diag_metadata_byte
991 
992  ! nc_diag_metadata - input integer(i_short)
993  ! Corresponding NetCDF4 type: short
994  subroutine nc_diag_metadata_short(metadata_name, metadata_value)
995  character(len=*), intent(in) :: metadata_name
996  integer(i_short), intent(in) :: metadata_value
997 
998  integer(i_long) :: var_index
999 
1000 #ifdef ENABLE_ACTION_MSGS
1001  character(len=1000) :: action_str
1002 
1003  if (nclayer_enable_action) then
1004  write(action_str, "(A, I0, A)") "nc_diag_metadata_short(metadata_name = " // metadata_name // ", metadata_value = ", metadata_value, ")"
1005  call nclayer_actionm(trim(action_str))
1006  end if
1007 #endif
1008 
1009  if (diag_metadata_store%data_lock) then
1010  call nclayer_error("Can't add new data - data have already been written and locked!")
1011  end if
1012 
1013  var_index = nc_diag_metadata_lookup_var(metadata_name)
1014 
1015  if (var_index == -1) then
1016  ! First, check to make sure we can still define new variables.
1017  if (diag_metadata_store%def_lock) then
1018  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1019  end if
1020 
1022 
1023  diag_metadata_store%total = diag_metadata_store%total + 1
1024 
1025  diag_metadata_store%names(diag_metadata_store%total) = metadata_name
1027 
1028  var_index = diag_metadata_store%total
1029  end if
1030 
1031  ! We just need to add one entry...
1032  call nc_diag_metadata_resize_iarr(var_index, 1_i_llong)
1033  call nc_diag_metadata_resize_short(1_i_llong)
1034 
1035  ! Now add the actual entry!
1036  diag_metadata_store%m_short(diag_metadata_store%acount(2)) = metadata_value
1037  diag_metadata_store%stor_i_arr(var_index)%index_arr(diag_metadata_store%stor_i_arr(var_index)%icount) = &
1038  diag_metadata_store%acount(2)
1039  end subroutine nc_diag_metadata_short
1040 
1041  ! nc_diag_metadata - input integer(i_long)
1042  ! Corresponding NetCDF4 type: int (old: long)
1043  subroutine nc_diag_metadata_long(metadata_name, metadata_value)
1044  character(len=*), intent(in) :: metadata_name
1045  integer(i_long), intent(in) :: metadata_value
1046 
1047  integer(i_long) :: var_index
1048 
1049 #ifdef ENABLE_ACTION_MSGS
1050  character(len=1000) :: action_str
1051 
1052  if (nclayer_enable_action) then
1053  write(action_str, "(A, I0, A)") "nc_diag_metadata_long(metadata_name = " // metadata_name // ", metadata_value = ", metadata_value, ")"
1054  call nclayer_actionm(trim(action_str))
1055  end if
1056 #endif
1057 
1058  if (diag_metadata_store%data_lock) then
1059  call nclayer_error("Can't add new data - data have already been written and locked!")
1060  end if
1061 
1062  var_index = nc_diag_metadata_lookup_var(metadata_name)
1063 
1064  if (var_index == -1) then
1065  ! First, check to make sure we can still define new variables.
1066  if (diag_metadata_store%def_lock) then
1067  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1068  end if
1069 
1071 
1072  diag_metadata_store%total = diag_metadata_store%total + 1
1073 
1074  diag_metadata_store%names(diag_metadata_store%total) = metadata_name
1076 
1077  var_index = diag_metadata_store%total
1078  end if
1079 
1080 #ifdef _DEBUG_MEM_
1081  call nclayer_debug("Current total:")
1082  print *, diag_metadata_store%total
1083 #endif
1084 
1085  ! We just need to add one entry...
1086  call nc_diag_metadata_resize_iarr(var_index, 1_i_llong)
1087  call nc_diag_metadata_resize_long(1_i_llong)
1088 
1089  ! Now add the actual entry!
1090  diag_metadata_store%m_long(diag_metadata_store%acount(3)) = metadata_value
1091  diag_metadata_store%stor_i_arr(var_index)%index_arr(diag_metadata_store%stor_i_arr(var_index)%icount) = &
1092  diag_metadata_store%acount(3)
1093  end subroutine nc_diag_metadata_long
1094 
1095  ! nc_diag_metadata - input real(r_single)
1096  ! Corresponding NetCDF4 type: float (or real)
1097  subroutine nc_diag_metadata_rsingle(metadata_name, metadata_value)
1098  character(len=*), intent(in) :: metadata_name
1099  real(r_single), intent(in) :: metadata_value
1100 
1101  integer(i_long) :: var_index
1102 
1103 #ifdef ENABLE_ACTION_MSGS
1104  character(len=1000) :: action_str
1105 
1106  if (nclayer_enable_action) then
1107  write(action_str, "(A, F0.5, A)") "nc_diag_metadata_rsingle(metadata_name = " // metadata_name // ", metadata_value = ", metadata_value, ")"
1108  call nclayer_actionm(trim(action_str))
1109  end if
1110 #endif
1111 
1112  if (diag_metadata_store%data_lock) then
1113  call nclayer_error("Can't add new data - data have already been written and locked!")
1114  end if
1115 
1116  var_index = nc_diag_metadata_lookup_var(metadata_name)
1117 
1118  if (var_index == -1) then
1119  ! First, check to make sure we can still define new variables.
1120  if (diag_metadata_store%def_lock) then
1121  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1122  end if
1123 #ifdef _DEBUG_MEM_
1124  write (*, "(A, A, A, F)") "NEW METADATA: ", metadata_name, " | First value: ", metadata_value
1125 #endif
1127 
1128  diag_metadata_store%total = diag_metadata_store%total + 1
1129 
1130  diag_metadata_store%names(diag_metadata_store%total) = metadata_name
1132 
1133  var_index = diag_metadata_store%total
1134  end if
1135 
1136  ! We just need to add one entry...
1137  call nc_diag_metadata_resize_iarr(var_index, 1_i_llong)
1138  call nc_diag_metadata_resize_rsingle(1_i_llong)
1139 
1140  ! Now add the actual entry!
1141  diag_metadata_store%m_rsingle(diag_metadata_store%acount(4)) = metadata_value
1142  diag_metadata_store%stor_i_arr(var_index)%index_arr(diag_metadata_store%stor_i_arr(var_index)%icount) = &
1143  diag_metadata_store%acount(4)
1144  end subroutine nc_diag_metadata_rsingle
1145 
1146  ! nc_diag_metadata - input real(r_double)
1147  ! Corresponding NetCDF4 type: double
1148  subroutine nc_diag_metadata_rdouble(metadata_name, metadata_value)
1149  character(len=*), intent(in) :: metadata_name
1150  real(r_double), intent(in) :: metadata_value
1151 
1152  integer(i_long) :: var_index
1153 
1154 #ifdef ENABLE_ACTION_MSGS
1155  character(len=1000) :: action_str
1156 
1157  if (nclayer_enable_action) then
1158  write(action_str, "(A, F0.5, A)") "nc_diag_metadata_rdouble(metadata_name = " // metadata_name // ", metadata_value = ", metadata_value, ")"
1159  call nclayer_actionm(trim(action_str))
1160  end if
1161 #endif
1162 
1163  if (diag_metadata_store%data_lock) then
1164  call nclayer_error("Can't add new data - data have already been written and locked!")
1165  end if
1166 
1167  var_index = nc_diag_metadata_lookup_var(metadata_name)
1168 
1169  if (var_index == -1) then
1170  ! First, check to make sure we can still define new variables.
1171  if (diag_metadata_store%def_lock) then
1172  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1173  end if
1174 
1176 
1177  diag_metadata_store%total = diag_metadata_store%total + 1
1178 
1179  diag_metadata_store%names(diag_metadata_store%total) = metadata_name
1181 
1182  var_index = diag_metadata_store%total
1183  end if
1184 
1185  ! We just need to add one entry...
1186  call nc_diag_metadata_resize_iarr(var_index, 1_i_llong)
1187  call nc_diag_metadata_resize_rdouble(1_i_llong)
1188 
1189  ! Now add the actual entry!
1190  diag_metadata_store%m_rdouble(diag_metadata_store%acount(5)) = metadata_value
1191  diag_metadata_store%stor_i_arr(var_index)%index_arr(diag_metadata_store%stor_i_arr(var_index)%icount) = &
1192  diag_metadata_store%acount(5)
1193  end subroutine nc_diag_metadata_rdouble
1194 
1195  ! nc_diag_metadata - input character(len=*)
1196  ! Corresponding NetCDF4 type: string? char?
1197  subroutine nc_diag_metadata_string(metadata_name, metadata_value)
1198  character(len=*), intent(in) :: metadata_name
1199  character(len=*), intent(in) :: metadata_value
1200 
1201  integer(i_long) :: var_index
1202 
1203 #ifdef ENABLE_ACTION_MSGS
1204  character(len=1000) :: action_str
1205 
1206  if (nclayer_enable_action) then
1207  write(action_str, "(A)") "nc_diag_metadata_string(metadata_name = " // metadata_name // ", metadata_value = " // trim(metadata_value) // ")"
1208  call nclayer_actionm(trim(action_str))
1209  end if
1210 #endif
1211 
1212  if (diag_metadata_store%data_lock) then
1213  call nclayer_error("Can't add new data - data have already been written and locked!")
1214  end if
1215 
1216  var_index = nc_diag_metadata_lookup_var(metadata_name)
1217 
1218  if (var_index == -1) then
1219  ! First, check to make sure we can still define new variables.
1220  if (diag_metadata_store%def_lock) then
1221  call nclayer_error("Can't add new variable - definitions have already been written and locked!")
1222  end if
1223 
1225 
1226  diag_metadata_store%total = diag_metadata_store%total + 1
1227 
1228  diag_metadata_store%names(diag_metadata_store%total) = metadata_name
1230 
1231  var_index = diag_metadata_store%total
1232  else
1233  ! Check max string length
1234 #ifdef _DEBUG_MEM_
1235  print *, "len_trim(metadata_value) = ", len_trim(metadata_value)
1236  print *, "diag_metadata_store%max_str_lens(var_index) = ", diag_metadata_store%max_str_lens(var_index)
1237 #endif
1238  if ((diag_metadata_store%def_lock) .AND. &
1239  (len_trim(metadata_value) > diag_metadata_store%max_str_lens(var_index))) &
1240  call nclayer_error("Cannot expand variable string length after locking variable definitions!")
1241  end if
1242 
1243  ! We just need to add one entry...
1244  ! Strings can't be vectored (at least for attributes), so no 2nd argument
1245  ! here.
1246  call nc_diag_metadata_resize_iarr(var_index, 1_i_llong)
1247  call nc_diag_metadata_resize_string(1_i_llong)
1248 
1249  ! If trim isn't enabled, set our maximum string length here!
1250  if (.NOT. enable_trim) then
1251  if (diag_metadata_store%max_str_lens(var_index) == -1) then
1252  diag_metadata_store%max_str_lens(var_index) = len(metadata_value)
1253  else
1254  ! Validate that our non-first value isn't different from
1255  ! the initial string length
1256  if (diag_metadata_store%max_str_lens(var_index) /= len(metadata_value)) &
1257  call nclayer_error("Cannot change string size when trimming is disabled!")
1258  end if
1259  end if
1260 
1261  ! Now add the actual entry!
1262  diag_metadata_store%m_string(diag_metadata_store%acount(6)) = metadata_value
1263  diag_metadata_store%stor_i_arr(var_index)%index_arr(diag_metadata_store%stor_i_arr(var_index)%icount) = &
1264  diag_metadata_store%acount(6)
1265  end subroutine nc_diag_metadata_string
1266 end module ncdw_metadata
subroutine nc_diag_metadata_long(metadata_name, metadata_value)
subroutine nc_diag_metadata_resize_rdouble(addl_num_entries, update_acount_in)
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
type(diag_metadata), allocatable diag_metadata_store
Definition: ncdw_state.f90:17
integer, parameter, public i_long
Definition: ncd_kinds.F90:47
subroutine nc_diag_metadata_resize_rsingle(addl_num_entries, update_acount_in)
subroutine nc_diag_metadata_resize_string(addl_num_entries, update_acount_in)
subroutine nc_diag_metadata_prealloc_vars_storage_all(num_of_addl_slots)
subroutine nc_diag_metadata_prealloc_vars_storage(nclayer_type, num_of_addl_slots)
subroutine nc_diag_metadata_allocmulti(multiplier)
logical init_done
Definition: ncdw_state.f90:9
integer(i_byte), parameter nlayer_double
Definition: ncdw_types.F90:14
subroutine nc_diag_metadata_resize_long(addl_num_entries, update_acount_in)
subroutine nc_diag_metadata_set_strict(enable_strict)
logical enable_trim
Definition: ncdw_state.f90:12
subroutine nc_diag_metadata_write_data(flush_data_only)
integer function max_len_string_array(str_arr, arr_length)
type(diag_varattr), allocatable diag_varattr_store
Definition: ncdw_state.f90:19
integer(i_byte), parameter nlayer_string
Definition: ncdw_types.F90:15
integer(i_long) ncid
Definition: ncdw_state.f90:8
subroutine nc_diag_metadata_prealloc_vars(num_of_addl_vars)
subroutine nc_diag_metadata_rdouble(metadata_name, metadata_value)
subroutine nc_diag_metadata_resize_iarr_type(addl_num_entries)
subroutine nc_diag_metadata_resize_byte(addl_num_entries, update_acount_in)
integer(i_long), parameter nlayer_compression
Definition: ncdw_types.F90:21
subroutine nc_diag_metadata_expand
subroutine nc_diag_metadata_short(metadata_name, metadata_value)
logical append_only
Definition: ncdw_state.f90:10
integer function nc_diag_metadata_lookup_var(metadata_name)
integer, parameter, public i_short
Definition: ncd_kinds.F90:46
subroutine nc_diag_metadata_rsingle(metadata_name, metadata_value)
subroutine nc_diag_metadata_byte(metadata_name, metadata_value)
integer(i_byte), parameter nlayer_byte
Definition: ncdw_types.F90:10
subroutine nc_diag_metadata_resize_short(addl_num_entries, update_acount_in)
subroutine nc_diag_metadata_load_def
subroutine nc_diag_varattr_add_var(var_name, var_type, var_id)
integer(i_short), parameter nlayer_default_ent
Definition: ncdw_types.F90:18
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_metadata_write_def(internal)
integer, parameter, public i_llong
Definition: ncd_kinds.F90:49
subroutine nc_diag_metadata_resize_iarr(iarr_index, addl_num_entries, update_icount_in)
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_metadata_string(metadata_name, metadata_value)
subroutine nc_diag_varattr_make_nobs_dim