PartMC  2.3.0
bin_grid.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2015 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_bin_grid module.
7 
8 !> The bin_grid_t structure and associated subroutines.
10 
11  use pmc_constants
12  use pmc_util
13  use pmc_spec_file
14  use pmc_mpi
15  use pmc_netcdf
16 #ifdef PMC_USE_MPI
17  use mpi
18 #endif
19 
20  !> Invalid type of bin grid.
21  integer, parameter :: BIN_GRID_TYPE_INVALID = 0
22  !> Logarithmically spaced bin grid.
23  integer, parameter :: BIN_GRID_TYPE_LOG = 1
24  !> Linearly spaced bin grid.
25  integer, parameter :: BIN_GRID_TYPE_LINEAR = 2
26 
27  !> 1D grid, either logarithmic or linear.
28  !!
29  !! The grid of bins is logarithmically spaced in volume, an
30  !! assumption that is quite heavily incorporated into the code. At
31  !! some point in the future it would be nice to relax this
32  !! assumption.
34  !> Type of grid spacing (BIN_GRID_TYPE_LOG, etc).
35  integer :: type
36  !> Number of bins.
37  integer :: n_bin
38  !> Bin centers.
39  real(kind=dp), allocatable :: centers(:)
40  !> Bin edges.
41  real(kind=dp), allocatable :: edges(:)
42  !> Bin widths.
43  real(kind=dp), allocatable :: widths(:)
44  end type bin_grid_t
45 
46 contains
47 
48 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
49 
50  !> Allocates a bin_grid.
51  subroutine bin_grid_allocate(bin_grid)
52 
53  !> Bin grid.
54  type(bin_grid_t), intent(out) :: bin_grid
55 
56  bin_grid%type = bin_grid_type_invalid
57  bin_grid%n_bin = 0
58  allocate(bin_grid%centers(0))
59  allocate(bin_grid%edges(0))
60  allocate(bin_grid%widths(0))
61 
62  end subroutine bin_grid_allocate
63 
64 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
65 
66  !> Allocates a bin_grid.
67  subroutine bin_grid_allocate_size(bin_grid, n_bin)
68 
69  !> Bin grid.
70  type(bin_grid_t), intent(out) :: bin_grid
71  !> Number of bins.
72  integer, intent(in) :: n_bin
73 
74  bin_grid%type = bin_grid_type_invalid
75  bin_grid%n_bin = n_bin
76  allocate(bin_grid%centers(n_bin))
77  allocate(bin_grid%edges(n_bin + 1))
78  allocate(bin_grid%widths(n_bin))
79 
80  end subroutine bin_grid_allocate_size
81 
82 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
83 
84  !> Frees all memory.
85  subroutine bin_grid_deallocate(bin_grid)
86 
87  !> Bin_grid to free.
88  type(bin_grid_t), intent(inout) :: bin_grid
89 
90  deallocate(bin_grid%centers)
91  deallocate(bin_grid%edges)
92  deallocate(bin_grid%widths)
93 
94  end subroutine bin_grid_deallocate
95 
96 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
97 
98  !> Copies a bin grid.
99  subroutine bin_grid_copy(bin_grid_from, bin_grid_to)
100 
101  !> Bin_grid to copy from.
102  type(bin_grid_t), intent(in) :: bin_grid_from
103  !> Bin_grid to copy to.
104  type(bin_grid_t), intent(inout) :: bin_grid_to
105 
106  call bin_grid_deallocate(bin_grid_to)
107  call bin_grid_allocate_size(bin_grid_to, bin_grid_from%n_bin)
108  bin_grid_to%type = bin_grid_from%type
109  bin_grid_to%n_bin = bin_grid_from%n_bin
110  bin_grid_to%centers = bin_grid_from%centers
111  bin_grid_to%edges = bin_grid_from%edges
112  bin_grid_to%widths = bin_grid_from%widths
113 
114  end subroutine bin_grid_copy
115 
116 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
117 
118  !> Convert a concentration f(vol)d(vol) to f(ln(r))d(ln(r))
119  !> where vol = 4/3 pi r^3.
120  subroutine vol_to_lnr(r, f_vol, f_lnr)
121 
122  !> Radius (m).
123  real(kind=dp), intent(in) :: r
124  !> Concentration as a function of volume.
125  real(kind=dp), intent(in) :: f_vol
126  !> Concentration as a function of ln(r).
127  real(kind=dp), intent(out) :: f_lnr
128 
129  f_lnr = f_vol * 4d0 * const%pi * r**3
130 
131  end subroutine vol_to_lnr
132 
133 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
134 
135  !> Generates the bin grid given the range and number of bins.
136  subroutine bin_grid_make(bin_grid, type, n_bin, min, max)
137 
138  !> New bin grid.
139  type(bin_grid_t), intent(inout) :: bin_grid
140  !> Type of bin grid.
141  integer, intent(in) :: type
142  !> Number of bins.
143  integer, intent(in) :: n_bin
144  !> Minimum edge value.
145  real(kind=dp), intent(in) :: min
146  !> Minimum edge value.
147  real(kind=dp), intent(in) :: max
148 
149  real(kind=dp) :: c1, c2
150 
151  call assert_msg(538534122, n_bin >= 0, &
152  "bin_grid requires a non-negative n_bin, not: " &
153  // trim(integer_to_string(n_bin)))
154  if (n_bin > 0) then
155  if (type == bin_grid_type_log) then
156  call assert_msg(966541762, min > 0d0, &
157  "log bin_grid requires a positive min value, not: " &
158  // trim(real_to_string(min)))
159  end if
160  call assert_msg(711537859, min < max, &
161  "bin_grid requires min < max, not: " &
162  // trim(real_to_string(min)) // " and " &
163  // trim(real_to_string(max)))
164  end if
165  call bin_grid_deallocate(bin_grid)
166  call bin_grid_allocate_size(bin_grid, n_bin)
167  bin_grid%type = type
168  if (n_bin == 0) return
169  if (type == bin_grid_type_log) then
170  call logspace(min, max, bin_grid%edges)
171  c1 = exp(interp_linear_disc(log(min), log(max), 2 * n_bin + 1, 2))
172  c2 = exp(interp_linear_disc(log(min), log(max), 2 * n_bin + 1, &
173  2 * n_bin))
174  call logspace(c1, c2, bin_grid%centers)
175  bin_grid%widths = (log(max) - log(min)) / real(n_bin, kind=dp)
176  elseif (bin_grid%type == bin_grid_type_linear) then
177  call linspace(min, max, bin_grid%edges)
178  c1 = interp_linear_disc(min, max, 2 * n_bin + 1, 2)
179  c2 = interp_linear_disc(min, max, 2 * n_bin + 1, 2 * n_bin)
180  call linspace(c1, c2, bin_grid%centers)
181  bin_grid%widths = (max - min) / real(n_bin, kind=dp)
182  else
183  call die_msg(678302366, "unknown bin_grid type: " &
184  // trim(integer_to_string(bin_grid%type)))
185  end if
186 
187  end subroutine bin_grid_make
188 
189 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
190 
191  !> Find the bin number that contains a given value.
192  !!
193  !! If a particle is below the smallest bin then its bin number is
194  !! 0. If a particle is above the largest bin then its bin number is
195  !! <tt>n_bin + 1</tt>.
196  integer function bin_grid_find(bin_grid, val)
197 
198  !> Bin_grid.
199  type(bin_grid_t), intent(in) :: bin_grid
200  !> Value to locate bin for.
201  real(kind=dp), intent(in) :: val
202 
203  call assert(448215689, bin_grid%n_bin >= 0)
204  bin_grid_find = 0
205  if (bin_grid%n_bin == 0) return
206  if (bin_grid%type == bin_grid_type_log) then
207  bin_grid_find = logspace_find(bin_grid%edges(1), &
208  bin_grid%edges(bin_grid%n_bin + 1), bin_grid%n_bin + 1, val)
209  elseif (bin_grid%type == bin_grid_type_linear) then
210  bin_grid_find = linspace_find(bin_grid%edges(1), &
211  bin_grid%edges(bin_grid%n_bin + 1), bin_grid%n_bin + 1, val)
212  else
213  call die_msg(348908641, "unknown bin_grid type: " &
214  // trim(integer_to_string(bin_grid%type)))
215  end if
216 
217  end function bin_grid_find
218 
219 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
220 
221  !> Make a histogram with of the given weighted data, scaled by the
222  !> bin sizes.
223  function bin_grid_histogram_1d(x_bin_grid, x_data, weight_data)
224 
225  !> x-axis bin grid.
226  type(bin_grid_t), intent(in) :: x_bin_grid
227  !> Data values on the x-axis.
228  real(kind=dp), intent(in) :: x_data(:)
229  !> Data value weights.
230  real(kind=dp), intent(in) :: weight_data(size(x_data))
231 
232  !> Return histogram.
233  real(kind=dp) :: bin_grid_histogram_1d(x_bin_grid%n_bin)
234 
235  integer :: i_data, x_bin
236 
238  do i_data = 1,size(x_data)
239  x_bin = bin_grid_find(x_bin_grid, x_data(i_data))
240  if ((x_bin >= 1) .and. (x_bin <= x_bin_grid%n_bin)) then
242  + weight_data(i_data) / x_bin_grid%widths(x_bin)
243  end if
244  end do
245 
246  end function bin_grid_histogram_1d
247 
248 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
249 
250  !> Make a 2D histogram with of the given weighted data, scaled by
251  !> the bin sizes.
252  function bin_grid_histogram_2d(x_bin_grid, x_data, y_bin_grid, y_data, &
253  weight_data)
254 
255  !> x-axis bin grid.
256  type(bin_grid_t), intent(in) :: x_bin_grid
257  !> Data values on the x-axis.
258  real(kind=dp), intent(in) :: x_data(:)
259  !> y-axis bin grid.
260  type(bin_grid_t), intent(in) :: y_bin_grid
261  !> Data values on the y-axis.
262  real(kind=dp), intent(in) :: y_data(size(x_data))
263  !> Data value weights.
264  real(kind=dp), intent(in) :: weight_data(size(x_data))
265 
266  !> Return histogram.
267  real(kind=dp) :: bin_grid_histogram_2d(x_bin_grid%n_bin, y_bin_grid%n_bin)
268 
269  integer :: i_data, x_bin, y_bin
270 
272  do i_data = 1,size(x_data)
273  x_bin = bin_grid_find(x_bin_grid, x_data(i_data))
274  y_bin = bin_grid_find(y_bin_grid, y_data(i_data))
275  if ((x_bin >= 1) .and. (x_bin <= x_bin_grid%n_bin) &
276  .and. (y_bin >= 1) .and. (y_bin <= y_bin_grid%n_bin)) then
277  bin_grid_histogram_2d(x_bin, y_bin) &
278  = bin_grid_histogram_2d(x_bin, y_bin) + weight_data(i_data) &
279  / x_bin_grid%widths(x_bin) / y_bin_grid%widths(y_bin)
280  end if
281  end do
282 
283  end function bin_grid_histogram_2d
284 
285 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
286 
287  !> Read the specification for a radius bin_grid from a spec file.
288  subroutine spec_file_read_radius_bin_grid(file, bin_grid)
289 
290  !> Spec file.
291  type(spec_file_t), intent(inout) :: file
292  !> Radius bin grid.
293  type(bin_grid_t), intent(inout) :: bin_grid
294 
295  integer :: n_bin
296  real(kind=dp) :: d_min, d_max
297 
298  !> \page input_format_diam_bin_grid Input File Format: Diameter Axis Bin Grid
299  !!
300  !! The diameter bin grid is logarithmic, consisting of
301  !! \f$n_{\rm bin}\f$ bins with centers \f$c_i\f$ (\f$i =
302  !! 1,\ldots,n_{\rm bin}\f$) and edges \f$e_i\f$ (\f$i =
303  !! 1,\ldots,(n_{\rm bin} + 1)\f$) such that \f$e_{i+1}/e_i\f$ is a
304  !! constant and \f$c_i/e_i = \sqrt{e_{i+1}/e_i}\f$. That is,
305  !! \f$\ln(e_i)\f$ are uniformly spaced and \f$\ln(c_i)\f$ are the
306  !! arithmetic centers.
307  !!
308  !! The diameter axis bin grid is specified by the parameters:
309  !! - \b n_bin (integer): The number of bins \f$n_{\rm bin}\f$.
310  !! - \b d_min (real, unit m): The left edge of the left-most bin,
311  !! \f$e_1\f$.
312  !! - \b d_max (real, unit m): The right edge of the right-most bin,
313  !! \f$e_{n_{\rm bin} + 1}\f$.
314  !!
315  !! See also:
316  !! - \ref spec_file_format --- the input file text format
317  !! - \ref output_format_diam_bin_grid --- the corresponding output format
318 
319  call spec_file_read_integer(file, 'n_bin', n_bin)
320  call spec_file_read_real(file, 'd_min', d_min)
321  call spec_file_read_real(file, 'd_max', d_max)
322  call bin_grid_make(bin_grid, bin_grid_type_log, n_bin, diam2rad(d_min), &
323  diam2rad(d_max))
324 
325  end subroutine spec_file_read_radius_bin_grid
326 
327 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
328 
329  !> Determines the number of bytes required to pack the given value.
330  integer function pmc_mpi_pack_size_bin_grid(val)
331 
332  !> Value to pack.
333  type(bin_grid_t), intent(in) :: val
334 
336  pmc_mpi_pack_size_integer(val%type) &
337  + pmc_mpi_pack_size_integer(val%n_bin) &
338  + pmc_mpi_pack_size_real_array(val%centers) &
339  + pmc_mpi_pack_size_real_array(val%edges) &
340  + pmc_mpi_pack_size_real_array(val%widths)
341 
342  end function pmc_mpi_pack_size_bin_grid
343 
344 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
345 
346  !> Packs the given value into the buffer, advancing position.
347  subroutine pmc_mpi_pack_bin_grid(buffer, position, val)
348 
349  !> Memory buffer.
350  character, intent(inout) :: buffer(:)
351  !> Current buffer position.
352  integer, intent(inout) :: position
353  !> Value to pack.
354  type(bin_grid_t), intent(in) :: val
355 
356 #ifdef PMC_USE_MPI
357  integer :: prev_position
358 
359  prev_position = position
360  call pmc_mpi_pack_integer(buffer, position, val%type)
361  call pmc_mpi_pack_integer(buffer, position, val%n_bin)
362  call pmc_mpi_pack_real_array(buffer, position, val%centers)
363  call pmc_mpi_pack_real_array(buffer, position, val%edges)
364  call pmc_mpi_pack_real_array(buffer, position, val%widths)
365  call assert(385455586, &
366  position - prev_position <= pmc_mpi_pack_size_bin_grid(val))
367 #endif
368 
369  end subroutine pmc_mpi_pack_bin_grid
370 
371 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
372 
373  !> Unpacks the given value from the buffer, advancing position.
374  subroutine pmc_mpi_unpack_bin_grid(buffer, position, val)
375 
376  !> Memory buffer.
377  character, intent(inout) :: buffer(:)
378  !> Current buffer position.
379  integer, intent(inout) :: position
380  !> Value to pack.
381  type(bin_grid_t), intent(inout) :: val
382 
383 #ifdef PMC_USE_MPI
384  integer :: prev_position
385 
386  prev_position = position
387  call pmc_mpi_unpack_integer(buffer, position, val%type)
388  call pmc_mpi_unpack_integer(buffer, position, val%n_bin)
389  call pmc_mpi_unpack_real_array_alloc(buffer, position, val%centers)
390  call pmc_mpi_unpack_real_array_alloc(buffer, position, val%edges)
391  call pmc_mpi_unpack_real_array_alloc(buffer, position, val%widths)
392  call assert(741838730, &
393  position - prev_position <= pmc_mpi_pack_size_bin_grid(val))
394 #endif
395 
396  end subroutine pmc_mpi_unpack_bin_grid
397 
398 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
399 
400  !> Check whether all processors have the same value.
401  logical function pmc_mpi_allequal_bin_grid(val)
402 
403  !> Value to compare.
404  type(bin_grid_t), intent(inout) :: val
405 
406 #ifdef PMC_USE_MPI
407  if (.not. pmc_mpi_allequal_integer(val%type)) then
408  pmc_mpi_allequal_bin_grid = .false.
409  return
410  end if
411 
412  if (.not. pmc_mpi_allequal_integer(val%n_bin)) then
413  pmc_mpi_allequal_bin_grid = .false.
414  return
415  end if
416 
417  if (val%n_bin == 0) then
419  return
420  end if
421 
422  if (pmc_mpi_allequal_real(val%edges(1)) &
423  .and. pmc_mpi_allequal_real(val%edges(val%n_bin))) then
425  else
426  pmc_mpi_allequal_bin_grid = .false.
427  end if
428 #else
430 #endif
431 
432  end function pmc_mpi_allequal_bin_grid
433 
434 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
435 
436  !> Write a bin grid dimension to the given NetCDF file if it is
437  !> not already present and in any case return the associated dimid.
438  subroutine bin_grid_netcdf_dim(bin_grid, ncid, dim_name, unit, dimid, &
439  long_name, scale)
440 
441  !> Bin_grid structure.
442  type(bin_grid_t), intent(in) :: bin_grid
443  !> NetCDF file ID, in data mode.
444  integer, intent(in) :: ncid
445  !> Dimension name.
446  character(len=*), intent(in) :: dim_name
447  !> Units for the grid.
448  character(len=*), intent(in) :: unit
449  !> Dimid of the grid dimension.
450  integer, intent(out) :: dimid
451  !> Long dimension name to use.
452  character(len=*), intent(in), optional :: long_name
453  !> Factor to scale grid by before output.
454  real(kind=dp), intent(in), optional :: scale
455 
456  integer :: status, varid, dimid_edges, varid_edges, varid_widths, i
457  real(kind=dp) :: centers(bin_grid%n_bin), edges(bin_grid%n_bin + 1)
458  real(kind=dp) :: widths(bin_grid%n_bin)
459  character(len=(len_trim(dim_name)+10)) :: dim_name_edges
460  character(len=255) :: use_long_name
461 
462  status = nf90_inq_dimid(ncid, dim_name, dimid)
463  if (status == nf90_noerr) return
464  if (status /= nf90_ebaddim) call pmc_nc_check(status)
465 
466  ! dimension not defined, so define now define it
467 
468  dim_name_edges = trim(dim_name) // "_edges"
469  if (present(long_name)) then
470  call assert_msg(125084459, len_trim(long_name) <= len(use_long_name), &
471  "long_name is longer than " &
472  // trim(integer_to_string(len(use_long_name))))
473  use_long_name = trim(long_name)
474  else
475  call assert_msg(660927086, len_trim(dim_name) <= len(use_long_name), &
476  "dim_name is longer than " &
477  // trim(integer_to_string(len(use_long_name))))
478  use_long_name = trim(dim_name)
479  end if
480 
481  call pmc_nc_check(nf90_redef(ncid))
482  call pmc_nc_check(nf90_def_dim(ncid, dim_name, bin_grid%n_bin, dimid))
483  call pmc_nc_check(nf90_def_dim(ncid, dim_name_edges, bin_grid%n_bin + 1, &
484  dimid_edges))
485  call pmc_nc_check(nf90_enddef(ncid))
486 
487  centers = bin_grid%centers
488  edges = bin_grid%edges
489  widths = bin_grid%widths
490  if (bin_grid%type == bin_grid_type_log) then
491  if (present(scale)) then
492  centers = centers * scale
493  edges = edges * scale
494  end if
495  call pmc_nc_write_real_1d(ncid, centers, dim_name, (/ dimid /), &
496  unit=unit, long_name=(trim(use_long_name) // " grid centers"), &
497  description=("logarithmically spaced centers of " &
498  // trim(use_long_name) // " grid, so that " // trim(dim_name) &
499  // "(i) is the geometric mean of " // trim(dim_name_edges) &
500  // "(i) and " // trim(dim_name_edges) // "(i + 1)"))
501  call pmc_nc_write_real_1d(ncid, edges, dim_name_edges, &
502  (/ dimid_edges /), unit=unit, &
503  long_name=(trim(use_long_name) // " grid edges"), &
504  description=("logarithmically spaced edges of " &
505  // trim(use_long_name) // " grid, with one more edge than center"))
506  call pmc_nc_write_real_1d(ncid, widths, trim(dim_name) // "_widths", &
507  (/ dimid /), unit="1", &
508  long_name=(trim(use_long_name) // " grid widths"), &
509  description=("base-e logarithmic widths of " &
510  // trim(use_long_name) // " grid, with " // trim(dim_name) &
511  // "_widths(i) = ln(" // trim(dim_name_edges) // "(i + 1) / " &
512  // trim(dim_name_edges) // "(i))"))
513  elseif (bin_grid%type == bin_grid_type_linear) then
514  if (present(scale)) then
515  centers = centers * scale
516  edges = edges * scale
517  widths = widths * scale
518  end if
519  call pmc_nc_write_real_1d(ncid, centers, dim_name, (/ dimid /), &
520  unit=unit, long_name=(trim(use_long_name) // " grid centers"), &
521  description=("linearly spaced centers of " // trim(use_long_name) &
522  // " grid, so that " // trim(dim_name) // "(i) is the mean of " &
523  // trim(dim_name_edges) // "(i) and " // trim(dim_name_edges) &
524  // "(i + 1)"))
525  call pmc_nc_write_real_1d(ncid, edges, dim_name_edges, &
526  (/ dimid_edges /), unit=unit, &
527  long_name=(trim(use_long_name) // " grid edges"), &
528  description=("linearly spaced edges of " &
529  // trim(use_long_name) // " grid, with one more edge than center"))
530  call pmc_nc_write_real_1d(ncid, widths, trim(dim_name) // "_widths", &
531  (/ dimid /), unit=unit, &
532  long_name=(trim(use_long_name) // " grid widths"), &
533  description=("widths of " // trim(use_long_name) &
534  // " grid, with " // trim(dim_name) // "_widths(i) = " &
535  // trim(dim_name_edges) // "(i + 1) - " // trim(dim_name_edges) &
536  // "(i)"))
537  else
538  call die_msg(942560572, "unknown bin_grid type: " &
539  // trim(integer_to_string(bin_grid%type)))
540  end if
541 
542  end subroutine bin_grid_netcdf_dim
543 
544 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
545 
546  !> Write a bin grid to the given NetCDF file.
547  subroutine bin_grid_output_netcdf(bin_grid, ncid, dim_name, unit, &
548  long_name, scale)
549 
550  !> Bin_grid structure.
551  type(bin_grid_t), intent(in) :: bin_grid
552  !> NetCDF file ID, in data mode.
553  integer, intent(in) :: ncid
554  !> Dimension name.
555  character(len=*), intent(in) :: dim_name
556  !> Units for the grid.
557  character(len=*), intent(in) :: unit
558  !> Long dimension name to use.
559  character(len=*), intent(in), optional :: long_name
560  !> Factor to scale grid by before output.
561  real(kind=dp), intent(in), optional :: scale
562 
563  integer :: dimid
564 
565  call bin_grid_netcdf_dim(bin_grid, ncid, dim_name, unit, dimid, &
566  long_name, scale)
567 
568  end subroutine bin_grid_output_netcdf
569 
570 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
571 
572  !> Read full state.
573  subroutine bin_grid_input_netcdf(bin_grid, ncid, dim_name, scale)
574 
575  !> bin_grid to read.
576  type(bin_grid_t), intent(inout) :: bin_grid
577  !> NetCDF file ID, in data mode.
578  integer, intent(in) :: ncid
579  !> Dimension name.
580  character(len=*), intent(in) :: dim_name
581  !> Factor to scale grid by after input.
582  real(kind=dp), intent(in), optional :: scale
583 
584  integer :: dimid, varid, n_bin, type
585  character(len=1000) :: name, description
586  real(kind=dp), allocatable :: edges(:)
587 
588  call pmc_nc_check(nf90_inq_dimid(ncid, dim_name, dimid))
589  call pmc_nc_check(nf90_inquire_dimension(ncid, dimid, name, n_bin))
590  call pmc_nc_check(nf90_inq_varid(ncid, dim_name, varid))
591  call pmc_nc_check(nf90_get_att(ncid, varid, "description", description))
592 
593  allocate(edges(n_bin + 1))
594  call pmc_nc_read_real_1d(ncid, edges, dim_name // "_edges")
595 
596  if (starts_with(description, "logarithmically")) then
597  type = bin_grid_type_log
598  elseif (starts_with(description, "linearly")) then
599  type = bin_grid_type_linear
600  else
601  call die_msg(792158584, "cannot identify grid type for NetCDF " &
602  // "dimension: " // trim(dim_name))
603  end if
604 
605  if (present(scale)) then
606  call bin_grid_make(bin_grid, type, n_bin, scale * edges(1), &
607  scale * edges(n_bin + 1))
608  else
609  call bin_grid_make(bin_grid, type, n_bin, edges(1), edges(n_bin + 1))
610  end if
611 
612  deallocate(edges)
613 
614  end subroutine bin_grid_input_netcdf
615 
616 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
617 
618 end module pmc_bin_grid
real(kind=dp) function, dimension(x_bin_grid%n_bin, y_bin_grid%n_bin) bin_grid_histogram_2d(x_bin_grid, x_data, y_bin_grid, y_data, weight_data)
Make a 2D histogram with of the given weighted data, scaled by the bin sizes.
Definition: bin_grid.F90:252
subroutine pmc_nc_write_real_1d(ncid, var, name, dimids, dim_name, unit, long_name, standard_name, description)
Write a simple real array to a NetCDF file.
Definition: netcdf.F90:504
integer function pmc_mpi_pack_size_real_array(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:471
An input file with extra data for printing messages.
Definition: spec_file.F90:59
integer function pmc_mpi_pack_size_integer(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:347
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:133
subroutine bin_grid_allocate_size(bin_grid, n_bin)
Allocates a bin_grid.
Definition: bin_grid.F90:67
real(kind=dp) function, dimension(x_bin_grid%n_bin) bin_grid_histogram_1d(x_bin_grid, x_data, weight_data)
Make a histogram with of the given weighted data, scaled by the bin sizes.
Definition: bin_grid.F90:223
integer function linspace_find(min_x, max_x, n, x)
Find the position of a real number in a 1D linear array.
Definition: util.F90:517
subroutine bin_grid_output_netcdf(bin_grid, ncid, dim_name, unit, long_name, scale)
Write a bin grid to the given NetCDF file.
Definition: bin_grid.F90:547
logical function starts_with(string, start_string)
Checks whether a string starts with a given other string.
Definition: util.F90:1493
Physical constants.
Definition: constants.F90:9
subroutine pmc_mpi_unpack_real_array_alloc(buffer, position, val)
Unpacks the given value from the buffer to an allocated array, advancing position.
Definition: mpi.F90:1499
subroutine pmc_nc_check(status)
Check the status of a NetCDF function call.
Definition: netcdf.F90:22
subroutine bin_grid_make(bin_grid, type, n_bin, min, max)
Generates the bin grid given the range and number of bins.
Definition: bin_grid.F90:136
integer function bin_grid_find(bin_grid, val)
Find the bin number that contains a given value.
Definition: bin_grid.F90:196
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:76
subroutine logspace(min_x, max_x, x)
Makes a logarithmically spaced array of length n from min to max.
Definition: util.F90:481
logical function pmc_mpi_allequal_real(val)
Returns whether all processors have the same value.
Definition: mpi.F90:1382
subroutine bin_grid_netcdf_dim(bin_grid, ncid, dim_name, unit, dimid, long_name, scale)
Write a bin grid dimension to the given NetCDF file if it is not already present and in any case retu...
Definition: bin_grid.F90:438
real(kind=dp) function interp_linear_disc(x_1, x_n, n, i)
Linear interpolation over discrete indices.
Definition: util.F90:644
subroutine spec_file_read_integer(file, name, var)
Read an integer from a spec file that must have the given name.
Definition: spec_file.F90:555
subroutine bin_grid_allocate(bin_grid)
Allocates a bin_grid.
Definition: bin_grid.F90:51
subroutine pmc_mpi_pack_integer(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:534
subroutine bin_grid_input_netcdf(bin_grid, ncid, dim_name, scale)
Read full state.
Definition: bin_grid.F90:573
real(kind=dp) elemental function diam2rad(d)
Convert diameter (m) to radius (m).
Definition: util.F90:286
Common utility subroutines.
Definition: util.F90:9
integer function pmc_mpi_pack_size_bin_grid(val)
Determines the number of bytes required to pack the given value.
Definition: bin_grid.F90:330
logical function pmc_mpi_allequal_integer(val)
Returns whether all processors have the same value.
Definition: mpi.F90:1358
subroutine vol_to_lnr(r, f_vol, f_lnr)
Convert a concentration f(vol)d(vol) to f(ln(r))d(ln(r)) where vol = 4/3 pi r^3.
Definition: bin_grid.F90:120
subroutine pmc_nc_read_real_1d(ncid, var, name, must_be_present)
Read a simple real array from a NetCDF file.
Definition: netcdf.F90:218
subroutine pmc_mpi_unpack_integer(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:771
logical function pmc_mpi_allequal_bin_grid(val)
Check whether all processors have the same value.
Definition: bin_grid.F90:401
Wrapper functions for MPI.
Definition: mpi.F90:13
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:743
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
The bin_grid_t structure and associated subroutines.
Definition: bin_grid.F90:9
Reading formatted text input.
Definition: spec_file.F90:43
Wrapper functions for NetCDF. These all take a NetCDF ncid in data mode and return with it again in d...
Definition: netcdf.F90:11
subroutine pmc_mpi_pack_real_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:688
subroutine pmc_mpi_unpack_bin_grid(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: bin_grid.F90:374
character(len=pmc_util_convert_string_len) function real_to_string(val)
Convert a real to a string format.
Definition: util.F90:759
integer function logspace_find(min_x, max_x, n, x)
Find the position of a real number in a 1D logarithmic array.
Definition: util.F90:551
subroutine spec_file_read_real(file, name, var)
Read a real number from a spec file that must have the given name.
Definition: spec_file.F90:579
subroutine linspace(min_x, max_x, x)
Makes a linearly spaced array from min to max.
Definition: util.F90:453
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:102
subroutine bin_grid_copy(bin_grid_from, bin_grid_to)
Copies a bin grid.
Definition: bin_grid.F90:99
subroutine pmc_mpi_pack_bin_grid(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: bin_grid.F90:347
subroutine bin_grid_deallocate(bin_grid)
Frees all memory.
Definition: bin_grid.F90:85