PartMC  2.6.1
gas_state.F90
Go to the documentation of this file.
1 ! Copyright (C) 2007-2012 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_state module.
7 
8 !> The gas_state_t structure and associated subroutines.
10 
11  use pmc_util
12  use pmc_spec_file
13  use pmc_gas_data
14  use pmc_env_state
15  use pmc_mpi
16  use pmc_netcdf
17 #ifdef PMC_USE_CAMP
18  use camp_camp_state
19 #endif
20 #ifdef PMC_USE_MPI
21  use mpi
22 #endif
23 
24  !> Current state of the gas mixing ratios in the system.
25  !!
26  !! The gas species are defined by the gas_data_t structure, so that
27  !! \c gas_state%%mix_rat(i) is the current mixing ratio of the gas
28  !! with name \c gas_data%%name(i), etc.
29  !!
30  !! By convention, if gas_state_is_allocated() return \c .false.,
31  !! then the gas_state is treated as zero for all operations on
32  !! it. This will be the case for new \c gas_state_t structures.
34  !> Length n_spec, mixing ratio (ppb).
35  real(kind=dp), allocatable :: mix_rat(:)
36  end type gas_state_t
37 
38 contains
39 
40 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
41 
42  !> Determine whether the \c gas_state is correctly allocated.
43  logical function gas_state_is_allocated(gas_state)
44 
45  !> Gas state to check.
46  type(gas_state_t), intent(in) :: gas_state
47 
48  gas_state_is_allocated = allocated(gas_state%mix_rat)
49 
50  end function gas_state_is_allocated
51 
52 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
53 
54  !> Sets the sizes of the gas state.
55  subroutine gas_state_set_size(gas_state, n_spec)
56 
57  !> Gas state to be allocated.
58  type(gas_state_t), intent(inout) :: gas_state
59  !> Number of species.
60  integer, intent(in) :: n_spec
61 
62  if (allocated(gas_state%mix_rat)) deallocate(gas_state%mix_rat)
63  allocate(gas_state%mix_rat(n_spec))
64  call gas_state_zero(gas_state)
65 
66  end subroutine gas_state_set_size
67 
68 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
69 
70  !> Zeros the state.
71  subroutine gas_state_zero(gas_state)
72 
73  !> Gas state.
74  type(gas_state_t), intent(inout) :: gas_state
75 
76  if (gas_state_is_allocated(gas_state)) then
77  gas_state%mix_rat = 0d0
78  end if
79 
80  end subroutine gas_state_zero
81 
82 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
83 
84  !> Scale a gas state.
85  subroutine gas_state_scale(gas_state, alpha)
86 
87  !> Existing gas state.
88  type(gas_state_t), intent(inout) :: gas_state
89  !> Scale factor.
90  real(kind=dp), intent(in) :: alpha
91 
92  if (gas_state_is_allocated(gas_state)) then
93  gas_state%mix_rat = gas_state%mix_rat * alpha
94  end if
95 
96  end subroutine gas_state_scale
97 
98 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
99 
100  !> Adds the given gas_state_delta.
101  subroutine gas_state_add(gas_state, gas_state_delta)
102 
103  !> Existing gas state.
104  type(gas_state_t), intent(inout) :: gas_state
105  !> Incremental state.
106  type(gas_state_t), intent(in) :: gas_state_delta
107 
108  if (gas_state_is_allocated(gas_state_delta)) then
109  if (gas_state_is_allocated(gas_state)) then
110  gas_state%mix_rat = gas_state%mix_rat + gas_state_delta%mix_rat
111  else
112  gas_state%mix_rat = gas_state_delta%mix_rat
113  end if
114  end if
115 
116  end subroutine gas_state_add
117 
118 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
119 
120  !> Adds the given \c gas_state_delta scaled by \c alpha.
121  !!
122  !! Does gas_state += alpha * gas_state_delta.
123  subroutine gas_state_add_scaled(gas_state, gas_state_delta, alpha)
124 
125  !> Existing gas state.
126  type(gas_state_t), intent(inout) :: gas_state
127  !> Incremental state.
128  type(gas_state_t), intent(in) :: gas_state_delta
129  !> Scale factor.
130  real(kind=dp), intent(in) :: alpha
131 
132  if (gas_state_is_allocated(gas_state_delta)) then
133  if (gas_state_is_allocated(gas_state)) then
134  gas_state%mix_rat = gas_state%mix_rat &
135  + alpha * gas_state_delta%mix_rat
136  else
137  gas_state%mix_rat = alpha * gas_state_delta%mix_rat
138  end if
139  end if
140 
141  end subroutine gas_state_add_scaled
142 
143 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
144 
145  !> Subtracts the given gas_state_delta.
146  subroutine gas_state_sub(gas_state, gas_state_delta)
147 
148  !> Existing gas state.
149  type(gas_state_t), intent(inout) :: gas_state
150  !> Incremental state.
151  type(gas_state_t), intent(in) :: gas_state_delta
152 
153  if (gas_state_is_allocated(gas_state_delta)) then
154  if (gas_state_is_allocated(gas_state)) then
155  gas_state%mix_rat = gas_state%mix_rat - gas_state_delta%mix_rat
156  else
157  gas_state%mix_rat = - gas_state_delta%mix_rat
158  end if
159  end if
160 
161  end subroutine gas_state_sub
162 
163 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
164 
165  !> Set any negative values to zero.
166  subroutine gas_state_ensure_nonnegative(gas_state)
167 
168  !> Gas state.
169  type(gas_state_t) :: gas_state
170 
171  if (gas_state_is_allocated(gas_state)) then
172  gas_state%mix_rat = max(gas_state%mix_rat, 0d0)
173  end if
174 
175  end subroutine gas_state_ensure_nonnegative
176 
177 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
178 
179  !> Convert (mol m^{-3}) to (ppb).
180  subroutine gas_state_mole_dens_to_ppb(gas_state, env_state)
181 
182  !> Gas state.
183  type(gas_state_t), intent(inout) :: gas_state
184  !> Environment state.
185  type(env_state_t), intent(in) :: env_state
186 
187  call gas_state_scale(gas_state, 1d9 / env_state_air_molar_den(env_state))
188 
189  end subroutine gas_state_mole_dens_to_ppb
190 
191 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
192 
193 #ifdef PMC_USE_CAMP
194  !> Set CAMP gas-phase species concentrations
195  subroutine gas_state_set_camp_conc(gas_state, camp_state, gas_data)
196 
197  !> Gas state
198  class(gas_state_t), intent(in) :: gas_state
199  !> CAMP state
200  type(camp_state_t), intent(inout) :: camp_state
201  !> Gas data
202  type(gas_data_t), intent(in) :: gas_data
203 
204  real(kind=dp), parameter :: t_steam = 373.15 ! steam temperature (K)
205  real(kind=dp) :: a, water_vp
206 
207  camp_state%state_var(1:size(gas_state%mix_rat)) = gas_state%mix_rat(:) &
208  / 1000.0d0
209 
210  ! Convert relative humidity (1) to [H2O] (ppm)
211  ! From MOSAIC code - reference to Seinfeld & Pandis page 181
212  ! TODO Figure out how to have consistent RH<->ppm conversions
213  ! (There is only one environmental state for PartMC runs
214  call assert(590005048, associated(camp_state%env_states(1)%val))
215  a = 1.0 - t_steam / camp_state%env_states(1)%val%temp
216  a = (((-0.1299 * a - 0.6445) * a - 1.976) * a + 13.3185) * a
217  water_vp = 101325.0 * exp(a) ! (Pa)
218 
219  camp_state%state_var(gas_data%i_camp_water) = &
220  camp_state%env_states(1)%val%rel_humid * water_vp * 1.0e6 &
221  / camp_state%env_states(1)%val%pressure ! (ppm)
222 
223  end subroutine gas_state_set_camp_conc
224 
225 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
226 
227  !> Get CAMP gas-phase species concentrations
228  subroutine gas_state_get_camp_conc(gas_state, camp_state)
229 
230  !> Gas state
231  class(gas_state_t), intent(inout) :: gas_state
232  !> CAMP state
233  type(camp_state_t), intent(in) :: camp_state
234 
235  gas_state%mix_rat(:) = camp_state%state_var(1:size(gas_state%mix_rat)) &
236  * 1000.0d0
237 
238  end subroutine gas_state_get_camp_conc
239 #endif
240 
241 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
242 
243  !> Determine the current gas_state and rate by interpolating at the
244  !> current time with the lists of gas_states and rates.
245  subroutine gas_state_interp_1d(gas_state_list, time_list, &
246  rate_list, time, gas_state, rate)
247 
248  !> Gas states.
249  type(gas_state_t), intent(in) :: gas_state_list(:)
250  !> Times (s).
251  real(kind=dp), intent(in) :: time_list(size(gas_state_list))
252  !> Rates (s^{-1}).
253  real(kind=dp), intent(in) :: rate_list(size(gas_state_list))
254  !> Current time (s).
255  real(kind=dp), intent(in) :: time
256  !> Current gas state.
257  type(gas_state_t), intent(inout) :: gas_state
258  !> Current rate (s^{-1}).
259  real(kind=dp), intent(out) :: rate
260 
261  integer :: n, p
262  real(kind=dp) :: y, alpha
263 
264  n = size(gas_state_list)
265  p = find_1d(n, time_list, time)
266  if (p == 0) then
267  ! before the start, just use the first state and rate
268  gas_state = gas_state_list(1)
269  rate = rate_list(1)
270  elseif (p == n) then
271  ! after the end, just use the last state and rate
272  gas_state = gas_state_list(n)
273  rate = rate_list(n)
274  else
275  ! in the middle, use the previous state
276  gas_state = gas_state_list(p)
277  rate = rate_list(p)
278  end if
279 
280  end subroutine gas_state_interp_1d
281 
282 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
283 
284  !> Read gas state from the file named on the line read from file.
285  subroutine spec_file_read_gas_state(file, gas_data, gas_state)
286 
287  !> File to read gas state from.
288  type(spec_file_t), intent(inout) :: file
289  !> Gas data.
290  type(gas_data_t), intent(in) :: gas_data
291  !> Gas data to read.
292  type(gas_state_t), intent(inout) :: gas_state
293 
294  integer :: n_species, species, i
295  character(len=SPEC_LINE_MAX_VAR_LEN), allocatable :: species_name(:)
296  real(kind=dp), allocatable :: species_data(:,:)
297 
298  !> \page input_format_gas_state Input File Format: Gas State
299  !!
300  !! A gas state input file must consist of one line per gas
301  !! species, with each line having the species name followed by the
302  !! species mixing ratio in ppb (parts per billion). The valid
303  !! species names are those specfied by the \ref
304  !! input_format_gas_data file, but not all species have to be
305  !! listed. Any missing species will have mixing ratios of
306  !! zero. For example, a gas state file could contain:
307  !! <pre>
308  !! # gas mixing ratio (ppb)
309  !! H2SO4 0
310  !! HNO3 1
311  !! HCl 0.7
312  !! NH3 0.5
313  !! </pre>
314  !!
315  !! See also:
316  !! - \ref spec_file_format --- the input file text format
317  !! - \ref output_format_gas_state --- the corresponding output format
318  !! - \ref input_format_gas_data --- the gas species list and
319  !! material data
320 
321  call spec_file_read_real_named_array(file, 0, species_name, &
322  species_data)
323 
324  ! check the data size
325  n_species = size(species_data, 1)
326  if (.not. ((size(species_data, 2) == 1) .or. (n_species == 0))) then
327  call die_msg(686719840, 'each line in ' // trim(file%name) &
328  // ' must contain exactly one data value')
329  end if
330 
331  ! copy over the data
332  call gas_state_set_size(gas_state, gas_data_n_spec(gas_data))
333  do i = 1,n_species
334  species = gas_data_spec_by_name(gas_data, species_name(i))
335  if (species == 0) then
336  call die_msg(129794076, 'unknown species ' // &
337  trim(species_name(i)) // ' in file ' // trim(file%name))
338  end if
339  gas_state%mix_rat(species) = species_data(i,1)
340  end do
341 
342  end subroutine spec_file_read_gas_state
343 
344 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
345 
346  !> Read an array of gas states with associated times and rates from
347  !> the file named on the line read from the given file.
348  subroutine spec_file_read_gas_states_times_rates(file, gas_data, &
349  times, rates, gas_states)
350 
351  !> Spec file.
352  type(spec_file_t), intent(inout) :: file
353  !> Gas data.
354  type(gas_data_t), intent(in) :: gas_data
355  !> Times (s).
356  real(kind=dp), allocatable :: times(:)
357  !> Rates (s^{-1}).
358  real(kind=dp), allocatable :: rates(:)
359  !> Gas states.
360  type(gas_state_t), allocatable :: gas_states(:)
361 
362  integer :: n_lines, species, i, n_time, i_time
363  character(len=SPEC_LINE_MAX_VAR_LEN), allocatable :: species_name(:)
364  real(kind=dp), allocatable :: species_data(:,:)
365 
366  !> \page input_format_gas_profile Input File Format: Gas Profile
367  !!
368  !! A gas profile input file must consist of three or more
369  !! lines, consisting of:
370  !! - the first line must begin with \c time and should be followed
371  !! by \f$N\f$ space-separated real scalars, giving the times (in
372  !! s after the start of the simulation) of the gas set
373  !! points --- the times must be in increasing order
374  !! - the second line must begin with \c rate and should be
375  !! followed by \f$N\f$ space-separated real scalars, giving the
376  !! values at the corresponding times
377  !! - the third and subsequent lines specify gas species, one
378  !! species per line, with each line beginning with the species
379  !! name and followed by \f$N\f$ space-separated scalars giving
380  !! the gas state of that species at the corresponding times
381  !!
382  !! The units and meanings of the rate and species lines depends on
383  !! the type of gas profile:
384  !! - emissions gas profiles have dimensionless rates that are used
385  !! to scale the species rates and species giving emission rates
386  !! with units of mol/(m^2 s) --- the emission rate is divided by
387  !! the current mixing layer height to give a per-volume emission
388  !! rate
389  !! - background gas profiles have rates with units s^{-1} that are
390  !! dilution rates and species with units of ppb (parts per
391  !! billion) that are the background mixing ratios
392  !!
393  !! The species names must be those specified by the \ref
394  !! input_format_gas_data. Any species not listed are taken to be
395  !! zero.
396  !!
397  !! Between the specified times the gas profile is interpolated
398  !! step-wise and kept constant at its last value. That is, if the
399  !! times are \f$t_i\f$, the rates are \f$r_i\f$, and the gas
400  !! states are \f$g_i\f$ (all with \f$i = 1,\ldots,n\f$), then
401  !! between times \f$t_i\f$ and \f$t_{i+1}\f$ the gas state is
402  !! constant at \f$r_i g_i\f$. Before time \f$t_1\f$ the gas state
403  !! is \f$r_1 g_1\f$, while after time \f$t_n\f$ it is \f$r_n
404  !! g_n\f$.
405  !!
406  !! Example: an emissions gas profile could be:
407  !! <pre>
408  !! time 0 600 1800 # time (in s) after simulation start
409  !! rate 1 0.5 1 # scaling factor
410  !! H2SO4 0 0 0 # emission rate in mol/(m^2 s)
411  !! SO2 4e-9 5.6e-9 5e-9 # emission rate in mol/(m^2 s)
412  !! </pre>
413  !! Here there are no emissions of \f$\rm H_2SO_4\f$, while \f$\rm
414  !! SO_2\f$ starts out being emitted at \f$4\times
415  !! 10^{-9}\rm\ mol\ m^{-2}\ s^{-1}\f$ at the start of the simulation,
416  !! before falling to a rate of \f$2.8\times
417  !! 10^{-9}\rm\ mol\ m^{-2}\ s^{-1}\f$ at 10&nbsp;min (note the scaling
418  !! of 0.5), and then rising again to \f$5\times
419  !! 10^{-9}\rm\ mol\ m^{-2}\ s^{-1}\f$ after 30&nbsp;min. Between
420  !! 0&nbsp;min and 10&nbsp;min the emissions are the same as at
421  !! 0&nbsp;min, while between 10&nbsp;min and 30&nbsp;min they are
422  !! the same as at 10&nbsp;min. After 30&nbsp;min they are held
423  !! constant at their final value.
424  !!
425  !! See also:
426  !! - \ref spec_file_format --- the input file text format
427  !! - \ref input_format_gas_data --- the gas species list and
428  !! material data
429 
430  ! read the data from the file
431  call spec_file_read_real_named_array(file, 0, species_name, &
432  species_data)
433 
434  ! check the data size
435  n_lines = size(species_data, 1)
436  if (n_lines < 2) then
437  call die_msg(291542946, 'insufficient data lines in file ' &
438  // trim(file%name))
439  end if
440  if (trim(species_name(1)) /= 'time') then
441  call die_msg(525127793, 'row 1 in file ' &
442  // trim(file%name) // ' must start with: time')
443  end if
444  if (trim(species_name(2)) /= 'rate') then
445  call die_msg(506981322, 'row 2 in file ' &
446  // trim(file%name) // ' must start with: rate')
447  end if
448  n_time = size(species_data, 2)
449  if (n_time < 1) then
450  call die_msg(398532628, 'each line in file ' &
451  // trim(file%name) // ' must contain at least one data value')
452  end if
453 
454  ! copy over the data
455  times = species_data(1,:)
456  rates = species_data(2,:)
457  if (allocated(gas_states)) deallocate(gas_states)
458  allocate(gas_states(n_time))
459  do i_time = 1,n_time
460  call gas_state_set_size(gas_states(i_time), gas_data_n_spec(gas_data))
461  end do
462  do i = 3,n_lines
463  species = gas_data_spec_by_name(gas_data, species_name(i))
464  if (species == 0) then
465  call die_msg(806500079, 'unknown species ' &
466  // trim(species_name(i)) // ' in file ' &
467  // trim(file%name))
468  end if
469  do i_time = 1,n_time
470  gas_states(i_time)%mix_rat(species) = species_data(i,i_time)
471  end do
472  end do
473 
474  end subroutine spec_file_read_gas_states_times_rates
475 
476 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
477 
478  !> Average val over all processes.
479  subroutine gas_state_mix(val)
480 
481  !> Value to average.
482  type(gas_state_t), intent(inout) :: val
483 
484 #ifdef PMC_USE_MPI
485  type(gas_state_t) :: val_avg
486 
487  call gas_state_set_size(val_avg, size(val%mix_rat))
488  call pmc_mpi_allreduce_average_real_array(val%mix_rat, val_avg%mix_rat)
489  val%mix_rat = val_avg%mix_rat
490 #endif
491 
492  end subroutine gas_state_mix
493 
494 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
495 
496  !> Average val over all processes, with the result only on the root
497  !> process.
498  subroutine gas_state_reduce_avg(val)
499 
500  !> Value to average.
501  type(gas_state_t), intent(inout) :: val
502 
503 #ifdef PMC_USE_MPI
504  type(gas_state_t) :: val_avg
505 
506  call gas_state_set_size(val_avg, size(val%mix_rat))
507  call pmc_mpi_reduce_avg_real_array(val%mix_rat, val_avg%mix_rat)
508  if (pmc_mpi_rank() == 0) then
509  val%mix_rat = val_avg%mix_rat
510  end if
511 #endif
512 
513  end subroutine gas_state_reduce_avg
514 
515 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
516 
517  !> Determines the number of bytes required to pack the given value.
518  integer function pmc_mpi_pack_size_gas_state(val)
519 
520  !> Value to pack.
521  type(gas_state_t), intent(in) :: val
522 
524  + pmc_mpi_pack_size_real_array(val%mix_rat)
525 
526  end function pmc_mpi_pack_size_gas_state
527 
528 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
529 
530  !> Packs the given value into the buffer, advancing position.
531  subroutine pmc_mpi_pack_gas_state(buffer, position, val)
532 
533  !> Memory buffer.
534  character, intent(inout) :: buffer(:)
535  !> Current buffer position.
536  integer, intent(inout) :: position
537  !> Value to pack.
538  type(gas_state_t), intent(in) :: val
539 
540 #ifdef PMC_USE_MPI
541  integer :: prev_position
542 
543  prev_position = position
544  call pmc_mpi_pack_real_array(buffer, position, val%mix_rat)
545  call assert(655827004, &
546  position - prev_position <= pmc_mpi_pack_size_gas_state(val))
547 #endif
548 
549  end subroutine pmc_mpi_pack_gas_state
550 
551 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
552 
553  !> Unpacks the given value from the buffer, advancing position.
554  subroutine pmc_mpi_unpack_gas_state(buffer, position, val)
555 
556  !> Memory buffer.
557  character, intent(inout) :: buffer(:)
558  !> Current buffer position.
559  integer, intent(inout) :: position
560  !> Value to pack.
561  type(gas_state_t), intent(inout) :: val
562 
563 #ifdef PMC_USE_MPI
564  integer :: prev_position
565 
566  prev_position = position
567  call pmc_mpi_unpack_real_array(buffer, position, val%mix_rat)
568  call assert(520815247, &
569  position - prev_position <= pmc_mpi_pack_size_gas_state(val))
570 #endif
571 
572  end subroutine pmc_mpi_unpack_gas_state
573 
574 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
575 
576  !> Computes the average of val across all processes, storing the
577  !> result in val_avg on the root process.
578  subroutine pmc_mpi_reduce_avg_gas_state(val, val_avg)
579 
580  !> Value to average.
581  type(gas_state_t), intent(in) :: val
582  !> Result.
583  type(gas_state_t), intent(inout) :: val_avg
584 
585  call pmc_mpi_reduce_avg_real_array(val%mix_rat, val_avg%mix_rat)
586 
587  end subroutine pmc_mpi_reduce_avg_gas_state
588 
589 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
590 
591  !> Write full state.
592  subroutine gas_state_output_netcdf(gas_state, ncid, gas_data)
593 
594  !> Gas state to write.
595  type(gas_state_t), intent(in) :: gas_state
596  !> NetCDF file ID, in data mode.
597  integer, intent(in) :: ncid
598  !> Gas data.
599  type(gas_data_t), intent(in) :: gas_data
600 
601  integer :: dimid_gas_species
602 
603  !> \page output_format_gas_state Output File Format: Gas State
604  !!
605  !! The gas state uses the \c gas_species NetCDF dimension as specified
606  !! in the \ref output_format_gas_data section.
607  !!
608  !! The gas state NetCDF variables are:
609  !! - \b gas_mixing_ratio (unit ppb, dim \c gas_species): current mixing
610  !! ratios of each gas species
611  !!
612  !! See also:
613  !! - \ref output_format_gas_data --- the gas species list and material
614  !! data output format
615  !! - \ref input_format_gas_state --- the corresponding input format
616 
617  call gas_data_netcdf_dim_gas_species(gas_data, ncid, &
618  dimid_gas_species)
619 
620  call pmc_nc_write_real_1d(ncid, gas_state%mix_rat, &
621  "gas_mixing_ratio", (/ dimid_gas_species /), unit="ppb", &
622  long_name="mixing ratios of gas species")
623 
624  end subroutine gas_state_output_netcdf
625 
626 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
627 
628  !> Read full state.
629  subroutine gas_state_input_netcdf(gas_state, ncid, gas_data)
630 
631  !> Gas state to read.
632  type(gas_state_t), intent(inout) :: gas_state
633  !> NetCDF file ID, in data mode.
634  integer, intent(in) :: ncid
635  !> Gas data.
636  type(gas_data_t), intent(in) :: gas_data
637 
638  call pmc_nc_read_real_1d(ncid, gas_state%mix_rat, "gas_mixing_ratio")
639 
640  end subroutine gas_state_input_netcdf
641 
642 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
643 
644 end module pmc_gas_state
pmc_mpi::pmc_mpi_allreduce_average_real_array
subroutine pmc_mpi_allreduce_average_real_array(val, val_avg)
Computes the average of val across all processes, storing the result in val_avg on all processes.
Definition: mpi.F90:1305
pmc_gas_data::gas_data_t
Constant gas data.
Definition: gas_data.F90:35
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_gas_state::gas_state_mix
subroutine gas_state_mix(val)
Average val over all processes.
Definition: gas_state.F90:480
pmc_gas_state::gas_state_ensure_nonnegative
subroutine gas_state_ensure_nonnegative(gas_state)
Set any negative values to zero.
Definition: gas_state.F90:167
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_gas_state::pmc_mpi_reduce_avg_gas_state
subroutine pmc_mpi_reduce_avg_gas_state(val, val_avg)
Computes the average of val across all processes, storing the result in val_avg on the root process.
Definition: gas_state.F90:579
pmc_gas_state::gas_state_mole_dens_to_ppb
subroutine gas_state_mole_dens_to_ppb(gas_state, env_state)
Convert (mol m^{-3}) to (ppb).
Definition: gas_state.F90:181
pmc_gas_state::gas_state_set_size
subroutine gas_state_set_size(gas_state, n_spec)
Sets the sizes of the gas state.
Definition: gas_state.F90:56
pmc_gas_data
The gas_data_t structure and associated subroutines.
Definition: gas_data.F90:9
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_mpi::pmc_mpi_rank
integer function pmc_mpi_rank()
Returns the rank of the current process.
Definition: mpi.F90:117
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_state::gas_state_add_scaled
subroutine gas_state_add_scaled(gas_state, gas_state_delta, alpha)
Adds the given gas_state_delta scaled by alpha.
Definition: gas_state.F90:124
pmc_env_state::env_state_t
Current environment state.
Definition: env_state.F90:29
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:103
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_spec_file::spec_file_t
An input file with extra data for printing messages.
Definition: spec_file.F90:59
pmc_gas_state
The gas_state_t structure and associated subroutines.
Definition: gas_state.F90:9
pmc_gas_state::pmc_mpi_unpack_gas_state
subroutine pmc_mpi_unpack_gas_state(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: gas_state.F90:555
pmc_gas_state::gas_state_reduce_avg
subroutine gas_state_reduce_avg(val)
Average val over all processes, with the result only on the root process.
Definition: gas_state.F90:499
pmc_env_state::env_state_air_molar_den
real(kind=dp) function env_state_air_molar_den(env_state)
Air molar density (mol m^{-3}).
Definition: env_state.F90:167
pmc_gas_state::gas_state_sub
subroutine gas_state_sub(gas_state, gas_state_delta)
Subtracts the given gas_state_delta.
Definition: gas_state.F90:147
pmc_gas_state::gas_state_scale
subroutine gas_state_scale(gas_state, alpha)
Scale a gas state.
Definition: gas_state.F90:86
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_env_state
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
pmc_gas_state::gas_state_t
Current state of the gas mixing ratios in the system.
Definition: gas_state.F90:33
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_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_state::pmc_mpi_pack_gas_state
subroutine pmc_mpi_pack_gas_state(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: gas_state.F90:532
pmc_gas_state::gas_state_is_allocated
logical function gas_state_is_allocated(gas_state)
Determine whether the gas_state is correctly allocated.
Definition: gas_state.F90:44
pmc_gas_state::gas_state_interp_1d
subroutine gas_state_interp_1d(gas_state_list, time_list, rate_list, time, gas_state, rate)
Determine the current gas_state and rate by interpolating at the current time with the lists of gas_s...
Definition: gas_state.F90:247
pmc_gas_state::gas_state_add
subroutine gas_state_add(gas_state, gas_state_delta)
Adds the given gas_state_delta.
Definition: gas_state.F90:102
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_mpi::pmc_mpi_reduce_avg_real_array
subroutine pmc_mpi_reduce_avg_real_array(val, val_avg)
Computes the average of val across all processes, storing the result in val_avg on the root process.
Definition: mpi.F90:1226
pmc_gas_state::pmc_mpi_pack_size_gas_state
integer function pmc_mpi_pack_size_gas_state(val)
Determines the number of bytes required to pack the given value.
Definition: gas_state.F90:519
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_gas_state::gas_state_input_netcdf
subroutine gas_state_input_netcdf(gas_state, ncid, gas_data)
Read full state.
Definition: gas_state.F90:630
pmc_util::find_1d
integer function find_1d(n, x_vals, x)
Find the position of a real number in an arbitrary 1D array.
Definition: util.F90:592
pmc_gas_state::gas_state_zero
subroutine gas_state_zero(gas_state)
Zeros the state.
Definition: gas_state.F90:72