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