PartMC  2.6.1
gas_data.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2012, 2021 Nicole Riemer and Matthew West
2 ! Licensed under the GNU General Public License version 2 or (at your
3 ! option) any later version. See the file COPYING for details.
4 
5 !> \file
6 !> The pmc_gas_data module.
7 
8 !> The gas_data_t structure and associated subroutines.
10 
11  use pmc_spec_file
12  use pmc_mpi
13  use pmc_util
14  use pmc_netcdf
15 #ifdef PMC_USE_CAMP
16  use camp_camp_core
17  use camp_chem_spec_data
18  use camp_property
19  use camp_util, only: string_t, split_string, split_char
20 #endif
21 #ifdef PMC_USE_MPI
22  use mpi
23 #endif
24 
25  !> Maximum length of the name of a gas.
26  integer, parameter :: gas_name_len = 100
27 
28  !> Constant gas data.
29  !!
30  !! Each gas species is identified by an integer \c i between 1 and
31  !! \c gas_data_n_spec(gas_data). Species \c i has name \c gas_data%%name(i).
32  !! The variable gas data describing the current mixing ratios is stored
33  !! in the gas_state_t structure, so the mixing ratio of species \c i
34  !! is gas_state%%mix_rat(i).
36  !> Species name [length \c gas_data_n_spec(gas_data)].
37  character(len=GAS_NAME_LEN), allocatable :: name(:)
38  !> Index of the corresponding MOSAIC species [length \c
39  !> gas_data_n_spec(gas_data)]. \c to_mosaic(i) is the mosaic index of
40  !> species \c i, or 0 if there is no match.
41  integer, allocatable :: mosaic_index(:)
42  !> Index for gas-phase water in the CAMP state array
43  integer :: i_camp_water = 0
44  end type gas_data_t
45 
46 contains
47 
48 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
49 
50 #ifdef PMC_USE_CAMP
51  !> Initialize the gas_data_t instance from camp_core data
52  subroutine gas_data_initialize(gas_data, camp_core)
53 
54  !> Gas-phase species data
55  class(gas_data_t), intent(inout) :: gas_data
56  !> CAMP core
57  class(camp_core_t), intent(in) :: camp_core
58 
59  type(chem_spec_data_t), pointer :: chem_spec_data
60  integer :: i_spec
61  type(string_t), allocatable :: gas_spec_names(:)
62  type(property_t), pointer :: property_set
63  character(len=:), allocatable :: prop_name
64  logical :: bool_val
65 
66  ! Get the chemical species data
67  call assert_msg(139566827, &
68  camp_core%get_chem_spec_data(chem_spec_data), &
69  "No chemical species data in camp core.")
70 
71  ! Get the gas-phase species names
72  gas_spec_names = chem_spec_data%get_spec_names( &
73  spec_phase = chem_spec_gas_phase)
74 
75  ! Allocate space for the gas-phase species
76  allocate(gas_data%name(size(gas_spec_names)))
77 
78  ! Set the species names and locate gas-phase water
79  prop_name = "is gas-phase water"
80  do i_spec = 1, size(gas_spec_names)
81  gas_data%name(i_spec) = gas_spec_names(i_spec)%string
82  call assert_msg(990037352, &
83  chem_spec_data%get_property_set( &
84  gas_spec_names(i_spec)%string, &
85  property_set), &
86  "Missing property set for gas species "// &
87  gas_spec_names(i_spec)%string)
88  if (property_set%get_logical(prop_name, bool_val)) then
89  call assert_msg(423633615, gas_data%i_camp_water == 0, &
90  "More than one gas-phase water species specified")
91  gas_data%i_camp_water = i_spec
92  end if
93  end do
94 
95  call assert_msg(134440820, gas_data%i_camp_water /= 0, &
96  "No gas-phase water species specified.")
97 
98  ! Allocate the mosaic index array and set to zero
99  allocate(gas_data%mosaic_index(size(gas_spec_names)))
100  gas_data%mosaic_index(:) = 0
101 
102  end subroutine gas_data_initialize
103 #endif
104 
105 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
106 
107  !> Return the number of gas species.
108  elemental integer function gas_data_n_spec(gas_data)
109 
110  !> Aero data structure.
111  type(gas_data_t), intent(in) :: gas_data
112 
113  if (allocated(gas_data%name)) then
114  gas_data_n_spec = size(gas_data%name)
115  else
116  gas_data_n_spec = 0
117  end if
118 
119  end function gas_data_n_spec
120 
121 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
122 
123  !> Returns the number of the species in gas with the given name, or
124  !> returns 0 if there is no such species.
125  integer function gas_data_spec_by_name(gas_data, name)
126 
127  !> Gas data.
128  type(gas_data_t), intent(in) :: gas_data
129  !> Name of species to find.
130  character(len=*), intent(in) :: name
131 
132  integer i
133  logical found
134 
135  found = .false.
136  do i = 1,gas_data_n_spec(gas_data)
137  if (name == gas_data%name(i)) then
138  found = .true.
139  exit
140  end if
141  end do
142  if (found) then
144  else
146  end if
147 
148  end function gas_data_spec_by_name
149 
150 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
151 
152  !> Fills in gas_data%%mosaic_index.
153  subroutine gas_data_set_mosaic_map(gas_data)
154 
155  !> Gas data.
156  type(gas_data_t), intent(inout) :: gas_data
157 
158  integer, parameter :: n_mosaic_species = 77
159  character(len=GAS_NAME_LEN), parameter, dimension(n_mosaic_species) &
160  :: mosaic_species = [ &
161  "H2SO4 ", "HNO3 ", "HCl ", &
162  "NH3 ", "NO ", "NO2 ", &
163  "NO3 ", "N2O5 ", "HONO ", &
164  "HNO4 ", "O3 ", "O1D ", &
165  "O3P ", "OH ", "HO2 ", &
166  "H2O2 ", "CO ", "SO2 ", &
167  "CH4 ", "C2H6 ", "CH3O2 ", &
168  "ETHP ", "HCHO ", "CH3OH ", &
169  "ANOL ", "CH3OOH ", "ETHOOH ", &
170  "ALD2 ", "HCOOH ", "RCOOH ", &
171  "C2O3 ", "PAN ", "ARO1 ", &
172  "ARO2 ", "ALK1 ", "OLE1 ", &
173  "API1 ", "API2 ", "LIM1 ", &
174  "LIM2 ", "PAR ", "AONE ", &
175  "MGLY ", "ETH ", "OLET ", &
176  "OLEI ", "TOL ", "XYL ", &
177  "CRES ", "TO2 ", "CRO ", &
178  "OPEN ", "ONIT ", "ROOH ", &
179  "RO2 ", "ANO2 ", "NAP ", &
180  "XO2 ", "XPAR ", "ISOP ", &
181  "ISOPRD ", "ISOPP ", "ISOPN ", &
182  "ISOPO2 ", "API ", "LIM ", &
183  "DMS ", "MSA ", "DMSO ", &
184  "DMSO2 ", "CH3SO2H ", "CH3SCH2OO ", &
185  "CH3SO2 ", "CH3SO3 ", "CH3SO2OO ", &
186  "CH3SO2CH2OO ", "SULFHOX "]
187 
188  integer spec, mosaic_spec, i
189 
190  gas_data%mosaic_index = 0
191  do spec = 1,gas_data_n_spec(gas_data)
192  mosaic_spec = 0
193  do i = 1,n_mosaic_species
194  if (gas_data%name(spec) == mosaic_species(i)) then
195  mosaic_spec = i
196  end if
197  end do
198  if (mosaic_spec > 0) then
199  gas_data%mosaic_index(spec) = mosaic_spec
200  end if
201  end do
202 
203  end subroutine gas_data_set_mosaic_map
204 
205 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
206 
207  !> Read gas data from a .spec file.
208  subroutine spec_file_read_gas_data(file, gas_data)
209 
210  !> Spec file to read data from.
211  type(spec_file_t), intent(inout) :: file
212  !> Gas data.
213  type(gas_data_t), intent(inout) :: gas_data
214 
215  integer :: n_species, species, i
216  character(len=SPEC_LINE_MAX_VAR_LEN), allocatable :: species_name(:)
217  real(kind=dp), allocatable :: species_data(:,:)
218 
219  !> \page input_format_gas_data Input File Format: Gas Material Data
220  !!
221  !! A gas material data file must consist of one line per gas
222  !! species, with each line having the species name. This specifies
223  !! which species are to be recognized as gases. For example, a \c
224  !! gas_data file could contain:
225  !! <pre>
226  !! H2SO4
227  !! HNO3
228  !! HCl
229  !! NH3
230  !! </pre>
231  !!
232  !! See also:
233  !! - \ref spec_file_format --- the input file text format
234  !! - \ref output_format_gas_data --- the corresponding output format
235 
236  ! read the gas data from the specified file
237  call spec_file_read_real_named_array(file, 0, species_name, &
238  species_data)
239 
240  ! check the data size
241  if (size(species_data, 2) /= 0) then
242  call die_msg(614290516, 'each line in ' // trim(file%name) &
243  // ' must only contain the species name')
244  end if
245 
246  ! allocate and copy over the data
247  n_species = size(species_data, 1)
248  call ensure_string_array_size(gas_data%name, n_species)
249  do i = 1,n_species
250  gas_data%name(i) = species_name(i)(1:gas_name_len)
251  end do
252 
253  call ensure_integer_array_size(gas_data%mosaic_index, n_species)
254  call gas_data_set_mosaic_map(gas_data)
255 
256  end subroutine spec_file_read_gas_data
257 
258 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
259 
260  !> Determines the number of bytes required to pack the given value.
261  integer function pmc_mpi_pack_size_gas_data(val)
262 
263  !> Value to pack.
264  type(gas_data_t), intent(in) :: val
265 
268  + pmc_mpi_pack_size_integer_array(val%mosaic_index)
269 
270  end function pmc_mpi_pack_size_gas_data
271 
272 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
273 
274  !> Packs the given value into the buffer, advancing position.
275  subroutine pmc_mpi_pack_gas_data(buffer, position, val)
276 
277  !> Memory buffer.
278  character, intent(inout) :: buffer(:)
279  !> Current buffer position.
280  integer, intent(inout) :: position
281  !> Value to pack.
282  type(gas_data_t), intent(in) :: val
283 
284 #ifdef PMC_USE_MPI
285  integer :: prev_position
286 
287  prev_position = position
288  call pmc_mpi_pack_string_array(buffer, position, val%name)
289  call pmc_mpi_pack_integer_array(buffer, position, val%mosaic_index)
290  call assert(449872094, &
291  position - prev_position <= pmc_mpi_pack_size_gas_data(val))
292 #endif
293 
294  end subroutine pmc_mpi_pack_gas_data
295 
296 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
297 
298  !> Unpacks the given value from the buffer, advancing position.
299  subroutine pmc_mpi_unpack_gas_data(buffer, position, val)
300 
301  !> Memory buffer.
302  character, intent(inout) :: buffer(:)
303  !> Current buffer position.
304  integer, intent(inout) :: position
305  !> Value to pack.
306  type(gas_data_t), intent(inout) :: val
307 
308 #ifdef PMC_USE_MPI
309  integer :: prev_position
310 
311  prev_position = position
312  call pmc_mpi_unpack_string_array(buffer, position, val%name)
313  call pmc_mpi_unpack_integer_array(buffer, position, val%mosaic_index)
314  call assert(137879163, &
315  position - prev_position <= pmc_mpi_pack_size_gas_data(val))
316 #endif
317 
318  end subroutine pmc_mpi_unpack_gas_data
319 
320 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
321 
322  !> Write the gas species dimension to the given NetCDF file if it
323  !> is not already present and in any case return the associated
324  !> dimid.
325  subroutine gas_data_netcdf_dim_gas_species(gas_data, ncid, &
326  dimid_gas_species)
327 
328  !> Gas_data structure.
329  type(gas_data_t), intent(in) :: gas_data
330  !> NetCDF file ID, in data mode.
331  integer, intent(in) :: ncid
332  !> Dimid of the species dimension.
333  integer, intent(out) :: dimid_gas_species
334 
335  integer :: status, i_spec
336  integer :: varid_gas_species
337  integer :: gas_species_centers(gas_data_n_spec(gas_data))
338  character(len=((GAS_NAME_LEN + 2) * gas_data_n_spec(gas_data))) &
339  :: gas_species_names
340 
341  ! try to get the dimension ID
342  status = nf90_inq_dimid(ncid, "gas_species", dimid_gas_species)
343  if (status == nf90_noerr) return
344  if (status /= nf90_ebaddim) call pmc_nc_check(status)
345 
346  ! dimension not defined, so define now define it
347  call pmc_nc_check(nf90_redef(ncid))
348 
349  call pmc_nc_check(nf90_def_dim(ncid, "gas_species", &
350  gas_data_n_spec(gas_data), dimid_gas_species))
351  gas_species_names = ""
352  do i_spec = 1,gas_data_n_spec(gas_data)
353  gas_species_names((len_trim(gas_species_names) + 1):) &
354  = trim(gas_data%name(i_spec))
355  if (i_spec < gas_data_n_spec(gas_data)) then
356  gas_species_names((len_trim(gas_species_names) + 1):) = ","
357  end if
358  end do
359  call pmc_nc_check(nf90_def_var(ncid, "gas_species", nf90_int, &
360  dimid_gas_species, varid_gas_species))
361  call pmc_nc_check(nf90_put_att(ncid, varid_gas_species, "names", &
362  gas_species_names))
363  call pmc_nc_check(nf90_put_att(ncid, varid_gas_species, "description", &
364  "dummy dimension variable (no useful value) - read species names " &
365  // "as comma-separated values from the 'names' attribute"))
366 
367  call pmc_nc_check(nf90_enddef(ncid))
368 
369  do i_spec = 1,gas_data_n_spec(gas_data)
370  gas_species_centers(i_spec) = i_spec
371  end do
372  call pmc_nc_check(nf90_put_var(ncid, varid_gas_species, &
373  gas_species_centers))
374 
375  end subroutine gas_data_netcdf_dim_gas_species
376 
377 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
378 
379  !> Write full state.
380  subroutine gas_data_output_netcdf(gas_data, ncid)
381 
382  !> Gas_data to write.
383  type(gas_data_t), intent(in) :: gas_data
384  !> NetCDF file ID, in data mode.
385  integer, intent(in) :: ncid
386 
387  integer :: dimid_gas_species
388 
389  !> \page output_format_gas_data Output File Format: Gas Material Data
390  !!
391  !! The gas material data NetCDF dimensions are:
392  !! - \b gas_species: number of gas species
393  !!
394  !! The gas material data NetCDF variables are:
395  !! - \b gas_species (dim \c gas_species): dummy dimension variable
396  !! (no useful value) - read species names as comma-separated values
397  !! from the 'names' attribute
398  !! - \b gas_mosaic_index (dim \c gas_species): MOSAIC indices of
399  !! gas species
400  !!
401  !! See also:
402  !! - \ref input_format_gas_data --- the corresponding input format
403 
404  call gas_data_netcdf_dim_gas_species(gas_data, ncid, &
405  dimid_gas_species)
406 
407  call pmc_nc_write_integer_1d(ncid, gas_data%mosaic_index, &
408  "gas_mosaic_index", (/ dimid_gas_species /), &
409  long_name="MOSAIC indices of gas species")
410 
411  end subroutine gas_data_output_netcdf
412 
413 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
414 
415  !> Read full state.
416  subroutine gas_data_input_netcdf(gas_data, ncid)
417 
418  !> Gas_data to read.
419  type(gas_data_t), intent(inout) :: gas_data
420  !> NetCDF file ID, in data mode.
421  integer, intent(in) :: ncid
422 
423  character(len=1000) :: name
424  integer :: dimid_gas_species, n_spec, varid_gas_species, i_spec, i
425  character(len=:), allocatable :: gas_species_names
426 
427  call pmc_nc_check(nf90_inq_dimid(ncid, "gas_species", dimid_gas_species))
428  call pmc_nc_check(nf90_inquire_dimension(ncid, dimid_gas_species, name, &
429  n_spec))
430  call assert(719237193, n_spec < 1000)
431 
432  if (allocated(gas_data%name)) deallocate(gas_data%name)
433  if (allocated(gas_data%mosaic_index)) deallocate(gas_data%mosaic_index)
434  allocate(gas_data%name(n_spec))
435  allocate(gas_data%mosaic_index(n_spec))
436 
437  call pmc_nc_read_integer_1d(ncid, gas_data%mosaic_index, &
438  "gas_mosaic_index")
439 
440  call pmc_nc_check(nf90_inq_varid(ncid, "gas_species", varid_gas_species))
441  allocate(character(len=((GAS_NAME_LEN + 2) * 1000)) :: gas_species_names)
442  call pmc_nc_check(nf90_get_att(ncid, varid_gas_species, "names", &
443  gas_species_names))
444  ! gas_species_names are comma-separated, so unpack them
445  do i_spec = 1,gas_data_n_spec(gas_data)
446  i = 1
447  do while ((gas_species_names(i:i) /= " ") &
448  .and. (gas_species_names(i:i) /= ","))
449  i = i + 1
450  end do
451  call assert(173021381, i > 1)
452  gas_data%name(i_spec) = gas_species_names(1:(i-1))
453  gas_species_names = gas_species_names((i+1):)
454  end do
455  call assert(469721220, gas_species_names == "")
456 
457  end subroutine gas_data_input_netcdf
458 
459 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
460 
461 end module pmc_gas_data
pmc_mpi::pmc_mpi_pack_string_array
subroutine pmc_mpi_pack_string_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:752
pmc_netcdf::pmc_nc_read_integer_1d
subroutine pmc_nc_read_integer_1d(ncid, var, name, must_be_present)
Read a simple integer array from a NetCDF file.
Definition: netcdf.F90:260
pmc_gas_data::gas_data_t
Constant gas data.
Definition: gas_data.F90:35
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_mpi::pmc_mpi_unpack_string_array
subroutine pmc_mpi_unpack_string_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1010
pmc_mpi::pmc_mpi_unpack_integer_array
subroutine pmc_mpi_unpack_integer_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:946
pmc_gas_data::gas_data_netcdf_dim_gas_species
subroutine gas_data_netcdf_dim_gas_species(gas_data, ncid, dimid_gas_species)
Write the gas species dimension to the given NetCDF file if it is not already present and in any case...
Definition: gas_data.F90:327
pmc_mpi::pmc_mpi_pack_integer_array
subroutine pmc_mpi_pack_integer_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:688
pmc_gas_data
The gas_data_t structure and associated subroutines.
Definition: gas_data.F90:9
pmc_netcdf::pmc_nc_write_integer_1d
subroutine pmc_nc_write_integer_1d(ncid, var, name, dimids, dim_name, unit, long_name, standard_name, description)
Write a simple integer array to a NetCDF file.
Definition: netcdf.F90:585
pmc_util::die_msg
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:134
pmc_netcdf
Wrapper functions for NetCDF. These all take a NetCDF ncid in data mode and return with it again in d...
Definition: netcdf.F90:11
pmc_constants::dp
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
pmc_spec_file
Reading formatted text input.
Definition: spec_file.F90:43
pmc_gas_data::pmc_mpi_pack_gas_data
subroutine pmc_mpi_pack_gas_data(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: gas_data.F90:276
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:103
pmc_spec_file::spec_file_t
An input file with extra data for printing messages.
Definition: spec_file.F90:59
pmc_gas_data::gas_name_len
integer, parameter gas_name_len
Maximum length of the name of a gas.
Definition: gas_data.F90:26
pmc_util::assert_msg
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:77
pmc_spec_file::spec_file_read_real_named_array
subroutine spec_file_read_real_named_array(file, max_lines, names, vals)
Read an array of named lines with real data. All lines must have the same number of data elements.
Definition: spec_file.F90:650
pmc_util::ensure_string_array_size
subroutine ensure_string_array_size(x, n)
Allocate or reallocate the given array to ensure it is of the given size, without preserving data.
Definition: util.F90:1189
pmc_gas_data::gas_data_input_netcdf
subroutine gas_data_input_netcdf(gas_data, ncid)
Read full state.
Definition: gas_data.F90:417
pmc_gas_data::gas_data_spec_by_name
integer function gas_data_spec_by_name(gas_data, name)
Returns the number of the species in gas with the given name, or returns 0 if there is no such specie...
Definition: gas_data.F90:126
pmc_mpi::pmc_mpi_pack_size_integer_array
integer function pmc_mpi_pack_size_integer_array(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:447
pmc_gas_data::pmc_mpi_pack_size_gas_data
integer function pmc_mpi_pack_size_gas_data(val)
Determines the number of bytes required to pack the given value.
Definition: gas_data.F90:262
pmc_netcdf::pmc_nc_check
subroutine pmc_nc_check(status)
Check the status of a NetCDF function call.
Definition: netcdf.F90:23
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_gas_data::gas_data_n_spec
elemental integer function gas_data_n_spec(gas_data)
Return the number of gas species.
Definition: gas_data.F90:109
pmc_gas_data::pmc_mpi_unpack_gas_data
subroutine pmc_mpi_unpack_gas_data(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: gas_data.F90:300
pmc_mpi::pmc_mpi_pack_size_string_array
integer function pmc_mpi_pack_size_string_array(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:507
pmc_gas_data::gas_data_set_mosaic_map
subroutine gas_data_set_mosaic_map(gas_data)
Fills in gas_data%mosaic_index.
Definition: gas_data.F90:154
pmc_util::ensure_integer_array_size
subroutine ensure_integer_array_size(x, n, only_grow)
Allocate or reallocate the given array to ensure it is of the given size, preserving any data and/or ...
Definition: util.F90:1115